76 std::vector<long double> points =
77 Polynomials::jacobi_polynomial_roots<long double>(n, 0, 0);
79 for (
unsigned int i = 0; i < (points.size() + 1) / 2; ++i)
85 const long double pp =
88 const long double x = -1. + 2. * points[i];
89 const double w = 1. / ((1. - x * x) * pp * pp);
106 long double result = n - 1;
107 for (
int i = n - 2; i > 1; --i)
121 std::vector<long double>
126 const unsigned int q = x.size();
127 std::vector<long double> w(q);
129 const long double factor =
131 ((q - 1) *
gamma(q) *
gamma(alpha + beta + q + 1));
132 for (
unsigned int i = 0; i < q; ++i)
134 const long double s =
136 w[i] = factor / (s * s);
139 w[q - 1] *= (alpha + 1);
154 std::vector<long double> points =
155 Polynomials::jacobi_polynomial_roots<long double>(n - 2, 1, 1);
156 points.insert(points.begin(), 0);
157 points.push_back(1.);
158 std::vector<long double>
w =
162 for (
unsigned int i = 0; i < points.size(); ++i)
165 this->weights[i] = 0.5 *
w[i];
185 static const double xpts[] = {0.0, 1.0};
186 static const double wts[] = {0.5, 0.5};
188 for (
unsigned int i = 0; i < this->
size(); ++i)
201 static const double xpts[] = {0.0, 0.5, 1.0};
202 static const double wts[] = {1. / 6., 2. / 3., 1. / 6.};
204 for (
unsigned int i = 0; i < this->
size(); ++i)
217 static const double xpts[] = {0.0, .25, .5, .75, 1.0};
218 static const double wts[] = {
219 7. / 90., 32. / 90., 12. / 90., 32. / 90., 7. / 90.};
221 for (
unsigned int i = 0; i < this->
size(); ++i)
234 static const double xpts[] = {
235 0.0, 1. / 6., 1. / 3., .5, 2. / 3., 5. / 6., 1.0};
236 static const double wts[] = {41. / 840.,
244 for (
unsigned int i = 0; i < this->
size(); ++i)
259 for (
unsigned int i = 0; i < this->
size(); ++i)
266 this->
weights[i] = revert ? w[n - 1 - i] : w[i];
274 std::vector<double> q_points(n);
279 q_points[0] = 0.3333333333333333;
283 q_points[0] = 0.1120088061669761;
284 q_points[1] = 0.6022769081187381;
288 q_points[0] = 0.06389079308732544;
289 q_points[1] = 0.3689970637156184;
290 q_points[2] = 0.766880303938942;
294 q_points[0] = 0.04144848019938324;
295 q_points[1] = 0.2452749143206022;
296 q_points[2] = 0.5561654535602751;
297 q_points[3] = 0.848982394532986;
301 q_points[0] = 0.02913447215197205;
302 q_points[1] = 0.1739772133208974;
303 q_points[2] = 0.4117025202849029;
304 q_points[3] = 0.6773141745828183;
305 q_points[4] = 0.89477136103101;
309 q_points[0] = 0.02163400584411693;
310 q_points[1] = 0.1295833911549506;
311 q_points[2] = 0.3140204499147661;
312 q_points[3] = 0.5386572173517997;
313 q_points[4] = 0.7569153373774084;
314 q_points[5] = 0.922668851372116;
319 q_points[0] = 0.0167193554082585;
320 q_points[1] = 0.100185677915675;
321 q_points[2] = 0.2462942462079286;
322 q_points[3] = 0.4334634932570557;
323 q_points[4] = 0.6323509880476823;
324 q_points[5] = 0.81111862674023;
325 q_points[6] = 0.940848166743287;
329 q_points[0] = 0.01332024416089244;
330 q_points[1] = 0.07975042901389491;
331 q_points[2] = 0.1978710293261864;
332 q_points[3] = 0.354153994351925;
333 q_points[4] = 0.5294585752348643;
334 q_points[5] = 0.7018145299391673;
335 q_points[6] = 0.849379320441094;
336 q_points[7] = 0.953326450056343;
340 q_points[0] = 0.01086933608417545;
341 q_points[1] = 0.06498366633800794;
342 q_points[2] = 0.1622293980238825;
343 q_points[3] = 0.2937499039716641;
344 q_points[4] = 0.4466318819056009;
345 q_points[5] = 0.6054816627755208;
346 q_points[6] = 0.7541101371585467;
347 q_points[7] = 0.877265828834263;
348 q_points[8] = 0.96225055941096;
352 q_points[0] = 0.00904263096219963;
353 q_points[1] = 0.05397126622250072;
354 q_points[2] = 0.1353118246392511;
355 q_points[3] = 0.2470524162871565;
356 q_points[4] = 0.3802125396092744;
357 q_points[5] = 0.5237923179723384;
358 q_points[6] = 0.6657752055148032;
359 q_points[7] = 0.7941904160147613;
360 q_points[8] = 0.898161091216429;
361 q_points[9] = 0.9688479887196;
366 q_points[0] = 0.007643941174637681;
367 q_points[1] = 0.04554182825657903;
368 q_points[2] = 0.1145222974551244;
369 q_points[3] = 0.2103785812270227;
370 q_points[4] = 0.3266955532217897;
371 q_points[5] = 0.4554532469286375;
372 q_points[6] = 0.5876483563573721;
373 q_points[7] = 0.7139638500230458;
374 q_points[8] = 0.825453217777127;
375 q_points[9] = 0.914193921640008;
376 q_points[10] = 0.973860256264123;
380 q_points[0] = 0.006548722279080035;
381 q_points[1] = 0.03894680956045022;
382 q_points[2] = 0.0981502631060046;
383 q_points[3] = 0.1811385815906331;
384 q_points[4] = 0.2832200676673157;
385 q_points[5] = 0.398434435164983;
386 q_points[6] = 0.5199526267791299;
387 q_points[7] = 0.6405109167754819;
388 q_points[8] = 0.7528650118926111;
389 q_points[9] = 0.850240024421055;
390 q_points[10] = 0.926749682988251;
391 q_points[11] = 0.977756129778486;
407 std::vector<double> quadrature_weights(n);
412 quadrature_weights[0] = -1.0;
415 quadrature_weights[0] = -0.7185393190303845;
416 quadrature_weights[1] = -0.2814606809696154;
420 quadrature_weights[0] = -0.5134045522323634;
421 quadrature_weights[1] = -0.3919800412014877;
422 quadrature_weights[2] = -0.0946154065661483;
426 quadrature_weights[0] = -0.3834640681451353;
427 quadrature_weights[1] = -0.3868753177747627;
428 quadrature_weights[2] = -0.1904351269501432;
429 quadrature_weights[3] = -0.03922548712995894;
433 quadrature_weights[0] = -0.2978934717828955;
434 quadrature_weights[1] = -0.3497762265132236;
435 quadrature_weights[2] = -0.234488290044052;
436 quadrature_weights[3] = -0.0989304595166356;
437 quadrature_weights[4] = -0.01891155214319462;
441 quadrature_weights[0] = -0.2387636625785478;
442 quadrature_weights[1] = -0.3082865732739458;
443 quadrature_weights[2] = -0.2453174265632108;
444 quadrature_weights[3] = -0.1420087565664786;
445 quadrature_weights[4] = -0.05545462232488041;
446 quadrature_weights[5] = -0.01016895869293513;
451 quadrature_weights[0] = -0.1961693894252476;
452 quadrature_weights[1] = -0.2703026442472726;
453 quadrature_weights[2] = -0.239681873007687;
454 quadrature_weights[3] = -0.1657757748104267;
455 quadrature_weights[4] = -0.0889432271377365;
456 quadrature_weights[5] = -0.03319430435645653;
457 quadrature_weights[6] = -0.005932787015162054;
461 quadrature_weights[0] = -0.164416604728002;
462 quadrature_weights[1] = -0.2375256100233057;
463 quadrature_weights[2] = -0.2268419844319134;
464 quadrature_weights[3] = -0.1757540790060772;
465 quadrature_weights[4] = -0.1129240302467932;
466 quadrature_weights[5] = -0.05787221071771947;
467 quadrature_weights[6] = -0.02097907374214317;
468 quadrature_weights[7] = -0.003686407104036044;
472 quadrature_weights[0] = -0.1400684387481339;
473 quadrature_weights[1] = -0.2097722052010308;
474 quadrature_weights[2] = -0.211427149896601;
475 quadrature_weights[3] = -0.1771562339380667;
476 quadrature_weights[4] = -0.1277992280331758;
477 quadrature_weights[5] = -0.07847890261203835;
478 quadrature_weights[6] = -0.0390225049841783;
479 quadrature_weights[7] = -0.01386729555074604;
480 quadrature_weights[8] = -0.002408041036090773;
484 quadrature_weights[0] = -0.12095513195457;
485 quadrature_weights[1] = -0.1863635425640733;
486 quadrature_weights[2] = -0.1956608732777627;
487 quadrature_weights[3] = -0.1735771421828997;
488 quadrature_weights[4] = -0.135695672995467;
489 quadrature_weights[5] = -0.0936467585378491;
490 quadrature_weights[6] = -0.05578772735275126;
491 quadrature_weights[7] = -0.02715981089692378;
492 quadrature_weights[8] = -0.00951518260454442;
493 quadrature_weights[9] = -0.001638157633217673;
498 quadrature_weights[0] = -0.1056522560990997;
499 quadrature_weights[1] = -0.1665716806006314;
500 quadrature_weights[2] = -0.1805632182877528;
501 quadrature_weights[3] = -0.1672787367737502;
502 quadrature_weights[4] = -0.1386970574017174;
503 quadrature_weights[5] = -0.1038334333650771;
504 quadrature_weights[6] = -0.06953669788988512;
505 quadrature_weights[7] = -0.04054160079499477;
506 quadrature_weights[8] = -0.01943540249522013;
507 quadrature_weights[9] = -0.006737429326043388;
508 quadrature_weights[10] = -0.001152486965101561;
512 quadrature_weights[0] = -0.09319269144393;
513 quadrature_weights[1] = -0.1497518275763289;
514 quadrature_weights[2] = -0.166557454364573;
515 quadrature_weights[3] = -0.1596335594369941;
516 quadrature_weights[4] = -0.1384248318647479;
517 quadrature_weights[5] = -0.1100165706360573;
518 quadrature_weights[6] = -0.07996182177673273;
519 quadrature_weights[7] = -0.0524069547809709;
520 quadrature_weights[8] = -0.03007108900074863;
521 quadrature_weights[9] = -0.01424924540252916;
522 quadrature_weights[10] = -0.004899924710875609;
523 quadrature_weights[11] = -0.000834029009809656;
531 return quadrature_weights;
539 const bool factor_out_singularity)
541 ((origin[0] == 0) || (origin[0] == 1)) ? (alpha == 1 ? n : 2 * n) : 4 * n)
542 , fraction(((origin[0] == 0) || (origin[0] == 1.)) ? 1. : origin[0])
567 unsigned int ns_offset = (
fraction == 1) ? n : 2 * n;
569 for (
unsigned int i = 0, j = ns_offset; i < n; ++i, ++j)
577 if ((alpha != 1) || (
fraction != 1))
597 if (factor_out_singularity ==
true)
598 for (
unsigned int i = 0; i <
size(); ++i)
603 "The singularity cannot be on a Gauss point of the same order!"));
606 Assert(denominator != 0.0,
608 "The quadrature formula you are using does not allow to "
609 "factor out the singularity, which is zero at one point."));
610 this->
weights[i] /= denominator;
619 const double eps = 1e-8;
620 bool on_edge =
false;
621 for (
unsigned int d = 0; d < 2; ++d)
622 on_edge = on_edge || (
std::abs(singularity[d]) < eps ||
623 std::abs(singularity[d] - 1.0) < eps);
624 const bool on_vertex =
638 const bool factor_out_singularity)
647 std::vector<QGaussOneOverR<2>> quads;
648 std::vector<Point<2>> origins;
651 quads.emplace_back(n, 3, factor_out_singularity);
652 quads.emplace_back(n, 2, factor_out_singularity);
653 quads.emplace_back(n, 1, factor_out_singularity);
654 quads.emplace_back(n, 0, factor_out_singularity);
656 origins.emplace_back(0., 0.);
657 origins.emplace_back(singularity[0], 0.);
658 origins.emplace_back(0., singularity[1]);
659 origins.push_back(singularity);
664 unsigned int q_id = 0;
667 for (
unsigned int box = 0; box < 4; ++box)
671 double area = dist[0] * dist[1];
673 for (
unsigned int q = 0; q < quads[box].size(); ++q, ++q_id)
675 const Point<2> &qp = quads[box].point(q);
677 origins[box] +
Point<2>(dist[0] * qp[0], dist[1] * qp[1]);
678 this->
weights[q_id] = quads[box].weight(q) * area;
686 const unsigned int vertex_index,
687 const bool factor_out_singularity)
723 std::vector<double> & ws = this->
weights;
726 for (
unsigned int q = 0; q < gauss.
size(); ++q)
730 ps[q][1] = gp[0] *
std::tan(pi4 * gp[1]);
732 if (factor_out_singularity)
736 ws[gauss.
size() + q] = ws[q];
737 ps[gauss.
size() + q][0] = ps[q][1];
738 ps[gauss.
size() + q][1] = ps[q][0];
743 switch (vertex_index)
763 if (vertex_index != 0)
764 for (
unsigned int q = 0; q <
size(); ++q)
766 double x = ps[q][0] - .5, y = ps[q][1] - .5;
768 ps[q][0] = R00 * x + R01 * y + .5;
769 ps[q][1] = R10 * x + R11 * y + .5;
778 std::vector<unsigned int> permutation(quad.
size());
779 for (
unsigned int i = 0; i < quad.
size(); ++i)
782 std::sort(permutation.begin(),
784 [
this](
const unsigned int x,
const unsigned int y) {
785 return this->compare_weights(x, y);
795 for (
unsigned int i = 0; i < quad.
size(); ++i)
799 if (permutation[i] != i)
809 return (this->weights[a] < this->weights[b]);
895 const double eta_bar = singularity[0] * 2. - 1.;
896 const double eta_star = eta_bar * eta_bar - 1.;
900 std::vector<double> weights_dummy(
weights.size());
901 unsigned int cont = 0;
902 const double tol = 1e-10;
908 weights_dummy[d - cont] =
weights[d];
920 weights.resize(weights_dummy.size() - 1);
937 gamma_bar = (eta_bar * eta_star +
std::abs(eta_star)) /
941 (eta_bar * eta_star -
std::abs(eta_star)) /
950 double eta = (Utilities::fixed_power<3>(gamma - gamma_bar) +
951 gamma_bar * (gamma_bar * gamma_bar + 3)) /
952 (1 + 3 * gamma_bar * gamma_bar);
954 double J = 3 * ((gamma - gamma_bar) * (gamma - gamma_bar)) /
955 (1 + 3 * gamma_bar * gamma_bar);
972 std::vector<double> points(n);
974 for (
unsigned short i = 0; i < n; ++i)
981 (1. +
double(2 * i + 1) /
double(2 * (n - 1) + 2))));
994 std::vector<double> weights(n);
996 for (
unsigned short i = 0; i < n; ++i)
1012 Assert(n > 0,
ExcMessage(
"Need at least one point for the quadrature rule"));
1016 for (
unsigned int i = 0; i < this->
size(); ++i)
1039 std::vector<double> points(n);
1041 for (
unsigned short i = 0; i < n; ++i)
1047 case ::QGaussRadauChebyshev<1>::left:
1053 (1 + 2 *
double(i) / (2 *
double(n - 1) + 1.))));
1057 case ::QGaussRadauChebyshev<1>::right:
1062 (2 *
double(n - 1) + 1.))));
1070 "This constructor can only be called with either "
1071 "QGaussRadauChebyshev::left or QGaussRadauChebyshev::right as "
1072 "second argument."));
1085 std::vector<double> weights(n);
1087 for (
unsigned short i = 0; i < n; ++i)
1090 weights[i] = 2. *
numbers::PI / double(2 * (n - 1) + 1.);
1110 std::vector<double> p =
1112 std::vector<double> w =
1115 for (
unsigned int i = 0; i < this->
size(); ++i)
1142 std::vector<double> points(n);
1144 for (
unsigned short i = 0; i < n; ++i)
1159 std::vector<double> weights(n);
1161 for (
unsigned short i = 0; i < n; ++i)
1165 if (i == 0 || i == (n - 1))
1182 "Need at least two points for Gauss-Lobatto quadrature rule"));
1183 std::vector<double> p =
1185 std::vector<double> w =
1188 for (
unsigned int i = 0; i < this->
size(); ++i)
1206 std::vector<Point<dim>> qpoints;
1207 std::vector<double> weights;
1209 for (
unsigned int i = 0; i < quad.
size(); ++i)
1216 for (
int d = 0; d < dim; ++d)
1217 r += quad.
point(i)[d];
1220 this->quadrature_points.push_back(quad.
point(i));
1221 this->weights.push_back(quad.
weight(i));
1229template <
int spacedim>
1235 ExcMessage(
"Invalid combination of dim and spacedim ."));
1237 for (
unsigned int d = 0; d < dim; ++d)
1241 const double J =
std::abs(B.determinant());
1247 std::vector<Point<spacedim>> qp(this->size());
1248 std::vector<double> w(this->size());
1250 for (
unsigned int i = 0; i < this->size(); ++i)
1254 w[i] = J * this->weight(i);
1263template <
int spacedim>
1266 const std::vector<std::array<
Point<spacedim>, dim + 1>> &simplices)
const
1268 Assert(!(dim == 1 && spacedim == 1),
1269 ExcMessage(
"This function is not supposed to work in 1D-1d case."));
1271 ExcMessage(
"Invalid combination of dim and spacedim ."));
1273 std::vector<Point<spacedim>> qp;
1274 std::vector<double> ws;
1275 for (
const auto &simplex : simplices)
1277 const auto rule = this->compute_affine_transformation(simplex);
1278 std::transform(rule.get_points().begin(),
1279 rule.get_points().end(),
1280 std::back_inserter(qp),
1282 std::transform(rule.get_weights().begin(),
1283 rule.get_weights().end(),
1284 std::back_inserter(ws),
1285 [&](
const double w) { return w; });
1299 for (
unsigned int i = 0; i < base.
size(); ++i)
1301 const auto &q = base.
point(i);
1302 const auto w = base.
weight(i);
1304 const auto xhat = q[0];
1305 const auto yhat = q[1];
1311 const double r = xhat / (st + ct);
1313 const double J = pi * xhat / (2 * (
std::sin(pi * yhat) + 1));
1336 for (
unsigned int i = 0; i < base.
size(); ++i)
1338 const auto &q = base.
point(i);
1339 const auto w = base.
weight(i);
1341 const auto xhat = q[0];
1342 const auto yhat = q[1];
1344 const double x =
std::pow(xhat, beta) * (1 - yhat);
1345 const double y =
std::pow(xhat, beta) * yhat;
1347 const double J = beta *
std::pow(xhat, 2. * beta - 1.);
1367 "The split point should be inside the unit reference cell."));
1369 std::array<Point<dim>, dim + 1>
vertices;
1377 for (
unsigned int start = 0; start < (dim > 2 ? 2 : 1); ++start)
1379 for (
unsigned int i = 0; i < dim; ++i)
1385 this->quadrature_points.insert(this->quadrature_points.end(),
1386 quad.get_points().begin(),
1387 quad.get_points().end());
1388 this->weights.insert(this->weights.end(),
1389 quad.get_weights().begin(),
1390 quad.get_weights().end());
1402 if (dim == 0 || dim == 1)
1404 const ::QGauss<dim> quad(n_points_1D);
1407 this->
weights = quad.get_weights();
1411 if (n_points_1D == 1)
1413 const double p = 1.0 / 3.0;
1415 this->
weights.emplace_back(0.5);
1417 else if (n_points_1D == 2)
1423 const double Q12 = 1.0 / 2.0;
1425 0.1550510257216822);
1427 0.6449489742783178);
1429 0.1550510257216822);
1431 0.6449489742783178);
1433 this->
weights.emplace_back(0.31804138174397717 * Q12);
1434 this->
weights.emplace_back(0.18195861825602283 * Q12);
1435 this->
weights.emplace_back(0.31804138174397717 * Q12);
1436 this->
weights.emplace_back(0.18195861825602283 * Q12);
1438 else if (n_points_1D == 3)
1441 const double p0 = 2.0 / 7.0 -
std::sqrt(15.0) / 21.0;
1442 const double p1 = 2.0 / 7.0 +
std::sqrt(15.0) / 21.0;
1443 const double p2 = 3.0 / 7.0 - 2.0 *
std::sqrt(15.0) / 21.0;
1444 const double p3 = 3.0 / 7.0 + 2.0 *
std::sqrt(15.0) / 21.0;
1453 const double q12 = 0.5;
1454 const double w0 = 9.0 / 40.0;
1455 const double w1 = 31.0 / 240.0 -
std::sqrt(15.0) / 1200.0;
1456 const double w2 = 31.0 / 240.0 +
std::sqrt(15.0) / 1200.0;
1457 this->
weights.emplace_back(q12 * w0);
1458 this->
weights.emplace_back(q12 * w1);
1459 this->
weights.emplace_back(q12 * w1);
1460 this->
weights.emplace_back(q12 * w1);
1461 this->
weights.emplace_back(q12 * w2);
1462 this->
weights.emplace_back(q12 * w2);
1463 this->
weights.emplace_back(q12 * w2);
1465 else if (n_points_1D == 4)
1473 if (n_points_1D == 1)
1475 const double Q14 = 1.0 / 4.0;
1476 const double Q16 = 1.0 / 6.0;
1479 this->
weights.emplace_back(Q16);
1486 else if (n_points_1D == 2)
1488 const double Q16 = 1.0 / 6.0;
1489 this->
weights.emplace_back(0.1223220027573451 * Q16);
1490 this->
weights.emplace_back(0.1280664127107469 * Q16);
1491 this->
weights.emplace_back(0.1325680271444452 * Q16);
1492 this->
weights.emplace_back(0.1406244096604032 * Q16);
1493 this->
weights.emplace_back(0.2244151669175574 * Q16);
1494 this->
weights.emplace_back(0.2520039808095023 * Q16);
1498 0.01271836631368145);
1501 0.3621268299455338);
1503 0.01140332944455717,
1504 0.3586207204668839);
1507 0.6384932999617267);
1510 0.1308471689520965);
1513 0.1403728057942107);
1517 else if (n_points_1D == 3)
1522 else if (n_points_1D == 4)
1532 "QGaussSimplex is currently only implemented for "
1533 "n_points_1D = 1, 2, 3, and 4 while you are asking for "
1540 template <std::
size_t b_dim>
1541 std::vector<std::array<double, b_dim>>
1542 all_permutations(
const std::array<double, b_dim> &b_point)
1544 std::vector<std::array<double, b_dim>> output;
1549 std::array<double, b_dim> temp = b_point;
1550 std::sort(temp.begin(), temp.end());
1553 output.push_back(temp);
1555 while (std::next_permutation(temp.begin(), temp.end()));
1565 const unsigned int n_points_1D,
1566 const bool use_odd_order)
1578 std::array<double, dim + 1> centroid;
1579 std::fill(centroid.begin(), centroid.end(), 1.0 / (dim + 1.0));
1580 std::vector<std::vector<std::array<double, dim + 1>>> b_point_permutations;
1581 std::vector<double> b_weights;
1591 auto process_point_1 = [&](
const double a,
const double w) {
1592 const double b = 1.0 - dim * a;
1593 std::array<double, dim + 1> b_point;
1594 std::fill(b_point.begin(), b_point.begin() + dim, a);
1597 b_weights.push_back(w);
1598 b_point_permutations.push_back(all_permutations(b_point));
1603 auto process_point_2 = [&](
const double a,
const double w) {
1605 const double b = (1.0 - 2.0 * a) / 2.0;
1606 std::array<double, dim + 1> b_point;
1607 std::fill(b_point.begin(), b_point.begin() + dim - 1, a);
1608 b_point[dim - 1] = b;
1611 b_weights.push_back(w);
1612 b_point_permutations.push_back(all_permutations(b_point));
1618 auto process_point_3 = [&](
const double a,
const double b,
const double w) {
1619 const double c = 1.0 - (dim - 1.0) * a - b;
1620 std::array<double, dim + 1> b_point;
1621 std::fill(b_point.begin(), b_point.begin() + dim - 1, a);
1622 b_point[dim - 1] = b;
1625 b_weights.push_back(w);
1626 b_point_permutations.push_back(all_permutations(b_point));
1629 switch (n_points_1D)
1638 b_point_permutations.push_back({centroid});
1639 b_weights.push_back(1.0000000000000000e+00);
1644 process_point_1(1.6666666666666669e-01,
1645 3.3333333333333331e-01);
1652 b_point_permutations.push_back({centroid});
1653 b_weights.push_back(1.0000000000000000e+00);
1658 process_point_1(1.3819660112501050e-01,
1659 2.5000000000000000e-01);
1671 process_point_1(9.1576213509770743e-02, 1.0995174365532187e-01);
1672 process_point_1(4.4594849091596489e-01, 2.2338158967801147e-01);
1678 process_point_1(3.2816330251638171e-01,
1679 1.3621784253708741e-01);
1680 process_point_1(1.0804724989842859e-01,
1681 1.1378215746291261e-01);
1700 b_point_permutations.push_back({centroid});
1701 b_weights.push_back(2.2500000000000001e-01);
1702 process_point_1(1.0128650732345634e-01,
1703 1.2593918054482714e-01);
1704 process_point_1(4.7014206410511511e-01,
1705 1.3239415278850619e-01);
1710 process_point_1(6.3089014491502227e-02,
1711 5.0844906370206819e-02);
1712 process_point_1(2.4928674517091043e-01,
1713 1.1678627572637937e-01);
1714 process_point_3(5.3145049844816938e-02,
1715 3.1035245103378439e-01,
1716 8.2851075618373571e-02);
1723 process_point_1(3.1088591926330061e-01,
1724 1.1268792571801590e-01);
1725 process_point_1(9.2735250310891248e-02,
1726 7.3493043116361956e-02);
1727 process_point_2(4.5503704125649642e-02,
1728 4.2546020777081472e-02);
1733 process_point_1(4.0673958534611372e-02,
1734 1.0077211055320640e-02);
1735 process_point_1(3.2233789014227548e-01,
1736 5.5357181543654717e-02);
1737 process_point_1(2.1460287125915201e-01,
1738 3.9922750258167487e-02);
1739 process_point_3(6.3661001875017442e-02,
1740 6.0300566479164919e-01,
1741 4.8214285714285710e-02);
1755 process_point_1(3.3730648554587850e-02,
1756 1.6545050110792131e-02);
1757 process_point_1(4.7430969250471822e-01,
1758 7.7086646185986069e-02);
1759 process_point_1(2.4157738259540357e-01,
1760 1.2794417123015558e-01);
1761 process_point_3(4.7036644652595216e-02,
1762 1.9868331479735168e-01,
1763 5.5878732903199779e-02);
1768 b_point_permutations.push_back({centroid});
1769 b_weights.push_back(1.4431560767778717e-01);
1770 process_point_1(5.0547228317030957e-02,
1771 3.2458497623198079e-02);
1772 process_point_1(4.5929258829272313e-01,
1773 9.5091634267284619e-02);
1774 process_point_1(1.7056930775176021e-01,
1775 1.0321737053471824e-01);
1776 process_point_3(8.3947774099575878e-03,
1777 2.6311282963463811e-01,
1778 2.7230314174434993e-02);
1785 b_point_permutations.push_back({centroid});
1786 b_weights.push_back(9.5485289464130846e-02);
1787 process_point_1(3.1570114977820279e-01,
1788 4.2329581209967028e-02);
1789 process_point_2(5.0489822598396350e-02,
1790 3.1896927832857580e-02);
1791 process_point_3(1.8883383102600099e-01,
1792 5.7517163758699996e-01,
1793 3.7207130728334620e-02);
1794 process_point_3(2.1265472541483140e-02,
1795 8.1083024109854862e-01,
1796 8.1107708299033420e-03);
1801 process_point_1(1.0795272496221089e-01,
1802 2.6426650908408830e-02);
1803 process_point_1(1.8510948778258660e-01,
1804 5.2031747563738531e-02);
1805 process_point_1(4.2316543684767283e-02,
1806 7.5252561535401989e-03);
1807 process_point_1(3.1418170912403898e-01,
1808 4.1763782856934897e-02);
1809 process_point_2(4.3559132858383021e-01,
1810 3.6280930261308818e-02);
1811 process_point_3(2.1433930127130570e-02,
1812 7.1746406342630831e-01,
1813 7.1569028908444327e-03);
1814 process_point_3(2.0413933387602909e-01,
1815 5.8379737830214440e-01,
1816 1.5453486150960340e-02);
1830 b_point_permutations.push_back({centroid});
1831 b_weights.push_back(9.7135796282798836e-02);
1832 process_point_1(4.4729513394452691e-02,
1833 2.5577675658698031e-02);
1834 process_point_1(4.8968251919873762e-01,
1835 3.1334700227139071e-02);
1836 process_point_1(4.3708959149293664e-01,
1837 7.7827541004774278e-02);
1838 process_point_1(1.8820353561903275e-01,
1839 7.9647738927210249e-02);
1840 process_point_3(3.6838412054736258e-02,
1841 2.2196298916076568e-01,
1842 4.3283539377289376e-02);
1847 b_point_permutations.push_back({centroid});
1848 b_weights.push_back(8.1743329146285973e-02);
1849 process_point_1(3.2055373216943517e-02,
1850 1.3352968813149567e-02);
1851 process_point_1(1.4216110105656438e-01,
1852 4.5957963604744731e-02);
1853 process_point_3(2.8367665339938453e-02,
1854 1.6370173373718250e-01,
1855 2.5297757707288385e-02);
1856 process_point_3(2.9619889488729734e-02,
1857 3.6914678182781102e-01,
1858 3.4184648162959429e-02);
1859 process_point_3(1.4813288578382056e-01,
1860 3.2181299528883545e-01,
1861 6.3904906396424044e-02);
1868 b_point_permutations.push_back({centroid});
1869 b_weights.push_back(5.8010548912480253e-02);
1870 process_point_1(6.1981697552226933e-10,
1871 6.4319281759256394e-05);
1872 process_point_1(1.6077453539526160e-01,
1873 2.3173338462425461e-02);
1874 process_point_1(3.2227652182142102e-01,
1875 2.9562912335429289e-02);
1876 process_point_1(4.5108918345413578e-02,
1877 8.0639799796161822e-03);
1878 process_point_2(1.1229654600437609e-01,
1879 3.8134080103702457e-02);
1880 process_point_3(4.5887144875245922e-01,
1881 2.5545792330413102e-03,
1882 8.3844221982985519e-03);
1883 process_point_3(3.3775870685338598e-02,
1884 7.1835032644207453e-01,
1885 1.0234559352745330e-02);
1886 process_point_3(1.8364136980992790e-01,
1887 3.4415910578175279e-02,
1888 2.0524915967988139e-02);
1893 b_point_permutations.push_back({centroid});
1894 b_weights.push_back(4.7399773556020743e-02);
1895 process_point_1(3.1225006869518868e-01,
1896 2.6937059992268701e-02);
1897 process_point_1(1.1430965385734609e-01,
1898 9.8691597167933822e-03);
1899 process_point_3(4.1043073921896539e-01,
1900 1.6548602561961109e-01,
1901 1.1393881220195230e-02);
1902 process_point_3(6.1380088247906528e-03,
1903 9.4298876734520487e-01,
1904 3.6194434433925362e-04);
1905 process_point_3(1.2105018114558939e-01,
1906 4.7719037990428043e-01,
1907 2.5739731980456069e-02);
1908 process_point_3(3.2779468216442620e-02,
1909 5.9425626948000698e-01,
1910 1.0135871679755789e-02);
1911 process_point_3(3.2485281564823047e-02,
1912 8.0117728465834437e-01,
1913 6.5761472770359038e-03);
1914 process_point_3(1.7497934218393901e-01,
1915 6.2807184547536599e-01,
1916 1.2907035798861989e-02);
1929 b_point_permutations.push_back({centroid});
1930 b_weights.push_back(8.5761179732224219e-02);
1931 process_point_1(2.8485417614371900e-02, 1.0431870512894697e-02);
1932 process_point_1(4.9589190096589092e-01, 1.6606273054585369e-02);
1933 process_point_1(1.0263548271224643e-01, 3.8630759237019321e-02);
1934 process_point_1(4.3846592676435220e-01, 6.7316154079468296e-02);
1935 process_point_1(2.1021995670317828e-01, 7.0515684111716576e-02);
1936 process_point_3(7.3254276860644785e-03,
1937 1.4932478865208237e-01,
1938 1.0290289572953278e-02);
1939 process_point_3(4.6010500165429957e-02,
1940 2.8958112563770588e-01,
1941 4.0332476640500554e-02);
1946 process_point_1(2.4646363436335583e-02, 7.9316425099736389e-03);
1947 process_point_1(4.8820375094554153e-01, 2.4266838081452032e-02);
1948 process_point_1(1.0925782765935427e-01, 2.8486052068877544e-02);
1949 process_point_1(4.4011164865859309e-01, 4.9918334928060942e-02);
1950 process_point_1(2.7146250701492608e-01, 6.2541213195902765e-02);
1951 process_point_3(2.1382490256170616e-02,
1952 1.2727971723358933e-01,
1953 1.5083677576511438e-02);
1954 process_point_3(2.3034156355267121e-02,
1955 2.9165567973834094e-01,
1956 2.1783585038607559e-02);
1957 process_point_3(1.1629601967792658e-01,
1958 2.5545422863851736e-01,
1959 4.3227363659414209e-02);
1968 b_point_permutations.push_back({centroid});
1969 b_weights.push_back(6.7960036586831640e-02);
1970 process_point_1(2.1509681108843159e-02, 6.0523371035391717e-03);
1971 process_point_1(4.8907694645253935e-01, 2.3994401928894731e-02);
1972 process_point_1(4.2694141425980042e-01, 5.5601967530453329e-02);
1973 process_point_1(2.2137228629183292e-01, 5.8278485119199981e-02);
1974 process_point_3(5.1263891023823893e-03,
1975 2.7251581777342970e-01,
1976 9.5906810035432631e-03);
1977 process_point_3(2.4370186901093827e-02,
1978 1.1092204280346341e-01,
1979 1.4965401105165668e-02);
1980 process_point_3(8.7895483032197297e-02,
1981 1.6359740106785048e-01,
1982 2.4179039811593819e-02);
1983 process_point_3(6.8012243554206653e-02,
1984 3.0844176089211778e-01,
1985 3.4641276140848373e-02);
1990 process_point_1(1.9390961248701044e-02, 4.9234036024000819e-03);
1991 process_point_1(6.1799883090872587e-02, 1.4433699669776668e-02);
1992 process_point_1(4.8896391036217862e-01, 2.1883581369428889e-02);
1993 process_point_1(4.1764471934045394e-01, 3.2788353544125348e-02);
1994 process_point_1(1.7720553241254344e-01, 4.2162588736993016e-02);
1995 process_point_1(2.7347752830883865e-01, 5.1774104507291585e-02);
1996 process_point_3(1.2683309328720416e-03,
1997 1.1897449769695684e-01,
1998 5.0102288385006719e-03);
1999 process_point_3(1.4646950055654417e-02,
2000 2.9837288213625779e-01,
2001 1.4436308113533840e-02);
2002 process_point_3(5.7124757403647919e-02,
2003 1.7226668782135557e-01,
2004 2.4665753212563674e-02);
2005 process_point_3(9.2916249356971847e-02,
2006 3.3686145979634496e-01,
2007 3.8571510787060684e-02);
2015 for (
unsigned int permutation_n = 0; permutation_n < b_weights.size();
2018 for (
const std::array<double, dim + 1> &b_point :
2019 b_point_permutations[permutation_n])
2021 const double volume = (dim == 2 ? 1.0 / 2.0 : 1.0 / 6.0);
2022 this->
weights.emplace_back(volume * b_weights[permutation_n]);
2024 for (
int d = 0; d < dim; ++d)
2025 c_point[d] = b_point[d];
2047 const unsigned int n_copies)
2057 const unsigned int n_copies)
2063 setup_qiterated_1D(base_quad, n_copies);
2068 const auto n_refinements =
2069 static_cast<unsigned int>(std::round(std::log2(n_copies)));
2070 Assert((1u << n_refinements) == n_copies,
2071 ExcMessage(
"The number of copies must be a power of 2."));
2073 const auto reference_cell = ReferenceCells::get_simplex<dim>();
2077 reference_cell.template get_default_linear_mapping<dim>();
2084 std::vector<Point<dim>> points;
2085 std::vector<double> weights;
2089 for (
unsigned int qp = 0; qp < base_quad.
size(); ++qp)
2092 weights.push_back(fe_values.
JxW(qp));
2117 for (
unsigned int i = 0; i < quad_line.
size(); ++i)
2118 for (
unsigned int j = 0; j < quad_tri.
size(); ++j)
2121 quad_tri.
point(j)[1],
2122 quad_line.
point(i)[0]);
2139 if (n_points_1D == 1)
2141 const double Q14 = 1.0 / 4.0;
2142 const double Q43 = 4.0 / 3.0;
2145 this->
weights.emplace_back(Q43);
2147 else if (n_points_1D == 2)
2150 this->
quadrature_points.emplace_back(-0.26318405556971, -0.26318405556971, 0.54415184401122);
2151 this->
quadrature_points.emplace_back(-0.50661630334979, -0.50661630334979, 0.12251482265544);
2152 this->
quadrature_points.emplace_back(-0.26318405556971, +0.26318405556971, 0.54415184401122);
2153 this->
quadrature_points.emplace_back(-0.50661630334979, +0.50661630334979, 0.12251482265544);
2154 this->
quadrature_points.emplace_back(+0.26318405556971, -0.26318405556971, 0.54415184401122);
2155 this->
quadrature_points.emplace_back(+0.50661630334979, -0.50661630334979, 0.12251482265544);
2156 this->
quadrature_points.emplace_back(+0.26318405556971, +0.26318405556971, 0.54415184401122);
2157 this->
quadrature_points.emplace_back(+0.50661630334979, +0.50661630334979, 0.12251482265544);
2160 this->
weights.emplace_back(0.10078588207983);
2161 this->
weights.emplace_back(0.23254745125351);
2162 this->
weights.emplace_back(0.10078588207983);
2163 this->
weights.emplace_back(0.23254745125351);
2164 this->
weights.emplace_back(0.10078588207983);
2165 this->
weights.emplace_back(0.23254745125351);
2166 this->
weights.emplace_back(0.10078588207983);
2167 this->
weights.emplace_back(0.23254745125351);
2271 const std::vector<std::array<
Point<2>, 1 + 1>> &simplices)
const;
2275 const std::vector<std::array<
Point<3>, 1 + 1>> &simplices)
const;
2279 const std::vector<std::array<
Point<2>, 2 + 1>> &simplices)
const;
2283 const std::vector<std::array<
Point<3>, 2 + 1>> &simplices)
const;
2287 const std::vector<std::array<
Point<3>, 3 + 1>> &simplices)
const;
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
double JxW(const unsigned int q_point) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
Abstract base class for mapping classes.
QDuffy(const Quadrature< 1 > &radial_quadrature, const Quadrature< 1 > &angular_quadrature, const double beta=1.0)
QGaussChebyshev(const unsigned int n)
Generate a formula with n quadrature points.
QGaussLobattoChebyshev(const unsigned int n)
Generate a formula with n quadrature points.
QGaussLobatto(const unsigned int n)
QGaussLogR(const unsigned int n, const Point< dim > &x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false)
QGaussLog(const unsigned int n, const bool revert=false)
static std::vector< double > get_quadrature_points(const unsigned int n)
static std::vector< double > get_quadrature_weights(const unsigned int n)
QGaussOneOverR(const unsigned int n, const Point< dim > &singularity, const bool factor_out_singular_weight=false)
static unsigned int quad_size(const Point< dim > &singularity, const unsigned int n)
QGaussPyramid(const unsigned int n_points_1D)
QGaussRadauChebyshev(const unsigned int n, EndPoint ep=QGaussRadauChebyshev::left)
Generate a formula with n quadrature points.
QGaussSimplex(const unsigned int n_points_1D)
QGaussWedge(const unsigned int n_points_1D)
QGauss(const unsigned int n)
QIteratedSimplex(const Quadrature< dim > &base_quadrature, const unsigned int n_copies)
Quadrature< spacedim > compute_affine_transformation(const std::array< Point< spacedim >, dim+1 > &vertices) const
QSimplex(const Quadrature< dim > &quad)
Quadrature< spacedim > mapped_quadrature(const std::vector< std::array< Point< spacedim >, dim+1 > > &simplices) const
QSorted(const Quadrature< dim > &quad)
bool compare_weights(const unsigned int a, const unsigned int b) const
QSplit(const QSimplex< dim > &base, const Point< dim > &split_point)
QTelles(const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
QTrianglePolar(const Quadrature< 1 > &radial_quadrature, const Quadrature< 1 > &angular_quadrature)
QWitherdenVincentSimplex(const unsigned int n_points_1D, const bool use_odd_order=true)
std::vector< Point< dim > > quadrature_points
Quadrature & operator=(const Quadrature< dim > &)
const Point< dim > & point(const unsigned int i) const
bool is_tensor_product_flag
double weight(const unsigned int i) const
std::vector< double > weights
unsigned int size() const
void refine_global(const unsigned int times=1)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ update_JxW_values
Transformed quadrature weights.
@ update_quadrature_points
Transformed quadrature points.
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Number jacobi_polynomial_value(const unsigned int degree, const int alpha, const int beta, const Number x)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
std::vector< double > get_quadrature_weights(const unsigned int n)
std::vector< double > get_quadrature_points(const unsigned int n)
std::vector< double > get_quadrature_points(const unsigned int n)
std::vector< double > get_quadrature_weights(const unsigned int n)
std::vector< long double > compute_quadrature_weights(const std::vector< long double > &x, const int alpha, const int beta)
long double gamma(const unsigned int n)
std::vector< double > get_quadrature_points(const unsigned int n, ::QGaussRadauChebyshev< 1 >::EndPoint ep)
std::vector< double > get_quadrature_weights(const unsigned int n, ::QGaussRadauChebyshev< 1 >::EndPoint ep)
static constexpr double PI_2
static constexpr double PI
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > tan(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
static Point< dim > unit_cell_vertex(const unsigned int vertex)
const ::Triangulation< dim, spacedim > & tria