39 for (j = 0; j < ny; j++)
42 for (i = 0; i < nx; i++)
63 for (
int iz = 0; iz < nz; ++iz)
66 for (
int iy = 0; iy < ny; ++iy)
69 for (
int ix = 0; ix < nx; ++ix)
85 if (weights.Size() != GetNPoints())
88 for (
int i = 0; i < GetNPoints(); i++)
90 weights[i] = IntPoint(i).weight;
96 void IntegrationRule::SetPointIndices()
98 for (
int i = 0; i < Size(); i++)
100 IntPoint(i).index = i;
104 void IntegrationRule::GrundmannMollerSimplexRule(
int s,
int n)
108 const int d = 2*s + 1;
113 for (
int i = 1; i < fact.Size(); i++)
115 fact(i) = fact(i - 1)*i;
120 for (
int i = 0; i <= n; i++)
122 np *= (s + i + 1),
f *= (i + 1);
129 for (
int i = 0; i <=
s; i++)
133 weight = pow(2., -2*s)*pow(static_cast<double>(d + n - 2*i),
134 d)/fact(i)/fact(d + n - i);
146 IntegrationPoint &ip = IntPoint(pt++);
148 ip.x = double(2*beta[0] + 1)/(d + n - 2*i);
149 ip.y = double(2*beta[1] + 1)/(d + n - 2*i);
152 ip.z = double(2*beta[2] + 1)/(d + n - 2*i);
166 for (j--; j >= 0; j--)
180 class HP_Quadrature1D
183 mpfr_t pi, z, pp, p1, p2, p3, dz, w, rtol;
186 static const mpfr_rnd_t rnd = GMP_RNDN;
187 static const int default_prec = 128;
190 HP_Quadrature1D(
const int prec = default_prec)
192 mpfr_inits2(prec, pi, z, pp, p1, p2, p3, dz, w, rtol, (mpfr_ptr) 0);
193 mpfr_const_pi(pi, rnd);
194 mpfr_set_si_2exp(rtol, 1, -32, rnd);
200 void SetRelTol(
const int exponent = -32)
202 mpfr_set_si_2exp(rtol, 1, exponent, rnd);
208 void ComputeGaussLegendrePoint(
const int n,
const int k)
210 MFEM_ASSERT(n > 0 && 0 <= k && k < n,
"invalid n = " << n
211 <<
" and/or k = " << k);
213 int i = (k < (n+1)/2) ? k+1 : n-k;
218 mpfr_set_si(z, n+1-2*i, rnd);
219 mpfr_div_si(z, z, 2*n+1, rnd);
220 mpfr_mul(z, z, pi, rnd);
226 mpfr_set_si(p2, 1, rnd);
227 mpfr_set(p1, z, rnd);
228 for (
int j = 2; j <= n; j++)
230 mpfr_set(p3, p2, rnd);
231 mpfr_set(p2, p1, rnd);
233 mpfr_mul_si(p1, z, 2*j-1, rnd);
234 mpfr_mul_si(p3, p3, j-1, rnd);
235 mpfr_fms(p1, p1, p2, p3, rnd);
236 mpfr_div_si(p1, p1, j, rnd);
242 mpfr_fms(pp, z, p1, p2, rnd);
243 mpfr_mul_si(pp, pp, n, rnd);
244 mpfr_sqr(p2, z, rnd);
245 mpfr_sub_si(p2, p2, 1, rnd);
246 mpfr_div(pp, pp, p2, rnd);
251 mpfr_div(dz, p1, pp, rnd);
254 mpfr_si_sub(atol, 1, z, rnd);
255 mpfr_mul(atol, atol, rtol, rnd);
256 if (mpfr_cmpabs(dz, atol) <= 0)
262 mpfr_sub(z, z, dz, rnd);
266 mpfr_si_sub(z, 1, z, rnd);
267 mpfr_div_2si(z, z, 1, rnd);
270 mpfr_sqr(w, pp, rnd);
271 mpfr_mul_2si(w, w, 2, rnd);
272 mpfr_mul(w, w, z, rnd);
273 mpfr_si_sub(p1, 1, z, rnd);
274 mpfr_mul(w, w, p1, rnd);
275 mpfr_si_div(w, 1, w, rnd);
277 if (k >= (n+1)/2) { mpfr_swap(z, p1); }
283 void ComputeGaussLobattoPoint(
const int n,
const int k)
285 MFEM_ASSERT(n > 1 && 0 <= k && k < n,
"invalid n = " << n
286 <<
" and/or k = " << k);
288 int i = (k < (n+1)/2) ? k : n-1-k;
292 mpfr_set_si(z, 0, rnd);
293 mpfr_set_si(p1, 1, rnd);
294 mpfr_set_si(w, n*(n-1), rnd);
295 mpfr_si_div(w, 1, w, rnd);
300 mpfr_set_si(z, 2*i-n+1, rnd);
301 mpfr_div_si(z, z, 2*(n-1), rnd);
302 mpfr_mul(z, pi, z, rnd);
305 for (
int iter = 0 ; true ; ++iter)
308 mpfr_set_si(p1, 1, rnd);
309 mpfr_set(p2, z, rnd);
311 for (
int l = 1 ; l < (n-1) ; ++l)
314 mpfr_mul_si(p1, p1, l, rnd);
315 mpfr_mul_si(p3, z, 2*l+1, rnd);
316 mpfr_fms(p3, p3, p2, p1, rnd);
317 mpfr_div_si(p3, p3, l+1, rnd);
319 mpfr_set(p1, p2, rnd);
320 mpfr_set(p2, p3, rnd);
324 mpfr_fms(dz, z, p2, p1, rnd);
325 mpfr_mul_si(p3, p2, n, rnd);
326 mpfr_div(dz, dz, p3, rnd);
328 mpfr_sub(z, z, dz, rnd);
331 mpfr_add_si(atol, z, 1, rnd);
332 mpfr_mul(atol, atol, rtol, rnd);
334 if (mpfr_cmpabs(dz, atol) <= 0)
340 MFEM_VERIFY(iter < 8,
"n = " << n <<
", i = " << i
341 <<
", dz = " << mpfr_get_d(dz, rnd));
344 mpfr_add_si(z, z, 1, rnd);
345 mpfr_div_2si(z, z, 1, rnd);
347 mpfr_si_sub(p1, 1, z, rnd);
349 mpfr_sqr(w, p2, rnd);
350 mpfr_mul_si(w, w, n*(n-1), rnd);
351 mpfr_si_div(w, 1, w, rnd);
353 if (k >= (n+1)/2) { mpfr_swap(z, p1); }
356 double GetPoint()
const {
return mpfr_get_d(z, rnd); }
357 double GetSymmPoint()
const {
return mpfr_get_d(p1, rnd); }
358 double GetWeight()
const {
return mpfr_get_d(w, rnd); }
360 const mpfr_t &GetHPPoint()
const {
return z; }
361 const mpfr_t &GetHPSymmPoint()
const {
return p1; }
362 const mpfr_t &GetHPWeight()
const {
return w; }
366 mpfr_clears(pi, z, pp, p1, p2, p3, dz, w, rtol, (mpfr_ptr) 0);
371 #endif // MFEM_USE_MPFR
396 const int m = (n+1)/2;
398 #ifndef MFEM_USE_MPFR
400 for (
int i = 1; i <= m; i++)
402 double z = cos(M_PI * (i - 0.25) / (n + 0.5));
403 double pp, p1, dz, xi = 0.;
409 for (
int j = 2; j <= n; j++)
413 p1 = ((2 * j - 1) * z * p2 - (j - 1) * p3) / j;
417 pp = n * (z*p1-p2) / (z*z - 1);
421 if (fabs(dz) < 1e-16)
425 xi = ((1 - z) + dz)/2;
438 #else // MFEM_USE_MPFR is defined
440 HP_Quadrature1D hp_quad;
441 for (
int i = 1; i <= m; i++)
443 hp_quad.ComputeGaussLegendrePoint(n, i-1);
445 ir->
IntPoint(i-1).
x = hp_quad.GetPoint();
446 ir->
IntPoint(n-i).
x = hp_quad.GetSymmPoint();
450 #endif // MFEM_USE_MPFR
488 #ifndef MFEM_USE_MPFR
497 for (
int i = 1 ; i <= (np-1)/2 ; ++i)
501 double x_i = std::sin(M_PI * ((
double)(i)/(np-1) - 0.5));
502 double z_i = 0., p_l;
504 for (
int iter = 0 ; true ; ++iter)
510 for (
int l = 1 ; l < (np-1) ; ++l)
515 double p_lp1 = ( (2*l + 1)*x_i*p_l - l*p_lm1)/(l + 1);
533 double dx = (x_i*p_l - p_lm1) / (np*p_l);
534 if (std::abs(dx) < 1e-16)
538 z_i = ((1.0 + x_i) - dx)/2;
542 MFEM_VERIFY(iter < 8,
"np = " << np <<
", i = " << i
551 ip.
weight = (double)(1.0 / (np*(np-1)*p_l*p_l));
555 symm_ip.
x = 1.0 - z_i;
559 #else // MFEM_USE_MPFR is defined
561 HP_Quadrature1D hp_quad;
563 for (
int i = 0 ; i <= (np-1)/2 ; ++i)
565 hp_quad.ComputeGaussLobattoPoint(np, i);
567 ir->
IntPoint(np-1-i).
x = hp_quad.GetSymmPoint();
572 #endif // MFEM_USE_MPFR
584 for (
int i = 0; i < np ; ++i)
586 ir->
IntPoint(i).
x = double(i+1) / double(np + 1);
589 CalculateUniformWeights(ir, Quadrature1D::OpenUniform);
592 void QuadratureFunctions1D::ClosedUniform(
const int np,
603 for (
int i = 0; i < np ; ++i)
608 CalculateUniformWeights(ir, Quadrature1D::ClosedUniform);
617 for (
int i = 0; i < np ; ++i)
619 ir->
IntPoint(i).
x = double(2*i+1) / (2*np);
622 CalculateUniformWeights(ir, Quadrature1D::OpenHalfUniform);
635 GaussLegendre(np-1, &gl_ir);
637 for (
int i = 1; i < np-1; ++i)
643 CalculateUniformWeights(ir, Quadrature1D::ClosedGL);
646 void QuadratureFunctions1D::GivePolyPoints(
const int np,
double *
pts,
653 case Quadrature1D::GaussLegendre:
655 GaussLegendre(np,&ir);
658 case Quadrature1D::GaussLobatto:
660 GaussLobatto(np, &ir);
663 case Quadrature1D::OpenUniform:
668 case Quadrature1D::ClosedUniform:
670 ClosedUniform(np,&ir);
673 case Quadrature1D::OpenHalfUniform:
675 OpenHalfUniform(np, &ir);
678 case Quadrature1D::ClosedGL:
685 MFEM_ABORT(
"Asking for an unknown type of 1D Quadrature points, "
690 for (
int i = 0 ; i < np ; ++i)
696 void QuadratureFunctions1D::CalculateUniformWeights(
IntegrationRule *ir,
707 const int n = ir->
Size();
719 #ifndef MFEM_USE_MPFR
722 const IntegrationRule &glob_ir =
IntRules.
Get(Geometry::SEGMENT, n-1);
725 for (
int j = 0; j < n; j++)
729 Poly_1D::Basis basis(n-1, xv.GetData());
733 for (
int i = 0; i < m; i++)
735 const IntegrationPoint &ip = glob_ir.IntPoint(i);
736 basis.Eval(ip.x, xv);
737 w.Add(ip.weight, xv);
739 for (
int j = 0; j < n; j++)
744 #else // MFEM_USE_MPFR is defined
746 static const mpfr_rnd_t rnd = HP_Quadrature1D::rnd;
747 HP_Quadrature1D hp_quad;
748 mpfr_t l, lk, w0, wi, tmp, *weights;
749 mpfr_inits2(hp_quad.default_prec, l, lk, w0, wi, tmp, (mpfr_ptr) 0);
750 weights =
new mpfr_t[n];
751 for (
int i = 0; i < n; i++)
753 mpfr_init2(weights[i], hp_quad.default_prec);
754 mpfr_set_si(weights[i], 0, rnd);
756 hp_quad.SetRelTol(-48);
759 int hinv = 0, ihoffset = 0;
762 case Quadrature1D::ClosedUniform:
767 case Quadrature1D::OpenUniform:
772 case Quadrature1D::OpenHalfUniform:
778 MFEM_ABORT(
"invalid Quadrature1D type: " << type);
781 mpfr_fac_ui(w0, p, rnd);
782 mpfr_ui_pow_ui(tmp, hinv, p, rnd);
783 mpfr_div(w0, w0, tmp, rnd);
784 if (p%2) { mpfr_neg(w0, w0, rnd); }
786 for (
int j = 0; j < m; j++)
788 hp_quad.ComputeGaussLegendrePoint(m, j);
793 mpfr_mul_si(tmp, hp_quad.GetHPPoint(), hinv, rnd);
794 mpfr_sub_d(tmp, tmp, 0.5*ihoffset, rnd);
795 mpfr_round(tmp, tmp);
796 int k = min(max((
int)mpfr_get_si(tmp, rnd), 0), p);
797 mpfr_set_si(lk, 1, rnd);
798 for (
int i = 0; i <=
p; i++)
800 mpfr_set_si(tmp, 2*i+ihoffset, rnd);
801 mpfr_div_si(tmp, tmp, 2*hinv, rnd);
802 mpfr_sub(tmp, hp_quad.GetHPPoint(), tmp, rnd);
805 mpfr_mul(lk, lk, tmp, rnd);
809 mpfr_set(l, tmp, rnd);
812 mpfr_mul(l, l, lk, rnd);
813 mpfr_set(wi, w0, rnd);
814 for (
int i = 0;
true; i++)
819 mpfr_set_si(tmp, 2*i+ihoffset, rnd);
820 mpfr_div_si(tmp, tmp, 2*hinv, rnd);
821 mpfr_sub(tmp, hp_quad.GetHPPoint(), tmp, rnd);
822 mpfr_mul(tmp, tmp, wi, rnd);
823 mpfr_div(tmp, l, tmp, rnd);
828 mpfr_div(tmp, lk, wi, rnd);
831 mpfr_mul(tmp, tmp, hp_quad.GetHPWeight(), rnd);
832 mpfr_add(weights[i], weights[i], tmp, rnd);
834 if (i == p) {
break; }
837 mpfr_mul_si(wi, wi, i+1, rnd);
838 mpfr_div_si(wi, wi, i-p, rnd);
841 for (
int i = 0; i < n; i++)
844 mpfr_clear(weights[i]);
847 mpfr_clears(l, lk, w0, wi, tmp, (mpfr_ptr) 0);
849 #endif // MFEM_USE_MPFR
854 int Quadrature1D::CheckClosed(
int type)
866 int Quadrature1D::CheckOpen(
int type)
874 case OpenHalfUniform:
886 IntegrationRules::IntegrationRules(
int Ref,
int type_):
891 if (refined < 0) { own_rules = 0;
return; }
896 PointIntRules.SetSize(2, h_mt);
897 PointIntRules = NULL;
899 SegmentIntRules.SetSize(32, h_mt);
900 SegmentIntRules = NULL;
903 TriangleIntRules.SetSize(32, h_mt);
904 TriangleIntRules = NULL;
906 SquareIntRules.SetSize(32, h_mt);
907 SquareIntRules = NULL;
910 TetrahedronIntRules.SetSize(32, h_mt);
911 TetrahedronIntRules = NULL;
913 PrismIntRules.SetSize(32, h_mt);
914 PrismIntRules = NULL;
916 CubeIntRules.SetSize(32, h_mt);
934 mfem_error(
"IntegrationRules::Get(...) : Unknown geometry type!");
943 if (!HaveIntRule(*ir_array, Order))
945 #ifdef MFEM_USE_LEGACY_OPENMP
949 if (!HaveIntRule(*ir_array, Order))
952 int RealOrder = Order;
953 while (RealOrder+1 < ir_array->
Size() &&
954 (*ir_array)[RealOrder+1] == ir)
963 return *(*ir_array)[Order];
980 mfem_error(
"IntegrationRules::Set(...) : Unknown geometry type!");
984 if (HaveIntRule(*ir_array, Order))
986 MFEM_ABORT(
"Overwriting set rules is not supported!");
989 AllocIntRule(*ir_array, Order);
991 (*ir_array)[Order] = &IntRule;
1001 for (i = 0; i < ir_array.
Size(); i++)
1003 if (ir_array[i] != NULL && ir_array[i] != ir)
1013 if (!own_rules) {
return; }
1015 DeleteIntRuleArray(PointIntRules);
1016 DeleteIntRuleArray(SegmentIntRules);
1017 DeleteIntRuleArray(TriangleIntRules);
1018 DeleteIntRuleArray(SquareIntRules);
1019 DeleteIntRuleArray(TetrahedronIntRules);
1020 DeleteIntRuleArray(CubeIntRules);
1021 DeleteIntRuleArray(PrismIntRules);
1025 IntegrationRule *IntegrationRules::GenerateIntegrationRule(
int GeomType,
1031 return PointIntegrationRule(Order);
1033 return SegmentIntegrationRule(Order);
1035 return TriangleIntegrationRule(Order);
1037 return SquareIntegrationRule(Order);
1039 return TetrahedronIntegrationRule(Order);
1041 return CubeIntegrationRule(Order);
1043 return PrismIntegrationRule(Order);
1045 mfem_error(
"IntegrationRules::Set(...) : Unknown geometry type!");
1052 IntegrationRule *IntegrationRules::PointIntegrationRule(
int Order)
1056 mfem_error(
"Point Integration Rule of Order > 1 not defined");
1060 IntegrationRule *ir =
new IntegrationRule(1);
1061 ir->IntPoint(0).x = .0;
1062 ir->IntPoint(0).weight = 1.;
1064 PointIntRules[1] = PointIntRules[0] = ir;
1070 IntegrationRule *IntegrationRules::SegmentIntegrationRule(
int Order)
1072 int RealOrder = GetSegmentRealOrder(Order);
1074 AllocIntRule(SegmentIntRules, RealOrder);
1076 IntegrationRule tmp, *ir;
1077 ir = refined ? &tmp :
new IntegrationRule;
1121 MFEM_ABORT(
"unknown Quadrature1D type: " << quad_type);
1127 ir =
new IntegrationRule(2*n);
1128 for (
int j = 0; j < n; j++)
1130 ir->IntPoint(j).x = tmp.IntPoint(j).x/2.0;
1131 ir->IntPoint(j).weight = tmp.IntPoint(j).weight/2.0;
1132 ir->IntPoint(j+n).x = 0.5 + tmp.IntPoint(j).x/2.0;
1133 ir->IntPoint(j+n).weight = tmp.IntPoint(j).weight/2.0;
1136 SegmentIntRules[RealOrder-1] = SegmentIntRules[RealOrder] = ir;
1141 IntegrationRule *IntegrationRules::TriangleIntegrationRule(
int Order)
1143 IntegrationRule *ir = NULL;
1152 ir =
new IntegrationRule(1);
1153 ir->AddTriMidPoint(0, 0.5);
1154 TriangleIntRules[0] = TriangleIntRules[1] = ir;
1158 ir =
new IntegrationRule(3);
1159 ir->AddTriPoints3(0, 1./6., 1./6.);
1160 TriangleIntRules[2] = ir;
1165 ir =
new IntegrationRule(4);
1166 ir->AddTriMidPoint(0, -0.28125);
1167 ir->AddTriPoints3(1, 0.2, 25./96.);
1168 TriangleIntRules[3] = ir;
1172 ir =
new IntegrationRule(6);
1173 ir->AddTriPoints3(0, 0.091576213509770743460, 0.054975871827660933819);
1174 ir->AddTriPoints3(3, 0.44594849091596488632, 0.11169079483900573285);
1175 TriangleIntRules[4] = ir;
1179 ir =
new IntegrationRule(7);
1180 ir->AddTriMidPoint(0, 0.1125);
1181 ir->AddTriPoints3(1, 0.10128650732345633880, 0.062969590272413576298);
1182 ir->AddTriPoints3(4, 0.47014206410511508977, 0.066197076394253090369);
1183 TriangleIntRules[5] = ir;
1187 ir =
new IntegrationRule(12);
1188 ir->AddTriPoints3(0, 0.063089014491502228340, 0.025422453185103408460);
1189 ir->AddTriPoints3(3, 0.24928674517091042129, 0.058393137863189683013);
1190 ir->AddTriPoints6(6, 0.053145049844816947353, 0.31035245103378440542,
1191 0.041425537809186787597);
1192 TriangleIntRules[6] = ir;
1196 ir =
new IntegrationRule(12);
1197 ir->AddTriPoints3R(0, 0.062382265094402118174, 0.067517867073916085443,
1198 0.026517028157436251429);
1199 ir->AddTriPoints3R(3, 0.055225456656926611737, 0.32150249385198182267,
1200 0.043881408714446055037);
1202 ir->AddTriPoints3R(6, 0.034324302945097146470, 0.66094919618673565761,
1203 0.30472650086816719592, 0.028775042784981585738);
1204 ir->AddTriPoints3R(9, 0.51584233435359177926, 0.27771616697639178257,
1205 0.20644149867001643817, 0.067493187009802774463);
1206 TriangleIntRules[7] = ir;
1210 ir =
new IntegrationRule(16);
1211 ir->AddTriMidPoint(0, 0.0721578038388935841255455552445323);
1212 ir->AddTriPoints3(1, 0.170569307751760206622293501491464,
1213 0.0516086852673591251408957751460645);
1214 ir->AddTriPoints3(4, 0.0505472283170309754584235505965989,
1215 0.0162292488115990401554629641708902);
1216 ir->AddTriPoints3(7, 0.459292588292723156028815514494169,
1217 0.0475458171336423123969480521942921);
1218 ir->AddTriPoints6(10, 0.008394777409957605337213834539296,
1219 0.263112829634638113421785786284643,
1220 0.0136151570872174971324223450369544);
1221 TriangleIntRules[8] = ir;
1225 ir =
new IntegrationRule(19);
1226 ir->AddTriMidPoint(0, 0.0485678981413994169096209912536443);
1227 ir->AddTriPoints3b(1, 0.020634961602524744433,
1228 0.0156673501135695352684274156436046);
1229 ir->AddTriPoints3b(4, 0.12582081701412672546,
1230 0.0389137705023871396583696781497019);
1231 ir->AddTriPoints3(7, 0.188203535619032730240961280467335,
1232 0.0398238694636051265164458871320226);
1233 ir->AddTriPoints3(10, 0.0447295133944527098651065899662763,
1234 0.0127888378293490156308393992794999);
1235 ir->AddTriPoints6(13, 0.0368384120547362836348175987833851,
1236 0.2219629891607656956751025276931919,
1237 0.0216417696886446886446886446886446);
1238 TriangleIntRules[9] = ir;
1242 ir =
new IntegrationRule(25);
1243 ir->AddTriMidPoint(0, 0.0454089951913767900476432975500142);
1244 ir->AddTriPoints3b(1, 0.028844733232685245264984935583748,
1245 0.0183629788782333523585030359456832);
1246 ir->AddTriPoints3(4, 0.109481575485037054795458631340522,
1247 0.0226605297177639673913028223692986);
1248 ir->AddTriPoints6(7, 0.141707219414879954756683250476361,
1249 0.307939838764120950165155022930631,
1250 0.0363789584227100543021575883096803);
1251 ir->AddTriPoints6(13, 0.025003534762686386073988481007746,
1252 0.246672560639902693917276465411176,
1253 0.0141636212655287424183685307910495);
1254 ir->AddTriPoints6(19, 0.0095408154002994575801528096228873,
1255 0.0668032510122002657735402127620247,
1256 4.71083348186641172996373548344341E-03);
1257 TriangleIntRules[10] = ir;
1261 ir =
new IntegrationRule(28);
1262 ir->AddTriPoints6(0, 0.0,
1263 0.141129718717363295960826061941652,
1264 3.68119189165027713212944752369032E-03);
1265 ir->AddTriMidPoint(6, 0.0439886505811161193990465846607278);
1266 ir->AddTriPoints3(7, 0.0259891409282873952600324854988407,
1267 4.37215577686801152475821439991262E-03);
1268 ir->AddTriPoints3(10, 0.0942875026479224956305697762754049,
1269 0.0190407859969674687575121697178070);
1270 ir->AddTriPoints3b(13, 0.010726449965572372516734795387128,
1271 9.42772402806564602923839129555767E-03);
1272 ir->AddTriPoints3(16, 0.207343382614511333452934024112966,
1273 0.0360798487723697630620149942932315);
1274 ir->AddTriPoints3b(19, 0.122184388599015809877869236727746,
1275 0.0346645693527679499208828254519072);
1276 ir->AddTriPoints6(22, 0.0448416775891304433090523914688007,
1277 0.2772206675282791551488214673424523,
1278 0.0205281577146442833208261574536469);
1279 TriangleIntRules[11] = ir;
1283 ir =
new IntegrationRule(33);
1284 ir->AddTriPoints3b(0, 2.35652204523900E-02, 1.28655332202275E-02);
1285 ir->AddTriPoints3b(3, 1.20551215411079E-01, 2.18462722690190E-02);
1286 ir->AddTriPoints3(6, 2.71210385012116E-01, 3.14291121089425E-02);
1287 ir->AddTriPoints3(9, 1.27576145541586E-01, 1.73980564653545E-02);
1288 ir->AddTriPoints3(12, 2.13173504532100E-02, 3.08313052577950E-03);
1289 ir->AddTriPoints6(15, 1.15343494534698E-01, 2.75713269685514E-01,
1290 2.01857788831905E-02);
1291 ir->AddTriPoints6(21, 2.28383322222570E-02, 2.81325580989940E-01,
1292 1.11783866011515E-02);
1293 ir->AddTriPoints6(27, 2.57340505483300E-02, 1.16251915907597E-01,
1294 8.65811555432950E-03);
1295 TriangleIntRules[12] = ir;
1299 ir =
new IntegrationRule(37);
1300 ir->AddTriPoints3b(0, 0.0,
1301 2.67845189554543044455908674650066E-03);
1302 ir->AddTriMidPoint(3, 0.0293480398063595158995969648597808);
1303 ir->AddTriPoints3(4, 0.0246071886432302181878499494124643,
1304 3.92538414805004016372590903990464E-03);
1305 ir->AddTriPoints3b(7, 0.159382493797610632566158925635800,
1306 0.0253344765879434817105476355306468);
1307 ir->AddTriPoints3(10, 0.227900255506160619646298948153592,
1308 0.0250401630452545330803738542916538);
1309 ir->AddTriPoints3(13, 0.116213058883517905247155321839271,
1310 0.0158235572961491595176634480481793);
1311 ir->AddTriPoints3b(16, 0.046794039901841694097491569577008,
1312 0.0157462815379843978450278590138683);
1313 ir->AddTriPoints6(19, 0.0227978945382486125477207592747430,
1314 0.1254265183163409177176192369310890,
1315 7.90126610763037567956187298486575E-03);
1316 ir->AddTriPoints6(25, 0.0162757709910885409437036075960413,
1317 0.2909269114422506044621801030055257,
1318 7.99081889046420266145965132482933E-03);
1319 ir->AddTriPoints6(31, 0.0897330604516053590796290561145196,
1320 0.2723110556841851025078181617634414,
1321 0.0182757511120486476280967518782978);
1322 TriangleIntRules[13] = ir;
1326 ir =
new IntegrationRule(42);
1327 ir->AddTriPoints3b(0, 2.20721792756430E-02, 1.09417906847145E-02);
1328 ir->AddTriPoints3b(3, 1.64710561319092E-01, 1.63941767720625E-02);
1329 ir->AddTriPoints3(6, 2.73477528308839E-01, 2.58870522536460E-02);
1330 ir->AddTriPoints3(9, 1.77205532412543E-01, 2.10812943684965E-02);
1331 ir->AddTriPoints3(12, 6.17998830908730E-02, 7.21684983488850E-03);
1332 ir->AddTriPoints3(15, 1.93909612487010E-02, 2.46170180120000E-03);
1333 ir->AddTriPoints6(18, 5.71247574036480E-02, 1.72266687821356E-01,
1334 1.23328766062820E-02);
1335 ir->AddTriPoints6(24, 9.29162493569720E-02, 3.36861459796345E-01,
1336 1.92857553935305E-02);
1337 ir->AddTriPoints6(30, 1.46469500556540E-02, 2.98372882136258E-01,
1338 7.21815405676700E-03);
1339 ir->AddTriPoints6(36, 1.26833093287200E-03, 1.18974497696957E-01,
1340 2.50511441925050E-03);
1341 TriangleIntRules[14] = ir;
1345 ir =
new IntegrationRule(54);
1346 ir->AddTriPoints3b(0, 0.0834384072617499333, 0.016330909424402645);
1347 ir->AddTriPoints3b(3, 0.192779070841738867, 0.01370640901568218);
1348 ir->AddTriPoints3(6, 0.293197167913025367, 0.01325501829935165);
1349 ir->AddTriPoints3(9, 0.146467786942772933, 0.014607981068243055);
1350 ir->AddTriPoints3(12, 0.0563628676656034333, 0.005292304033121995);
1351 ir->AddTriPoints3(15, 0.0165751268583703333, 0.0018073215320460175);
1352 ir->AddTriPoints6(18, 0.0099122033092248, 0.239534554154794445,
1353 0.004263874050854718);
1354 ir->AddTriPoints6(24, 0.015803770630228, 0.404878807318339958,
1355 0.006958088258345965);
1356 ir->AddTriPoints6(30, 0.00514360881697066667, 0.0950021131130448885,
1357 0.0021459664703674175);
1358 ir->AddTriPoints6(36, 0.0489223257529888, 0.149753107322273969,
1359 0.008117664640887445);
1360 ir->AddTriPoints6(42, 0.0687687486325192, 0.286919612441334979,
1361 0.012803670460631195);
1362 ir->AddTriPoints6(48, 0.1684044181246992, 0.281835668099084562,
1363 0.016544097765822835);
1364 TriangleIntRules[15] = ir;
1369 ir =
new IntegrationRule(61);
1370 ir->AddTriMidPoint(0, 1.67185996454015E-02);
1371 ir->AddTriPoints3b(1, 5.65891888645200E-03, 2.54670772025350E-03);
1372 ir->AddTriPoints3b(4, 3.56473547507510E-02, 7.33543226381900E-03);
1373 ir->AddTriPoints3b(7, 9.95200619584370E-02, 1.21754391768360E-02);
1374 ir->AddTriPoints3b(10, 1.99467521245206E-01, 1.55537754344845E-02);
1375 ir->AddTriPoints3 (13, 2.52141267970953E-01, 1.56285556093100E-02);
1376 ir->AddTriPoints3 (16, 1.62047004658461E-01, 1.24078271698325E-02);
1377 ir->AddTriPoints3 (19, 7.58758822607460E-02, 7.02803653527850E-03);
1378 ir->AddTriPoints3 (22, 1.56547269678220E-02, 1.59733808688950E-03);
1379 ir->AddTriPoints6 (25, 1.01869288269190E-02, 3.34319867363658E-01,
1380 4.05982765949650E-03);
1381 ir->AddTriPoints6 (31, 1.35440871671036E-01, 2.92221537796944E-01,
1382 1.34028711415815E-02);
1383 ir->AddTriPoints6 (37, 5.44239242905830E-02, 3.19574885423190E-01,
1384 9.22999660541100E-03);
1385 ir->AddTriPoints6 (43, 1.28685608336370E-02, 1.90704224192292E-01,
1386 4.23843426716400E-03);
1387 ir->AddTriPoints6 (49, 6.71657824135240E-02, 1.80483211648746E-01,
1388 9.14639838501250E-03);
1389 ir->AddTriPoints6 (55, 1.46631822248280E-02, 8.07113136795640E-02,
1390 3.33281600208250E-03);
1391 TriangleIntRules[16] = TriangleIntRules[17] = ir;
1396 ir =
new IntegrationRule(73);
1397 ir->AddTriMidPoint(0, 0.0164531656944595);
1398 ir->AddTriPoints3b(1, 0.020780025853987, 0.005165365945636);
1399 ir->AddTriPoints3b(4, 0.090926214604215, 0.011193623631508);
1400 ir->AddTriPoints3b(7, 0.197166638701138, 0.015133062934734);
1401 ir->AddTriPoints3 (10, 0.255551654403098, 0.015245483901099);
1402 ir->AddTriPoints3 (13, 0.17707794215213, 0.0120796063708205);
1403 ir->AddTriPoints3 (16, 0.110061053227952, 0.0080254017934005);
1404 ir->AddTriPoints3 (19, 0.05552862425184, 0.004042290130892);
1405 ir->AddTriPoints3 (22, 0.012621863777229, 0.0010396810137425);
1406 ir->AddTriPoints6 (25, 0.003611417848412, 0.395754787356943,
1407 0.0019424384524905);
1408 ir->AddTriPoints6 (31, 0.13446675453078, 0.307929983880436,
1410 ir->AddTriPoints6 (37, 0.014446025776115, 0.26456694840652,
1412 ir->AddTriPoints6 (43, 0.046933578838178, 0.358539352205951,
1413 0.0080622733808655);
1414 ir->AddTriPoints6 (49, 0.002861120350567, 0.157807405968595,
1415 0.0012459709087455);
1416 ir->AddTriPoints6 (55, 0.075050596975911, 0.223861424097916,
1417 0.0091214200594755);
1418 ir->AddTriPoints6 (61, 0.03464707481676, 0.142421601113383,
1419 0.0051292818680995);
1420 ir->AddTriPoints6 (67, 0.065494628082938, 0.010161119296278,
1422 TriangleIntRules[18] = TriangleIntRules[19] = ir;
1426 ir =
new IntegrationRule(85);
1427 ir->AddTriMidPoint(0, 0.01380521349884976);
1428 ir->AddTriPoints3b(1, 0.001500649324429, 0.00088951477366337);
1429 ir->AddTriPoints3b(4, 0.0941397519389508667, 0.010056199056980585);
1430 ir->AddTriPoints3b(7, 0.2044721240895264, 0.013408923629665785);
1431 ir->AddTriPoints3(10, 0.264500202532787333, 0.012261566900751005);
1432 ir->AddTriPoints3(13, 0.211018964092076767, 0.008197289205347695);
1433 ir->AddTriPoints3(16, 0.107735607171271333, 0.0073979536993248);
1434 ir->AddTriPoints3(19, 0.0390690878378026667, 0.0022896411388521255);
1435 ir->AddTriPoints3(22, 0.0111743797293296333, 0.0008259132577881085);
1436 ir->AddTriPoints6(25, 0.00534961818733726667, 0.0635496659083522206,
1437 0.001174585454287792);
1438 ir->AddTriPoints6(31, 0.00795481706619893333, 0.157106918940706982,
1439 0.0022329628770908965);
1440 ir->AddTriPoints6(37, 0.0104223982812638, 0.395642114364374018,
1441 0.003049783403953986);
1442 ir->AddTriPoints6(43, 0.0109644147961233333, 0.273167570712910522,
1443 0.0034455406635941015);
1444 ir->AddTriPoints6(49, 0.0385667120854623333, 0.101785382485017108,
1445 0.0039987375362390815);
1446 ir->AddTriPoints6(55, 0.0355805078172182, 0.446658549176413815,
1447 0.003693067142668012);
1448 ir->AddTriPoints6(61, 0.0496708163627641333, 0.199010794149503095,
1449 0.00639966593932413);
1450 ir->AddTriPoints6(67, 0.0585197250843317333, 0.3242611836922827,
1451 0.008629035587848275);
1452 ir->AddTriPoints6(73, 0.121497787004394267, 0.208531363210132855,
1453 0.009336472951467735);
1454 ir->AddTriPoints6(79, 0.140710844943938733, 0.323170566536257485,
1455 0.01140911202919763);
1456 TriangleIntRules[20] = ir;
1464 ir =
new IntegrationRule(126);
1465 ir->AddTriPoints3b(0, 0.0279464830731742, 0.0040027909400102085);
1466 ir->AddTriPoints3b(3, 0.131178601327651467, 0.00797353841619525);
1467 ir->AddTriPoints3b(6, 0.220221729512072267, 0.006554570615397765);
1468 ir->AddTriPoints3 (9, 0.298443234019804467, 0.00979150048281781);
1469 ir->AddTriPoints3(12, 0.2340441723373718, 0.008235442720768635);
1470 ir->AddTriPoints3(15, 0.151468334609017567, 0.00427363953704605);
1471 ir->AddTriPoints3(18, 0.112733893545993667, 0.004080942928613246);
1472 ir->AddTriPoints3(21, 0.0777156920915263, 0.0030605732699918895);
1473 ir->AddTriPoints3(24, 0.034893093614297, 0.0014542491324683325);
1474 ir->AddTriPoints3(27, 0.00725818462093236667, 0.00034613762283099815);
1475 ir->AddTriPoints6(30, 0.0012923527044422, 0.227214452153364077,
1476 0.0006241445996386985);
1477 ir->AddTriPoints6(36, 0.0053997012721162, 0.435010554853571706,
1478 0.001702376454401511);
1479 ir->AddTriPoints6(42, 0.006384003033975, 0.320309599272204437,
1480 0.0016798271630320255);
1481 ir->AddTriPoints6(48, 0.00502821150199306667, 0.0917503222800051889,
1482 0.000858078269748377);
1483 ir->AddTriPoints6(54, 0.00682675862178186667, 0.0380108358587243835,
1484 0.000740428158357803);
1485 ir->AddTriPoints6(60, 0.0100161996399295333, 0.157425218485311668,
1486 0.0017556563053643425);
1487 ir->AddTriPoints6(66, 0.02575781317339, 0.239889659778533193,
1488 0.003696775074853242);
1489 ir->AddTriPoints6(72, 0.0302278981199158, 0.361943118126060531,
1490 0.003991543738688279);
1491 ir->AddTriPoints6(78, 0.0305049901071620667, 0.0835519609548285602,
1492 0.0021779813065790205);
1493 ir->AddTriPoints6(84, 0.0459565473625693333, 0.148443220732418205,
1494 0.003682528350708916);
1495 ir->AddTriPoints6(90, 0.0674428005402775333, 0.283739708727534955,
1496 0.005481786423209775);
1497 ir->AddTriPoints6(96, 0.0700450914159106, 0.406899375118787573,
1498 0.00587498087177056);
1499 ir->AddTriPoints6(102, 0.0839115246401166, 0.194113987024892542,
1500 0.005007800356899285);
1501 ir->AddTriPoints6(108, 0.120375535677152667, 0.32413434700070316,
1502 0.00665482039381434);
1503 ir->AddTriPoints6(114, 0.148066899157366667, 0.229277483555980969,
1504 0.00707722325261307);
1505 ir->AddTriPoints6(120, 0.191771865867325067, 0.325618122595983752,
1506 0.007440689780584005);
1507 TriangleIntRules[21] =
1508 TriangleIntRules[22] =
1509 TriangleIntRules[23] =
1510 TriangleIntRules[24] =
1511 TriangleIntRules[25] = ir;
1516 int i = (Order / 2) * 2 + 1;
1517 AllocIntRule(TriangleIntRules, i);
1518 ir =
new IntegrationRule;
1519 ir->GrundmannMollerSimplexRule(i/2,2);
1520 TriangleIntRules[i-1] = TriangleIntRules[i] = ir;
1526 IntegrationRule *IntegrationRules::SquareIntegrationRule(
int Order)
1528 int RealOrder = GetSegmentRealOrder(Order);
1530 if (!HaveIntRule(SegmentIntRules, RealOrder))
1532 SegmentIntegrationRule(RealOrder);
1534 AllocIntRule(SquareIntRules, RealOrder);
1535 SquareIntRules[RealOrder-1] =
1536 SquareIntRules[RealOrder] =
1537 new IntegrationRule(*SegmentIntRules[RealOrder],
1538 *SegmentIntRules[RealOrder]);
1539 return SquareIntRules[Order];
1544 IntegrationRule *IntegrationRules::TetrahedronIntegrationRule(
int Order)
1546 IntegrationRule *ir;
1555 ir =
new IntegrationRule(1);
1556 ir->AddTetMidPoint(0, 1./6.);
1557 TetrahedronIntRules[0] = TetrahedronIntRules[1] = ir;
1561 ir =
new IntegrationRule(4);
1563 ir->AddTetPoints4b(0, 0.58541019662496845446, 1./24.);
1564 TetrahedronIntRules[2] = ir;
1568 ir =
new IntegrationRule(5);
1569 ir->AddTetMidPoint(0, -2./15.);
1570 ir->AddTetPoints4b(1, 0.5, 0.075);
1571 TetrahedronIntRules[3] = ir;
1575 ir =
new IntegrationRule(11);
1576 ir->AddTetPoints4(0, 1./14., 343./45000.);
1577 ir->AddTetMidPoint(4, -74./5625.);
1578 ir->AddTetPoints6(5, 0.10059642383320079500, 28./1125.);
1579 TetrahedronIntRules[4] = ir;
1583 ir =
new IntegrationRule(14);
1584 ir->AddTetPoints6(0, 0.045503704125649649492,
1585 7.0910034628469110730E-03);
1586 ir->AddTetPoints4(6, 0.092735250310891226402, 0.012248840519393658257);
1587 ir->AddTetPoints4b(10, 0.067342242210098170608,
1588 0.018781320953002641800);
1589 TetrahedronIntRules[5] = ir;
1593 ir =
new IntegrationRule(24);
1594 ir->AddTetPoints4(0, 0.21460287125915202929,
1595 6.6537917096945820166E-03);
1596 ir->AddTetPoints4(4, 0.040673958534611353116,
1597 1.6795351758867738247E-03);
1598 ir->AddTetPoints4b(8, 0.032986329573173468968,
1599 9.2261969239424536825E-03);
1600 ir->AddTetPoints12(12, 0.063661001875017525299, 0.26967233145831580803,
1601 8.0357142857142857143E-03);
1602 TetrahedronIntRules[6] = ir;
1606 ir =
new IntegrationRule(31);
1607 ir->AddTetPoints6(0, 0.0, 9.7001763668430335097E-04);
1608 ir->AddTetMidPoint(6, 0.018264223466108820291);
1609 ir->AddTetPoints4(7, 0.078213192330318064374, 0.010599941524413686916);
1610 ir->AddTetPoints4(11, 0.12184321666390517465,
1611 -0.062517740114331851691);
1612 ir->AddTetPoints4b(15, 2.3825066607381275412E-03,
1613 4.8914252630734993858E-03);
1614 ir->AddTetPoints12(19, 0.1, 0.2, 0.027557319223985890653);
1615 TetrahedronIntRules[7] = ir;
1619 ir =
new IntegrationRule(43);
1620 ir->AddTetPoints4(0, 5.7819505051979972532E-03,
1621 1.6983410909288737984E-04);
1622 ir->AddTetPoints4(4, 0.082103588310546723091,
1623 1.9670333131339009876E-03);
1624 ir->AddTetPoints12(8, 0.036607749553197423679, 0.19048604193463345570,
1625 2.1405191411620925965E-03);
1626 ir->AddTetPoints6(20, 0.050532740018894224426,
1627 4.5796838244672818007E-03);
1628 ir->AddTetPoints12(26, 0.22906653611681113960, 0.035639582788534043717,
1629 5.7044858086819185068E-03);
1630 ir->AddTetPoints4(38, 0.20682993161067320408, 0.014250305822866901248);
1631 ir->AddTetMidPoint(42, -0.020500188658639915841);
1632 TetrahedronIntRules[8] = ir;
1636 ir =
new IntegrationRule;
1637 ir->GrundmannMollerSimplexRule(4,3);
1638 TetrahedronIntRules[9] = ir;
1642 int i = (Order / 2) * 2 + 1;
1643 AllocIntRule(TetrahedronIntRules, i);
1644 ir =
new IntegrationRule;
1645 ir->GrundmannMollerSimplexRule(i/2,3);
1646 TetrahedronIntRules[i-1] = TetrahedronIntRules[i] = ir;
1652 IntegrationRule *IntegrationRules::PrismIntegrationRule(
int Order)
1656 int nt = irt.GetNPoints();
1657 int ns = irs.GetNPoints();
1658 AllocIntRule(PrismIntRules, Order);
1659 PrismIntRules[Order] =
new IntegrationRule(nt * ns);
1661 for (
int ks=0; ks<ns; ks++)
1663 const IntegrationPoint & ips = irs.IntPoint(ks);
1664 for (
int kt=0; kt<nt; kt++)
1666 int kp = ks * nt + kt;
1667 const IntegrationPoint & ipt = irt.IntPoint(kt);
1668 IntegrationPoint & ipp = PrismIntRules[Order]->IntPoint(kp);
1672 ipp.weight = ipt.weight * ips.weight;
1675 return PrismIntRules[Order];
1679 IntegrationRule *IntegrationRules::CubeIntegrationRule(
int Order)
1681 int RealOrder = GetSegmentRealOrder(Order);
1682 if (!HaveIntRule(SegmentIntRules, RealOrder))
1684 SegmentIntegrationRule(RealOrder);
1686 AllocIntRule(CubeIntRules, RealOrder);
1687 CubeIntRules[RealOrder-1] =
1688 CubeIntRules[RealOrder] =
1689 new IntegrationRule(*SegmentIntRules[RealOrder],
1690 *SegmentIntRules[RealOrder],
1691 *SegmentIntRules[RealOrder]);
1692 return CubeIntRules[Order];
int GetNPoints() const
Returns the number of the points in the integration rule.
int Size() const
Return the logical size of the array.
Class for an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
static void ClosedUniform(const int np, IntegrationRule *ir)
~IntegrationRules()
Destroys an IntegrationRules object.
Container class for integration rules.
static void OpenUniform(const int np, IntegrationRule *ir)
static void GaussLobatto(const int np, IntegrationRule *ir)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void SetPointIndices()
Sets the indices of each quadrature point on initialization.
double p(const Vector &x, double t)
static void OpenHalfUniform(const int np, IntegrationRule *ir)
MemoryType
Memory types supported by MFEM.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void SetOrder(const int order)
Sets the order of the integration rule. This is only for keeping order information, it does not alter any data in the IntegrationRule.
Class for integration point with weight.
Host memory; using new[] and delete[].
aka "open half" Newton-Cotes
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
void pts(int iphi, int t, double x[])
void Set(int GeomType, int Order, IntegrationRule &IntRule)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
void Set1w(const double x1, const double w)
static void GaussLegendre(const int np, IntegrationRule *ir)
double f(const Vector &p)