16 #include <deal.II/base/quadrature_lib.h> 17 #include <deal.II/base/geometry_info.h> 25 DEAL_II_NAMESPACE_OPEN
70 const unsigned int m = (n+1)/2;
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());
120 volatile long double runtime_one = 1.0;
121 const long double tolerance
122 = (runtime_one + long_double_eps != runtime_one
124 std::max (double_eps / 100, long_double_eps * 5)
130 for (
unsigned int i=1; i<=m; ++i)
132 long double z = std::cos(
numbers::PI * (i-.25)/(n+.5));
135 long double p1, p2, p3;
143 for (
unsigned int j=0; j<n; ++j)
147 p1 = ((2.*j+1.)*z*p2-j*p3)/(j+1);
149 pp = n*(z*p1-p2)/(z*z-1);
152 while (std::abs(p1/pp) > tolerance);
158 double w = 1./((1.-z*z)*pp*pp);
170 Assert (n >= 2, ExcNotImplemented());
172 std::vector<long double> points = compute_quadrature_points(n, 1, 1);
173 std::vector<long double> w = compute_quadrature_weights(points, 0, 0);
177 for (
unsigned int i=0; i<points.size(); ++i)
190 const int beta)
const 192 const unsigned int m = q-2;
193 std::vector<long double> x(m);
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());
207 volatile long double runtime_one = 1.0;
208 const long double tolerance
209 = (runtime_one + long_double_eps != runtime_one
211 std::max (double_eps / 100, long_double_eps * 5)
227 for (
unsigned int i=0; i<m; ++i)
228 x[i] = - std::cos( (
long double) (2*i+1)/(2*m) *
numbers::PI );
230 long double r, s, J_x, f, delta;
232 for (
unsigned int k=0; k<m; ++k)
241 for (
unsigned int i=0; i<k; ++i)
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);
249 while (std::fabs(delta) >= tolerance);
255 x.insert(x.begin(), -1.L);
267 const int beta)
const 269 const unsigned int q = x.size();
270 std::vector<long double> w(q);
273 const long double factor = std::pow(2., alpha+beta+1) *
276 ((q-1)*gamma(q)*gamma(alpha+beta+q+1));
277 for (
unsigned int i=0; i<q; ++i)
279 s = JacobiP(x[i], alpha, beta, q-1);
283 w[q-1] *= (alpha + 1);
294 const unsigned int n)
const 298 std::vector<long double> p(n+1);
299 int v, a1, a2, a3, a4;
303 if (n==0)
return p[0];
304 p[1] = ((alpha+beta+2)*x + (alpha-beta))/2;
305 if (n==1)
return p[1];
307 for (
unsigned int i=1; i<=(n-1); ++i)
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);
315 p[i+1] =
static_cast<long double>( (a2 + a3*x)*p[i] - a4*p[i-1])/a1;
325 long double result = n - 1;
326 for (
int i=n-2; i>1; --i)
349 static const double xpts[] = { 0.0, 1.0 };
350 static const double wts[] = { 0.5, 0.5 };
352 for (
unsigned int i=0; i<this->
size(); ++i)
366 static const double xpts[] = { 0.0, 0.5, 1.0 };
367 static const double wts[] = { 1./6., 2./3., 1./6. };
369 for (
unsigned int i=0; i<this->
size(); ++i)
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. };
386 for (
unsigned int i=0; i<this->
size(); ++i)
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.
405 for (
unsigned int i=0; i<this->
size(); ++i)
420 std::vector<double> p=set_quadrature_points(n);
421 std::vector<double> w=set_quadrature_weights(n);
423 for (
unsigned int i=0; i<this->
size(); ++i)
429 this->
weights[i] = revert ? w[n-1-i] : w[i];
439 std::vector<double> points(n);
444 points[0] = 0.3333333333333333;
448 points[0] = 0.1120088061669761;
449 points[1] = 0.6022769081187381;
453 points[0] = 0.06389079308732544;
454 points[1] = 0.3689970637156184;
455 points[2] = 0.766880303938942;
459 points[0] = 0.04144848019938324;
460 points[1] = 0.2452749143206022;
461 points[2] = 0.5561654535602751;
462 points[3] = 0.848982394532986;
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;
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;
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;
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;
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;
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;
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;
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;
560 Assert(
false, ExcNotImplemented());
573 std::vector<double>
weights(n);
581 weights[0] = -0.7185393190303845;
582 weights[1] = -0.2814606809696154;
586 weights[0] = -0.5134045522323634;
587 weights[1] = -0.3919800412014877;
588 weights[2] = -0.0946154065661483;
592 weights[0] =-0.3834640681451353;
593 weights[1] =-0.3868753177747627;
594 weights[2] =-0.1904351269501432;
595 weights[3] =-0.03922548712995894;
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;
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;
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;
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;
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;
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;
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;
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;
693 Assert(
false, ExcNotImplemented());
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] )
729 Assert( (fraction >= 0) && (fraction <= 1),
734 unsigned int ns_offset = (fraction == 1) ? n : 2*n;
736 for (
unsigned int i=0, j=ns_offset; i<n; ++i, ++j)
741 this->
weights[i] = quad1.weight(i)*fraction;
744 if ( (alpha != 1) || (fraction != 1) )
747 this->
weights[j] = -std::log(alpha/fraction)*quad.weight(i)*fraction;
753 this->
weights[i+n] = quad2.weight(i)*(1-fraction);
757 this->
weights[j+n] = -std::log(alpha/(1-fraction))*quad.weight(i)*(1-fraction);
760 if (factor_out_singularity ==
true)
761 for (
unsigned int i=0; i<
size(); ++i)
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;
776 const unsigned int n)
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 ) )
785 if (on_edge && (std::abs( (singularity-
Point<2>(.5, .5)).norm_square()-.5)
788 if (on_vertex)
return (2*n*n);
789 if (on_edge)
return (4*n*n);
796 const bool factor_out_singularity) :
805 std::vector<QGaussOneOverR<2> > quads;
806 std::vector<Point<2> > origins;
815 origins.push_back(
Point<2>(singularity[0],0.));
816 origins.push_back(
Point<2>(0.,singularity[1]));
817 origins.push_back(singularity);
822 unsigned int q_id = 0;
826 for (
unsigned int box=0; box<4; ++box)
829 dist =
Point<2>(std::abs(dist[0]), std::abs(dist[1]));
830 area = dist[0]*dist[1];
832 for (
unsigned int q=0; q<quads[box].size(); ++q, ++q_id)
834 const Point<2> &qp = quads[box].point(q);
837 Point<2>(dist[0]*qp[0], dist[1]*qp[1]);
838 this->
weights[q_id] = quads[box].weight(q)*area;
846 const unsigned int vertex_index,
847 const bool factor_out_singularity) :
856 Assert(vertex_index <4, ExcIndexRange(vertex_index, 0, 4));
862 Assert(gauss.size() == n*n, ExcInternalError());
868 Assert(vertex_index < 4, ExcIndexRange(vertex_index, 0, 4));
883 std::vector<double> &ws = this->
weights;
884 double pi4 = numbers::PI/4;
886 for (
unsigned int q=0; q<gauss.size(); ++q)
888 const Point<2> &gp = gauss.point(q);
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)
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];
903 switch (vertex_index)
910 theta = numbers::PI/2;
913 theta = -numbers::PI/2;
920 double R00 = std::cos(theta), R01 = -std::sin(theta);
921 double R10 = std::sin(theta), R11 = std::cos(theta);
923 if (vertex_index != 0)
924 for (
unsigned int q=0; q<
size(); ++q)
926 double x = ps[q][0]-.5, y = ps[q][1]-.5;
928 ps[q][0] = R00*x + R01*y + .5;
929 ps[q][1] = R10*x + R11*y + .5;
938 std::vector< std::pair<double, Point<dim> > > wp;
939 for (
unsigned int i=0; i<quad.
size(); ++i)
942 sort(wp.begin(), wp.end(), *
this);
943 for (
unsigned int i=0; i<quad.
size(); ++i)
945 this->
weights[i] = wp[i].first;
955 return (a.first < b.first);
1038 const unsigned int n,
const Point<dim> &singularity)
1062 const double eta_bar = singularity[0] * 2. - 1.;
1063 const double eta_star = eta_bar * eta_bar - 1.;
1067 std::vector<double> weights_dummy(
weights.size());
1068 unsigned int cont = 0;
1069 const double tol = 1e-10;
1075 weights_dummy[d-cont] =
weights[d];
1088 weights.resize(weights_dummy.size()-1);
1092 weights[d] = weights_dummy[d];
1096 if (std::abs(eta_star) <= tol)
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)
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)
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);
1117 double J = 3 * ((gamma - gamma_bar) *(gamma - gamma_bar))
1118 / (1 + 3 * gamma_bar * gamma_bar);
1133 std::vector<double> points(n);
1135 for (
unsigned short i=0; i<n; ++i)
1139 points[i] = 1./2.*(1.+std::cos(numbers::PI*(1.+
double(2*i+1)/
double(2*(n-1)+2))));
1150 std::vector<double>
weights(n);
1152 for (
unsigned short i=0; i<n; ++i)
1155 weights[i] = numbers::PI/double(n);
1170 std::vector<double> p=get_quadrature_points(n);
1171 std::vector<double> w=get_quadrature_weights(n);
1173 for (
unsigned int i=0; i<this->
size(); ++i)
1198 std::vector<double> points(n);
1200 for (
unsigned short i=0; i<n; ++i)
1204 if (ep == QGaussRadauChebyshev::left)
1205 points[i] = 1./2.*(1.-std::cos(numbers::PI*(1+2*
double(i)/(2*
double(n-1)+1.))));
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.))));
1222 std::vector<double>
weights(n);
1224 for (
unsigned short i=0; i<n; ++i)
1227 weights[i] = 2.*numbers::PI/double(2*(n-1)+1.);
1228 if (ep==left && i==0)
1230 else if (ep==right && i==(n-1))
1241 QGaussRadauChebyshev<1>::EndPoint ep)
1251 for (
unsigned int i=0; i<this->
size(); ++i)
1284 std::vector<double> points(n);
1286 for (
unsigned short i=0; i<n; ++i)
1290 points[i] = 1./2.*(1.+std::cos(numbers::PI*(1+
double(i)/
double(n-1))));
1301 std::vector<double>
weights(n);
1303 for (
unsigned short i=0; i<n; ++i)
1306 weights[i] = numbers::PI/double((n-1));
1307 if (i==0 || i==(n-1))
1322 Assert(n>1,
ExcMessage(
"Need at least two points for Gauss-Lobatto quadrature rule"));
1326 for (
unsigned int i=0; i<this->
size(); ++i)
1379 DEAL_II_NAMESPACE_CLOSE
QGaussLog(const unsigned int n, const bool revert=false)
std::vector< double > weights
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.
std::vector< double > set_quadrature_points(const unsigned int n) const
#define Assert(cond, exc)
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
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)
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 >)