Reference documentation for deal.II version 8.4.2
quadrature_lib.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/quadrature_lib.h>
17 #include <deal.II/base/geometry_info.h>
18 
19 
20 #include <cmath>
21 #include <limits>
22 #include <algorithm>
23 
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 // please note: for a given dimension, we need the quadrature formulae
29 // for all lower dimensions as well. That is why in this file the check
30 // is for deal_II_dimension >= any_number and not for ==
31 
32 
33 
34 template <>
35 QGauss<0>::QGauss (const unsigned int)
36  :
37  // there are n_q^dim == 1
38  // points
39  Quadrature<0> (1)
40 {
41  // the single quadrature point gets unit
42  // weight
43  this->weights[0] = 1;
44 }
45 
46 
47 
48 template <>
49 QGaussLobatto<0>::QGaussLobatto (const unsigned int)
50  :
51  // there are n_q^dim == 1
52  // points
53  Quadrature<0> (1)
54 {
55  // the single quadrature point gets unit
56  // weight
57  this->weights[0] = 1;
58 }
59 
60 
61 
62 template <>
63 QGauss<1>::QGauss (const unsigned int n)
64  :
65  Quadrature<1> (n)
66 {
67  if (n == 0)
68  return;
69 
70  const unsigned int m = (n+1)/2;
71 
72  // tolerance for the Newton
73  // iteration below. we need to make
74  // it adaptive since on some
75  // machines (for example PowerPC)
76  // long double is the same as
77  // double -- in that case we can
78  // only get to a certain multiple
79  // of the accuracy of double there,
80  // while on other machines we'd
81  // like to go further down
82  //
83  // the situation is complicated by
84  // the fact that even if long
85  // double exists and is described
86  // by std::numeric_limits, we may
87  // not actually get the additional
88  // precision. One case where this
89  // happens is on x86, where one can
90  // set hardware flags that disable
91  // long double precision even for
92  // long double variables. these
93  // flags are not usually set, but
94  // for example matlab sets them and
95  // this then breaks deal.II code
96  // that is run as a subroutine to
97  // matlab...
98  //
99  // a similar situation exists, btw,
100  // when running programs under
101  // valgrind up to and including at
102  // least version 3.3: valgrind's
103  // emulator only supports 64 bit
104  // arithmetic, even for 80 bit long
105  // doubles.
106  const long double
107  long_double_eps = static_cast<long double>(std::numeric_limits<long double>::epsilon()),
108  double_eps = static_cast<long double>(std::numeric_limits<double>::epsilon());
109 
110  // now check whether long double is more
111  // accurate than double, and set
112  // tolerances accordingly. generate a one
113  // that really is generated at run-time
114  // and is not optimized away by the
115  // compiler. that makes sure that the
116  // tolerance is set at run-time with the
117  // current behavior, not at compile-time
118  // (not doing so leads to trouble with
119  // valgrind for example).
120  volatile long double runtime_one = 1.0;
121  const long double tolerance
122  = (runtime_one + long_double_eps != runtime_one
123  ?
124  std::max (double_eps / 100, long_double_eps * 5)
125  :
126  double_eps * 5
127  );
128 
129 
130  for (unsigned int i=1; i<=m; ++i)
131  {
132  long double z = std::cos(numbers::PI * (i-.25)/(n+.5));
133 
134  long double pp;
135  long double p1, p2, p3;
136 
137  // Newton iteration
138  do
139  {
140  // compute L_n (z)
141  p1 = 1.;
142  p2 = 0.;
143  for (unsigned int j=0; j<n; ++j)
144  {
145  p3 = p2;
146  p2 = p1;
147  p1 = ((2.*j+1.)*z*p2-j*p3)/(j+1);
148  }
149  pp = n*(z*p1-p2)/(z*z-1);
150  z = z-p1/pp;
151  }
152  while (std::abs(p1/pp) > tolerance);
153 
154  double x = .5*z;
155  this->quadrature_points[i-1] = Point<1>(.5-x);
156  this->quadrature_points[n-i] = Point<1>(.5+x);
157 
158  double w = 1./((1.-z*z)*pp*pp);
159  this->weights[i-1] = w;
160  this->weights[n-i] = w;
161  }
162 }
163 
164 
165 template <>
166 QGaussLobatto<1>::QGaussLobatto (const unsigned int n)
167  :
168  Quadrature<1> (n)
169 {
170  Assert (n >= 2, ExcNotImplemented());
171 
172  std::vector<long double> points = compute_quadrature_points(n, 1, 1);
173  std::vector<long double> w = compute_quadrature_weights(points, 0, 0);
174 
175  // scale points to the interval
176  // [0.0, 1.0]:
177  for (unsigned int i=0; i<points.size(); ++i)
178  {
179  this->quadrature_points[i] = Point<1>(0.5 + 0.5*static_cast<double>(points[i]));
180  this->weights[i] = 0.5*w[i];
181  }
182 }
183 
184 
185 
186 template <>
187 std::vector<long double> QGaussLobatto<1>::
188 compute_quadrature_points(const unsigned int q,
189  const int alpha,
190  const int beta) const
191 {
192  const unsigned int m = q-2; // no. of inner points
193  std::vector<long double> x(m);
194 
195  // compute quadrature points with
196  // a Newton algorithm.
197 
198  // Set tolerance. See class QGauss
199  // for detailed explanation.
200  const long double
201  long_double_eps = static_cast<long double>(std::numeric_limits<long double>::epsilon()),
202  double_eps = static_cast<long double>(std::numeric_limits<double>::epsilon());
203 
204  // check whether long double is
205  // more accurate than double, and
206  // set tolerances accordingly
207  volatile long double runtime_one = 1.0;
208  const long double tolerance
209  = (runtime_one + long_double_eps != runtime_one
210  ?
211  std::max (double_eps / 100, long_double_eps * 5)
212  :
213  double_eps * 5
214  );
215 
216  // The following implementation
217  // follows closely the one given in
218  // the appendix of the book by
219  // Karniadakis and Sherwin:
220  // Spectral/hp element methods for
221  // computational fluid dynamics
222  // (Oxford University Press, 2005)
223 
224  // we take the zeros of the Chebyshev
225  // polynomial (alpha=beta=-0.5) as
226  // initial values:
227  for (unsigned int i=0; i<m; ++i)
228  x[i] = - std::cos( (long double) (2*i+1)/(2*m) * numbers::PI );
229 
230  long double r, s, J_x, f, delta;
231 
232  for (unsigned int k=0; k<m; ++k)
233  {
234  r = x[k];
235  if (k>0)
236  r = (r + x[k-1])/2;
237 
238  do
239  {
240  s = 0.;
241  for (unsigned int i=0; i<k; ++i)
242  s += 1./(r - x[i]);
243 
244  J_x = 0.5*(alpha + beta + m + 1)*JacobiP(r, alpha+1, beta+1, m-1);
245  f = JacobiP(r, alpha, beta, m);
246  delta = f/(f*s- J_x);
247  r += delta;
248  }
249  while (std::fabs(delta) >= tolerance);
250 
251  x[k] = r;
252  } // for
253 
254  // add boundary points:
255  x.insert(x.begin(), -1.L);
256  x.push_back(+1.L);
257 
258  return x;
259 }
260 
261 
262 
263 template <>
264 std::vector<long double> QGaussLobatto<1>::
265 compute_quadrature_weights(const std::vector<long double> &x,
266  const int alpha,
267  const int beta) const
268 {
269  const unsigned int q = x.size();
270  std::vector<long double> w(q);
271  long double s = 0.L;
272 
273  const long double factor = std::pow(2., alpha+beta+1) *
274  gamma(alpha+q) *
275  gamma(beta+q) /
276  ((q-1)*gamma(q)*gamma(alpha+beta+q+1));
277  for (unsigned int i=0; i<q; ++i)
278  {
279  s = JacobiP(x[i], alpha, beta, q-1);
280  w[i] = factor/(s*s);
281  }
282  w[0] *= (beta + 1);
283  w[q-1] *= (alpha + 1);
284 
285  return w;
286 }
287 
288 
289 
290 template <>
291 long double QGaussLobatto<1>::JacobiP(const long double x,
292  const int alpha,
293  const int beta,
294  const unsigned int n) const
295 {
296  // the Jacobi polynomial is evaluated
297  // using a recursion formula.
298  std::vector<long double> p(n+1);
299  int v, a1, a2, a3, a4;
300 
301  // initial values P_0(x), P_1(x):
302  p[0] = 1.0L;
303  if (n==0) return p[0];
304  p[1] = ((alpha+beta+2)*x + (alpha-beta))/2;
305  if (n==1) return p[1];
306 
307  for (unsigned int i=1; i<=(n-1); ++i)
308  {
309  v = 2*i + alpha + beta;
310  a1 = 2*(i+1)*(i + alpha + beta + 1)*v;
311  a2 = (v + 1)*(alpha*alpha - beta*beta);
312  a3 = v*(v + 1)*(v + 2);
313  a4 = 2*(i+alpha)*(i+beta)*(v + 2);
314 
315  p[i+1] = static_cast<long double>( (a2 + a3*x)*p[i] - a4*p[i-1])/a1;
316  } // for
317  return p[n];
318 }
319 
320 
321 
322 template <>
323 long double QGaussLobatto<1>::gamma(const unsigned int n) const
324 {
325  long double result = n - 1;
326  for (int i=n-2; i>1; --i)
327  result *= i;
328  return result;
329 }
330 
331 
332 
333 template <>
335  :
336  Quadrature<1>(1)
337 {
338  this->quadrature_points[0] = Point<1>(0.5);
339  this->weights[0] = 1.0;
340 }
341 
342 
343 
344 template <>
346  :
347  Quadrature<1> (2)
348 {
349  static const double xpts[] = { 0.0, 1.0 };
350  static const double wts[] = { 0.5, 0.5 };
351 
352  for (unsigned int i=0; i<this->size(); ++i)
353  {
354  this->quadrature_points[i] = Point<1>(xpts[i]);
355  this->weights[i] = wts[i];
356  };
357 }
358 
359 
360 
361 template <>
363  :
364  Quadrature<1> (3)
365 {
366  static const double xpts[] = { 0.0, 0.5, 1.0 };
367  static const double wts[] = { 1./6., 2./3., 1./6. };
368 
369  for (unsigned int i=0; i<this->size(); ++i)
370  {
371  this->quadrature_points[i] = Point<1>(xpts[i]);
372  this->weights[i] = wts[i];
373  };
374 }
375 
376 
377 
378 template <>
380  :
381  Quadrature<1> (5)
382 {
383  static const double xpts[] = { 0.0, .25, .5, .75, 1.0 };
384  static const double wts[] = { 7./90., 32./90., 12./90., 32./90., 7./90. };
385 
386  for (unsigned int i=0; i<this->size(); ++i)
387  {
388  this->quadrature_points[i] = Point<1>(xpts[i]);
389  this->weights[i] = wts[i];
390  };
391 }
392 
393 
394 
395 template <>
397  :
398  Quadrature<1> (7)
399 {
400  static const double xpts[] = { 0.0, 1./6., 1./3., .5, 2./3., 5./6., 1.0 };
401  static const double wts[] = { 41./840., 216./840., 27./840., 272./840.,
402  27./840., 216./840., 41./840.
403  };
404 
405  for (unsigned int i=0; i<this->size(); ++i)
406  {
407  this->quadrature_points[i] = Point<1>(xpts[i]);
408  this->weights[i] = wts[i];
409  };
410 }
411 
412 
413 template <>
414 QGaussLog<1>::QGaussLog(const unsigned int n,
415  const bool revert)
416  :
417  Quadrature<1> (n)
418 {
419 
420  std::vector<double> p=set_quadrature_points(n);
421  std::vector<double> w=set_quadrature_weights(n);
422 
423  for (unsigned int i=0; i<this->size(); ++i)
424  {
425  // Using the change of variables x=1-t, it's possible to show
426  // that int f(x)ln|1-x| = int f(1-t) ln|t|, which implies that
427  // we can use this quadrature formula also with weight ln|1-x|.
428  this->quadrature_points[i] = revert ? Point<1>(1-p[n-1-i]) : Point<1>(p[i]);
429  this->weights[i] = revert ? w[n-1-i] : w[i];
430  }
431 
432 }
433 
434 template <>
435 std::vector<double>
436 QGaussLog<1>::set_quadrature_points(const unsigned int n) const
437 {
438 
439  std::vector<double> points(n);
440 
441  switch (n)
442  {
443  case 1:
444  points[0] = 0.3333333333333333;
445  break;
446 
447  case 2:
448  points[0] = 0.1120088061669761;
449  points[1] = 0.6022769081187381;
450  break;
451 
452  case 3:
453  points[0] = 0.06389079308732544;
454  points[1] = 0.3689970637156184;
455  points[2] = 0.766880303938942;
456  break;
457 
458  case 4:
459  points[0] = 0.04144848019938324;
460  points[1] = 0.2452749143206022;
461  points[2] = 0.5561654535602751;
462  points[3] = 0.848982394532986;
463  break;
464 
465  case 5:
466  points[0] = 0.02913447215197205;
467  points[1] = 0.1739772133208974;
468  points[2] = 0.4117025202849029;
469  points[3] = 0.6773141745828183;
470  points[4] = 0.89477136103101;
471  break;
472 
473  case 6:
474  points[0] = 0.02163400584411693;
475  points[1] = 0.1295833911549506;
476  points[2] = 0.3140204499147661;
477  points[3] = 0.5386572173517997;
478  points[4] = 0.7569153373774084;
479  points[5] = 0.922668851372116;
480  break;
481 
482 
483  case 7:
484  points[0] = 0.0167193554082585;
485  points[1] = 0.100185677915675;
486  points[2] = 0.2462942462079286;
487  points[3] = 0.4334634932570557;
488  points[4] = 0.6323509880476823;
489  points[5] = 0.81111862674023;
490  points[6] = 0.940848166743287;
491  break;
492 
493  case 8:
494  points[0] = 0.01332024416089244;
495  points[1] = 0.07975042901389491;
496  points[2] = 0.1978710293261864;
497  points[3] = 0.354153994351925;
498  points[4] = 0.5294585752348643;
499  points[5] = 0.7018145299391673;
500  points[6] = 0.849379320441094;
501  points[7] = 0.953326450056343;
502  break;
503 
504  case 9:
505  points[0] = 0.01086933608417545;
506  points[1] = 0.06498366633800794;
507  points[2] = 0.1622293980238825;
508  points[3] = 0.2937499039716641;
509  points[4] = 0.4466318819056009;
510  points[5] = 0.6054816627755208;
511  points[6] = 0.7541101371585467;
512  points[7] = 0.877265828834263;
513  points[8] = 0.96225055941096;
514  break;
515 
516  case 10:
517  points[0] = 0.00904263096219963;
518  points[1] = 0.05397126622250072;
519  points[2] = 0.1353118246392511;
520  points[3] = 0.2470524162871565;
521  points[4] = 0.3802125396092744;
522  points[5] = 0.5237923179723384;
523  points[6] = 0.6657752055148032;
524  points[7] = 0.7941904160147613;
525  points[8] = 0.898161091216429;
526  points[9] = 0.9688479887196;
527  break;
528 
529 
530  case 11:
531  points[0] = 0.007643941174637681;
532  points[1] = 0.04554182825657903;
533  points[2] = 0.1145222974551244;
534  points[3] = 0.2103785812270227;
535  points[4] = 0.3266955532217897;
536  points[5] = 0.4554532469286375;
537  points[6] = 0.5876483563573721;
538  points[7] = 0.7139638500230458;
539  points[8] = 0.825453217777127;
540  points[9] = 0.914193921640008;
541  points[10] = 0.973860256264123;
542  break;
543 
544  case 12:
545  points[0] = 0.006548722279080035;
546  points[1] = 0.03894680956045022;
547  points[2] = 0.0981502631060046;
548  points[3] = 0.1811385815906331;
549  points[4] = 0.2832200676673157;
550  points[5] = 0.398434435164983;
551  points[6] = 0.5199526267791299;
552  points[7] = 0.6405109167754819;
553  points[8] = 0.7528650118926111;
554  points[9] = 0.850240024421055;
555  points[10] = 0.926749682988251;
556  points[11] = 0.977756129778486;
557  break;
558 
559  default:
560  Assert(false, ExcNotImplemented());
561  break;
562  }
563 
564  return points;
565 }
566 
567 
568 template <>
569 std::vector<double>
570 QGaussLog<1>::set_quadrature_weights(const unsigned int n) const
571 {
572 
573  std::vector<double> weights(n);
574 
575  switch (n)
576  {
577  case 1:
578  weights[0] = -1.0;
579  break;
580  case 2:
581  weights[0] = -0.7185393190303845;
582  weights[1] = -0.2814606809696154;
583  break;
584 
585  case 3:
586  weights[0] = -0.5134045522323634;
587  weights[1] = -0.3919800412014877;
588  weights[2] = -0.0946154065661483;
589  break;
590 
591  case 4:
592  weights[0] =-0.3834640681451353;
593  weights[1] =-0.3868753177747627;
594  weights[2] =-0.1904351269501432;
595  weights[3] =-0.03922548712995894;
596  break;
597 
598  case 5:
599  weights[0] =-0.2978934717828955;
600  weights[1] =-0.3497762265132236;
601  weights[2] =-0.234488290044052;
602  weights[3] =-0.0989304595166356;
603  weights[4] =-0.01891155214319462;
604  break;
605 
606  case 6:
607  weights[0] = -0.2387636625785478;
608  weights[1] = -0.3082865732739458;
609  weights[2] = -0.2453174265632108;
610  weights[3] = -0.1420087565664786;
611  weights[4] = -0.05545462232488041;
612  weights[5] = -0.01016895869293513;
613  break;
614 
615 
616  case 7:
617  weights[0] = -0.1961693894252476;
618  weights[1] = -0.2703026442472726;
619  weights[2] = -0.239681873007687;
620  weights[3] = -0.1657757748104267;
621  weights[4] = -0.0889432271377365;
622  weights[5] = -0.03319430435645653;
623  weights[6] = -0.005932787015162054;
624  break;
625 
626  case 8:
627  weights[0] = -0.164416604728002;
628  weights[1] = -0.2375256100233057;
629  weights[2] = -0.2268419844319134;
630  weights[3] = -0.1757540790060772;
631  weights[4] = -0.1129240302467932;
632  weights[5] = -0.05787221071771947;
633  weights[6] = -0.02097907374214317;
634  weights[7] = -0.003686407104036044;
635  break;
636 
637  case 9:
638  weights[0] = -0.1400684387481339;
639  weights[1] = -0.2097722052010308;
640  weights[2] = -0.211427149896601;
641  weights[3] = -0.1771562339380667;
642  weights[4] = -0.1277992280331758;
643  weights[5] = -0.07847890261203835;
644  weights[6] = -0.0390225049841783;
645  weights[7] = -0.01386729555074604;
646  weights[8] = -0.002408041036090773;
647  break;
648 
649  case 10:
650  weights[0] = -0.12095513195457;
651  weights[1] = -0.1863635425640733;
652  weights[2] = -0.1956608732777627;
653  weights[3] = -0.1735771421828997;
654  weights[4] = -0.135695672995467;
655  weights[5] = -0.0936467585378491;
656  weights[6] = -0.05578772735275126;
657  weights[7] = -0.02715981089692378;
658  weights[8] = -0.00951518260454442;
659  weights[9] = -0.001638157633217673;
660  break;
661 
662 
663  case 11:
664  weights[0] = -0.1056522560990997;
665  weights[1] = -0.1665716806006314;
666  weights[2] = -0.1805632182877528;
667  weights[3] = -0.1672787367737502;
668  weights[4] = -0.1386970574017174;
669  weights[5] = -0.1038334333650771;
670  weights[6] = -0.06953669788988512;
671  weights[7] = -0.04054160079499477;
672  weights[8] = -0.01943540249522013;
673  weights[9] = -0.006737429326043388;
674  weights[10] = -0.001152486965101561;
675  break;
676 
677  case 12:
678  weights[0] = -0.09319269144393;
679  weights[1] = -0.1497518275763289;
680  weights[2] = -0.166557454364573;
681  weights[3] = -0.1596335594369941;
682  weights[4] = -0.1384248318647479;
683  weights[5] = -0.1100165706360573;
684  weights[6] = -0.07996182177673273;
685  weights[7] = -0.0524069547809709;
686  weights[8] = -0.03007108900074863;
687  weights[9] = -0.01424924540252916;
688  weights[10] = -0.004899924710875609;
689  weights[11] = -0.000834029009809656;
690  break;
691 
692  default:
693  Assert(false, ExcNotImplemented());
694  break;
695  }
696 
697  return weights;
698 
699 }
700 
701 
702 template<>
703 QGaussLogR<1>::QGaussLogR(const unsigned int n,
704  const Point<1> origin,
705  const double alpha,
706  const bool factor_out_singularity) :
707  Quadrature<1>( ( (origin[0] == 0) || (origin[0] == 1) ) ?
708  (alpha == 1 ? n : 2*n ) : 4*n ),
709  fraction( ( (origin[0] == 0) || (origin[0] == 1.) ) ? 1. : origin[0] )
710 {
711  // The three quadrature formulas that make this one up. There are
712  // at most two when the origin is one of the extremes, and there is
713  // only one if the origin is one of the extremes and alpha is
714  // equal to one.
715  //
716  // If alpha is different from one, then we need a correction which
717  // is performed with a standard Gauss quadrature rule on each
718  // segment. This is not needed in the standard case where alpha is
719  // equal to one and the origin is on one of the extremes. We
720  // integrate with weight ln|(x-o)/alpha|. In the easy cases, we
721  // only need n quadrature points. In the most difficult one, we
722  // need 2*n points for the first segment, and 2*n points for the
723  // second segment.
724  QGaussLog<1> quad1(n, origin[0] != 0);
725  QGaussLog<1> quad2(n);
726  QGauss<1> quad(n);
727 
728  // Check that the origin is inside 0,1
729  Assert( (fraction >= 0) && (fraction <= 1),
730  ExcMessage("Origin is outside [0,1]."));
731 
732  // Non singular offset. This is the start of non singular quad
733  // points.
734  unsigned int ns_offset = (fraction == 1) ? n : 2*n;
735 
736  for (unsigned int i=0, j=ns_offset; i<n; ++i, ++j)
737  {
738  // The first i quadrature points are the same as quad1, and
739  // are by default singular.
740  this->quadrature_points[i] = quad1.point(i)*fraction;
741  this->weights[i] = quad1.weight(i)*fraction;
742 
743  // We need to scale with -log|fraction*alpha|
744  if ( (alpha != 1) || (fraction != 1) )
745  {
746  this->quadrature_points[j] = quad.point(i)*fraction;
747  this->weights[j] = -std::log(alpha/fraction)*quad.weight(i)*fraction;
748  }
749  // In case we need the second quadrature as well, do it now.
750  if (fraction != 1)
751  {
752  this->quadrature_points[i+n] = quad2.point(i)*(1-fraction)+Point<1>(fraction);
753  this->weights[i+n] = quad2.weight(i)*(1-fraction);
754 
755  // We need to scale with -log|fraction*alpha|
756  this->quadrature_points[j+n] = quad.point(i)*(1-fraction)+Point<1>(fraction);
757  this->weights[j+n] = -std::log(alpha/(1-fraction))*quad.weight(i)*(1-fraction);
758  }
759  }
760  if (factor_out_singularity == true)
761  for (unsigned int i=0; i<size(); ++i)
762  {
763  Assert( this->quadrature_points[i] != origin,
764  ExcMessage("The singularity cannot be on a Gauss point of the same order!") );
765  double denominator = std::log(std::abs( (this->quadrature_points[i]-origin)[0] )/alpha);
766  Assert( denominator != 0.0,
767  ExcMessage("The quadrature formula you are using does not allow to "
768  "factor out the singularity, which is zero at one point."));
769  this->weights[i] /= denominator;
770  }
771 }
772 
773 
774 template<>
775 unsigned int QGaussOneOverR<2>::quad_size(const Point<2> singularity,
776  const unsigned int n)
777 {
778  double eps=1e-8;
779  bool on_edge=false;
780  bool on_vertex=false;
781  for (unsigned int i=0; i<2; ++i)
782  if ( ( std::abs(singularity[i] ) < eps ) ||
783  ( std::abs(singularity[i]-1) < eps ) )
784  on_edge = true;
785  if (on_edge && (std::abs( (singularity-Point<2>(.5, .5)).norm_square()-.5)
786  < eps) )
787  on_vertex = true;
788  if (on_vertex) return (2*n*n);
789  if (on_edge) return (4*n*n);
790  return (8*n*n);
791 }
792 
793 template<>
794 QGaussOneOverR<2>::QGaussOneOverR(const unsigned int n,
795  const Point<2> singularity,
796  const bool factor_out_singularity) :
797  Quadrature<2>(quad_size(singularity, n))
798 {
799  // We treat all the cases in the
800  // same way. Split the element in 4
801  // pieces, measure the area, if
802  // it's relevant, add the
803  // quadrature connected to that
804  // singularity.
805  std::vector<QGaussOneOverR<2> > quads;
806  std::vector<Point<2> > origins;
807  // Id of the corner with a
808  // singularity
809  quads.push_back(QGaussOneOverR(n, 3, factor_out_singularity));
810  quads.push_back(QGaussOneOverR(n, 2, factor_out_singularity));
811  quads.push_back(QGaussOneOverR(n, 1, factor_out_singularity));
812  quads.push_back(QGaussOneOverR(n, 0, factor_out_singularity));
813 
814  origins.push_back(Point<2>(0.,0.));
815  origins.push_back(Point<2>(singularity[0],0.));
816  origins.push_back(Point<2>(0.,singularity[1]));
817  origins.push_back(singularity);
818 
819  // Lexicographical ordering.
820 
821  double eps = 1e-8;
822  unsigned int q_id = 0; // Current quad point index.
823  double area = 0;
824  Tensor<1,2> dist;
825 
826  for (unsigned int box=0; box<4; ++box)
827  {
828  dist = (singularity-GeometryInfo<2>::unit_cell_vertex(box));
829  dist = Point<2>(std::abs(dist[0]), std::abs(dist[1]));
830  area = dist[0]*dist[1];
831  if (area > eps)
832  for (unsigned int q=0; q<quads[box].size(); ++q, ++q_id)
833  {
834  const Point<2> &qp = quads[box].point(q);
835  this->quadrature_points[q_id] =
836  origins[box]+
837  Point<2>(dist[0]*qp[0], dist[1]*qp[1]);
838  this->weights[q_id] = quads[box].weight(q)*area;
839  }
840  }
841 }
842 
843 
844 template<>
845 QGaussOneOverR<2>::QGaussOneOverR(const unsigned int n,
846  const unsigned int vertex_index,
847  const bool factor_out_singularity) :
848  Quadrature<2>(2*n *n)
849 {
850  // This version of the constructor
851  // works only for the 4
852  // vertices. If you need a more
853  // general one, you should use the
854  // one with the Point<2> in the
855  // constructor.
856  Assert(vertex_index <4, ExcIndexRange(vertex_index, 0, 4));
857 
858  // Start with the gauss quadrature formula on the (u,v) reference
859  // element.
860  QGauss<2> gauss(n);
861 
862  Assert(gauss.size() == n*n, ExcInternalError());
863 
864  // For the moment we only implemented this for the vertices of a
865  // quadrilateral. We are planning to do this also for the support
866  // points of arbitrary FE_Q elements, to allow the use of this
867  // class in boundary element programs with higher order mappings.
868  Assert(vertex_index < 4, ExcIndexRange(vertex_index, 0, 4));
869 
870  // We create only the first one. All other pieces are rotation of
871  // this one.
872  // In this case the transformation is
873  //
874  // (x,y) = (u, u tan(pi/4 v))
875  //
876  // with Jacobian
877  //
878  // J = pi/4 R / cos(pi/4 v)
879  //
880  // And we get rid of R to take into account the singularity,
881  // unless specified differently in the constructor.
882  std::vector<Point<2> > &ps = this->quadrature_points;
883  std::vector<double> &ws = this->weights;
884  double pi4 = numbers::PI/4;
885 
886  for (unsigned int q=0; q<gauss.size(); ++q)
887  {
888  const Point<2> &gp = gauss.point(q);
889  ps[q][0] = gp[0];
890  ps[q][1] = gp[0] * std::tan(pi4*gp[1]);
891  ws[q] = gauss.weight(q)*pi4/std::cos(pi4 *gp[1]);
892  if (factor_out_singularity)
893  ws[q] *= (ps[q]-GeometryInfo<2>::unit_cell_vertex(0)).norm();
894  // The other half of the quadrilateral is symmetric with
895  // respect to xy plane.
896  ws[gauss.size()+q] = ws[q];
897  ps[gauss.size()+q][0] = ps[q][1];
898  ps[gauss.size()+q][1] = ps[q][0];
899  }
900 
901  // Now we distribute these vertices in the correct manner
902  double theta = 0;
903  switch (vertex_index)
904  {
905  case 0:
906  theta = 0;
907  break;
908  case 1:
909  //
910  theta = numbers::PI/2;
911  break;
912  case 2:
913  theta = -numbers::PI/2;
914  break;
915  case 3:
916  theta = numbers::PI;
917  break;
918  }
919 
920  double R00 = std::cos(theta), R01 = -std::sin(theta);
921  double R10 = std::sin(theta), R11 = std::cos(theta);
922 
923  if (vertex_index != 0)
924  for (unsigned int q=0; q<size(); ++q)
925  {
926  double x = ps[q][0]-.5, y = ps[q][1]-.5;
927 
928  ps[q][0] = R00*x + R01*y + .5;
929  ps[q][1] = R10*x + R11*y + .5;
930  }
931 }
932 
933 
934 template <int dim>
936  Quadrature<dim>(quad.size())
937 {
938  std::vector< std::pair<double, Point<dim> > > wp;
939  for (unsigned int i=0; i<quad.size(); ++i)
940  wp.push_back(std::pair<double, Point<dim> >(quad.weight(i),
941  quad.point(i)));
942  sort(wp.begin(), wp.end(), *this);
943  for (unsigned int i=0; i<quad.size(); ++i)
944  {
945  this->weights[i] = wp[i].first;
946  this->quadrature_points[i] = wp[i].second;
947  }
948 }
949 
950 
951 template <int dim>
952 bool QSorted<dim>::operator()(const std::pair<double, Point<dim> > &a,
953  const std::pair<double, Point<dim> > &b)
954 {
955  return (a.first < b.first);
956 }
957 
958 
959 // construct the quadrature formulae in higher dimensions by
960 // tensor product of lower dimensions
961 
962 template <int dim>
963 QGauss<dim>::QGauss (const unsigned int n)
964  : Quadrature<dim> (QGauss<dim-1>(n), QGauss<1>(n))
965 {}
966 
967 
968 
969 template <int dim>
970 QGaussLobatto<dim>::QGaussLobatto (const unsigned int n)
971  : Quadrature<dim> (QGaussLobatto<dim-1>(n), QGaussLobatto<1>(n))
972 {}
973 
974 
975 
976 template <int dim>
978  :
980 {}
981 
982 
983 
984 template <int dim>
986  :
987  Quadrature<dim> (QTrapez<dim-1>(), QTrapez<1>())
988 {}
989 
990 
991 
992 template <int dim>
994  :
995  Quadrature<dim> (QSimpson<dim-1>(), QSimpson<1>())
996 {}
997 
998 
999 
1000 template <int dim>
1002  :
1003  Quadrature<dim> (QMilne<dim-1>(), QMilne<1>())
1004 {}
1005 
1006 
1007 template <int dim>
1009  :
1010  Quadrature<dim> (QWeddle<dim-1>(), QWeddle<1>())
1011 {}
1012 
1013 template <int dim>
1015  const Quadrature<1> &base_quad, const Point<dim> &singularity)
1016  :
1022  Quadrature<dim>(
1023  dim == 2 ?
1024  QAnisotropic<dim>(
1025  QTelles<1>(base_quad, Point<1>(singularity[0])),
1026  QTelles<1>(base_quad, Point<1>(singularity[1]))) :
1027  dim == 3 ?
1028  QAnisotropic<dim>(
1029  QTelles<1>(base_quad, Point<1>(singularity[0])),
1030  QTelles<1>(base_quad, Point<1>(singularity[1])),
1031  QTelles<1>(base_quad, Point<1>(singularity[2]))) :
1032  Quadrature<dim>())
1033 {
1034 }
1035 
1036 template <int dim>
1038  const unsigned int n, const Point<dim> &singularity)
1039  :
1044  Quadrature<dim>(QTelles<dim>(QGauss<1>(n), singularity))
1045 {}
1046 
1047 
1048 
1049 template <>
1051  const Quadrature<1> &base_quad, const Point<1> &singularity)
1052  :
1056  Quadrature<1>(base_quad)
1057 {
1062  const double eta_bar = singularity[0] * 2. - 1.;
1063  const double eta_star = eta_bar * eta_bar - 1.;
1064  double gamma_bar;
1065 
1066  std::vector<Point<1> > quadrature_points_dummy(quadrature_points.size());
1067  std::vector<double> weights_dummy(weights.size());
1068  unsigned int cont = 0;
1069  const double tol = 1e-10;
1070  for (unsigned int d = 0; d < quadrature_points.size(); ++d)
1071  {
1072  if (std::abs(quadrature_points[d][0] - singularity[0]) > tol)
1073  {
1074  quadrature_points_dummy[d-cont] = quadrature_points[d];
1075  weights_dummy[d-cont] = weights[d];
1076  }
1077  else
1078  {
1079  // We need to remove the singularity point from the quadrature point
1080  // list. To do so we use the variable cont.
1081  cont = 1;
1082  }
1083 
1084  }
1085  if (cont == 1)
1086  {
1087  quadrature_points.resize(quadrature_points_dummy.size()-1);
1088  weights.resize(weights_dummy.size()-1);
1089  for (unsigned int d = 0; d < quadrature_points.size()-1; ++d)
1090  {
1091  quadrature_points[d] = quadrature_points_dummy[d];
1092  weights[d] = weights_dummy[d];
1093  }
1094  }
1095  // We need to check if the singularity is at the boundary of the interval.
1096  if (std::abs(eta_star) <= tol)
1097  {
1098  gamma_bar = std::pow((eta_bar * eta_star + std::abs(eta_star)),1.0 / 3.0)
1099  + std::pow((eta_bar * eta_star - std::abs(eta_star)), 1.0 / 3.0)
1100  + eta_bar;
1101  }
1102  else
1103  {
1104  gamma_bar = (eta_bar * eta_star + std::abs(eta_star))/std::abs(eta_bar * eta_star + std::abs(eta_star))*
1105  std::pow(std::abs(eta_bar * eta_star + std::abs(eta_star)),1.0 / 3.0)
1106  + (eta_bar * eta_star - std::abs(eta_star))/std::abs(eta_bar * eta_star - std::abs(eta_star))*
1107  std::pow(std::abs(eta_bar * eta_star - std::abs(eta_star)), 1.0 / 3.0)
1108  + eta_bar;
1109  }
1110  for (unsigned int q = 0; q < quadrature_points.size(); ++q)
1111  {
1112  double gamma = quadrature_points[q][0] * 2 - 1;
1113  double eta = (std::pow(gamma - gamma_bar, 3.0)
1114  + gamma_bar * (gamma_bar * gamma_bar + 3))
1115  / (1 + 3 * gamma_bar * gamma_bar);
1116 
1117  double J = 3 * ((gamma - gamma_bar) *(gamma - gamma_bar))
1118  / (1 + 3 * gamma_bar * gamma_bar);
1119 
1120  quadrature_points[q][0] = (eta + 1) / 2.0;
1121  weights[q] = J * weights[q];
1122 
1123  }
1124 }
1125 
1126 
1127 
1128 template <>
1129 std::vector<double>
1130 QGaussChebyshev<1>::get_quadrature_points(const unsigned int n)
1131 {
1132 
1133  std::vector<double> points(n);
1134  // n point quadrature: index from 0 to n-1
1135  for (unsigned short i=0; i<n; ++i)
1136  // would be cos((2i+1)Pi/(2N+2))
1137  // put + Pi so we start from the smallest point
1138  // then map from [-1,1] to [0,1]
1139  points[i] = 1./2.*(1.+std::cos(numbers::PI*(1.+double(2*i+1)/double(2*(n-1)+2))));
1140 
1141  return points;
1142 }
1143 
1144 
1145 template <>
1146 std::vector<double>
1147 QGaussChebyshev<1>::get_quadrature_weights(const unsigned int n)
1148 {
1149 
1150  std::vector<double> weights(n);
1151 
1152  for (unsigned short i=0; i<n; ++i)
1153  {
1154  // same weights as on [-1,1]
1155  weights[i] = numbers::PI/double(n);
1156  }
1157 
1158  return weights;
1159 
1160 }
1161 
1162 
1163 template <>
1164 QGaussChebyshev<1>::QGaussChebyshev(const unsigned int n)
1165  :
1166  Quadrature<1> (n)
1167 {
1168 
1169  Assert(n>0,ExcMessage("Need at least one point for the quadrature rule"));
1170  std::vector<double> p=get_quadrature_points(n);
1171  std::vector<double> w=get_quadrature_weights(n);
1172 
1173  for (unsigned int i=0; i<this->size(); ++i)
1174  {
1175  this->quadrature_points[i] = Point<1>(p[i]);
1176  this->weights[i] = w[i];
1177  }
1178 
1179 }
1180 
1181 
1182 template <int dim>
1184  :
1185  Quadrature<dim> (QGaussChebyshev<dim-1>(n), QGaussChebyshev<1>(n))
1186 {}
1187 
1188 
1189 
1190 
1191 
1192 template <>
1193 std::vector<double>
1195  EndPoint ep)
1196 {
1197 
1198  std::vector<double> points(n);
1199  // n point quadrature: index from 0 to n-1
1200  for (unsigned short i=0; i<n; ++i)
1201  // would be -cos(2i Pi/(2N+1))
1202  // put + Pi so we start from the smallest point
1203  // then map from [-1,1] to [0,1]
1204  if (ep == QGaussRadauChebyshev::left)
1205  points[i] = 1./2.*(1.-std::cos(numbers::PI*(1+2*double(i)/(2*double(n-1)+1.))));
1206  else
1207  {
1208  Assert(ep==QGaussRadauChebyshev::right,ExcInvalidConstructorCall());
1209  points[i] = 1./2.*(1.-std::cos(numbers::PI*(2*double(n-1-i)/(2*double(n-1)+1.))));
1210  }
1211 
1212  return points;
1213 }
1214 
1215 
1216 template <>
1217 std::vector<double>
1219  EndPoint ep)
1220 {
1221 
1222  std::vector<double> weights(n);
1223 
1224  for (unsigned short i=0; i<n; ++i)
1225  {
1226  // same weights as on [-1,1]
1227  weights[i] = 2.*numbers::PI/double(2*(n-1)+1.);
1228  if (ep==left && i==0)
1229  weights[i] /= 2.;
1230  else if (ep==right && i==(n-1))
1231  weights[i] /= 2.;
1232  }
1233 
1234  return weights;
1235 
1236 }
1237 
1238 
1239 template <>
1241  QGaussRadauChebyshev<1>::EndPoint ep)
1242  :
1243  Quadrature<1> (n),
1244  ep (ep)
1245 {
1246 
1247  Assert(n>0,ExcMessage("Need at least one point for quadrature rules"));
1248  std::vector<double> p=get_quadrature_points(n,ep);
1249  std::vector<double> w=get_quadrature_weights(n,ep);
1250 
1251  for (unsigned int i=0; i<this->size(); ++i)
1252  {
1253  this->quadrature_points[i] = Point<1>(p[i]);
1254  this->weights[i] = w[i];
1255  }
1256 }
1257 
1258 
1259 template <>
1260 QGaussRadauChebyshev<2>::QGaussRadauChebyshev (const unsigned int n,
1261  EndPoint ep)
1262  :
1263  Quadrature<2> (QGaussRadauChebyshev<1>(n, static_cast<QGaussRadauChebyshev<1>::EndPoint>(ep)),
1264  QGaussRadauChebyshev<1>(n, static_cast<QGaussRadauChebyshev<1>::EndPoint>(ep))),
1265  ep (ep)
1266 {}
1267 
1268 
1269 template <int dim>
1271  EndPoint ep)
1272  :
1273  Quadrature<dim> (QGaussRadauChebyshev<dim-1>(n,static_cast<typename QGaussRadauChebyshev<dim-1>::EndPoint>(ep)),
1274  QGaussRadauChebyshev<1>(n,static_cast<QGaussRadauChebyshev<1>::EndPoint>(ep))),
1275  ep (ep)
1276 {}
1277 
1278 
1279 template <>
1280 std::vector<double>
1282 {
1283 
1284  std::vector<double> points(n);
1285  // n point quadrature: index from 0 to n-1
1286  for (unsigned short i=0; i<n; ++i)
1287  // would be cos(i Pi/N)
1288  // put + Pi so we start from the smallest point
1289  // then map from [-1,1] to [0,1]
1290  points[i] = 1./2.*(1.+std::cos(numbers::PI*(1+double(i)/double(n-1))));
1291 
1292  return points;
1293 }
1294 
1295 
1296 template <>
1297 std::vector<double>
1299 {
1300 
1301  std::vector<double> weights(n);
1302 
1303  for (unsigned short i=0; i<n; ++i)
1304  {
1305  // same weights as on [-1,1]
1306  weights[i] = numbers::PI/double((n-1));
1307  if (i==0 || i==(n-1))
1308  weights[i] /= 2.;
1309  }
1310 
1311  return weights;
1312 
1313 }
1314 
1315 
1316 template <>
1318  :
1319  Quadrature<1> (n)
1320 {
1321 
1322  Assert(n>1,ExcMessage("Need at least two points for Gauss-Lobatto quadrature rule"));
1323  std::vector<double> p=get_quadrature_points(n);
1324  std::vector<double> w=get_quadrature_weights(n);
1325 
1326  for (unsigned int i=0; i<this->size(); ++i)
1327  {
1328  this->quadrature_points[i] = Point<1>(p[i]);
1329  this->weights[i] = w[i];
1330  }
1331 
1332 }
1333 
1334 
1335 template <int dim>
1337  :
1339 {}
1340 
1341 // explicit specialization
1342 // note that 1d formulae are specialized by implementation above
1343 template class QGauss<2>;
1344 template class QGaussLobatto<2>;
1345 template class QMidpoint<2>;
1346 template class QTrapez<2>;
1347 template class QSimpson<2>;
1348 template class QMilne<2>;
1349 template class QWeddle<2>;
1350 
1351 template class QGauss<3>;
1352 template class QGaussLobatto<3>;
1353 template class QMidpoint<3>;
1354 template class QTrapez<3>;
1355 template class QSimpson<3>;
1356 template class QMilne<3>;
1357 template class QWeddle<3>;
1358 
1359 template class QSorted<1>;
1360 template class QSorted<2>;
1361 template class QSorted<3>;
1362 
1363 template class QTelles<1> ;
1364 template class QTelles<2> ;
1365 template class QTelles<3> ;
1366 
1367 template class QGaussChebyshev<1>;
1368 template class QGaussChebyshev<2>;
1369 template class QGaussChebyshev<3>;
1370 
1371 template class QGaussRadauChebyshev<1>;
1372 template class QGaussRadauChebyshev<2>;
1373 template class QGaussRadauChebyshev<3>;
1374 
1375 template class QGaussLobattoChebyshev<1>;
1376 template class QGaussLobattoChebyshev<2>;
1377 template class QGaussLobattoChebyshev<3>;
1378 
1379 DEAL_II_NAMESPACE_CLOSE
QGaussLog(const unsigned int n, const bool revert=false)
std::vector< double > weights
Definition: quadrature.h:223
QGaussOneOverR(const unsigned int n, const Point< dim > singularity, const bool factor_out_singular_weight=false)
bool operator()(const std::pair< double, Point< dim > > &a, const std::pair< double, Point< dim > > &b)
static std::vector< double > get_quadrature_weights(const unsigned int n)
Computes the weights of the quadrature formula.
long double gamma(const unsigned int n) const
static std::vector< double > get_quadrature_weights(const unsigned int n, EndPoint ep)
Computes the weights of the quadrature formula.
std::vector< long double > compute_quadrature_points(const unsigned int q, const int alpha, const int beta) const
::ExceptionBase & ExcMessage(std::string arg1)
QGaussChebyshev(const unsigned int n)
Generate a formula with n quadrature points.
const Point< dim > & point(const unsigned int i) const
static Point< dim > unit_cell_vertex(const unsigned int vertex)
QGauss(const unsigned int n)
static std::vector< double > get_quadrature_weights(const unsigned int n)
Computes the weights of the quadrature formula.
Definition: point.h:89
static const double PI
Definition: numbers.h:74
std::vector< double > set_quadrature_points(const unsigned int n) const
#define Assert(cond, exc)
Definition: exceptions.h:294
std::vector< long double > compute_quadrature_weights(const std::vector< long double > &x, const int alpha, const int beta) const
QGaussRadauChebyshev(const unsigned int n, EndPoint ep=QGaussRadauChebyshev::left)
Generate a formula with n quadrature points.
static std::vector< double > get_quadrature_points(const unsigned int n)
Computes the points of the quadrature formula.
QGaussLobatto(const unsigned int n)
std::vector< double > set_quadrature_weights(const unsigned int n) const
long double JacobiP(const long double x, const int alpha, const int beta, const unsigned int n) const
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:217
unsigned int size() const
QTelles(const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
static std::vector< double > get_quadrature_points(const unsigned int n, EndPoint ep)
Computes the points of the quadrature formula.
double norm(const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Du)
Definition: divergence.h:604
Definition: mpi.h:48
static std::vector< double > get_quadrature_points(const unsigned int n)
Computes the points of the quadrature formula.
QGaussLobattoChebyshev(const unsigned int n)
Generate a formula with n quadrature points.
QGaussLogR(const unsigned int n, const Point< dim > x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false)
double weight(const unsigned int i) const
static unsigned int quad_size(const Point< dim > singularity, const unsigned int n)
QSorted(const Quadrature< dim >)