MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
intrules.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 // Implementation of IntegrationRule(s) classes
13 
14 // Acknowledgment: Some of the high-precision triangular and tetrahedral
15 // quadrature rules below were obtained from the Encyclopaedia of Cubature
16 // Formulas at http://nines.cs.kuleuven.be/research/ecf/ecf.html
17 
18 #include "fem.hpp"
19 #include <cmath>
20 
21 #ifdef MFEM_USE_MPFR
22 #include <mpfr.h>
23 #endif
24 
25 using namespace std;
26 
27 namespace mfem
28 {
29 
30 IntegrationRule::IntegrationRule(IntegrationRule &irx, IntegrationRule &iry)
31 {
32  int i, j, nx, ny;
33 
34  nx = irx.GetNPoints();
35  ny = iry.GetNPoints();
36  SetSize(nx * ny);
37 
38  for (j = 0; j < ny; j++)
39  {
40  IntegrationPoint &ipy = iry.IntPoint(j);
41  for (i = 0; i < nx; i++)
42  {
43  IntegrationPoint &ipx = irx.IntPoint(i);
44  IntegrationPoint &ip = IntPoint(j*nx+i);
45 
46  ip.x = ipx.x;
47  ip.y = ipy.x;
48  ip.weight = ipx.weight * ipy.weight;
49  }
50  }
51 
52  SetPointIndices();
53 }
54 
55 IntegrationRule::IntegrationRule(IntegrationRule &irx, IntegrationRule &iry,
56  IntegrationRule &irz)
57 {
58  const int nx = irx.GetNPoints();
59  const int ny = iry.GetNPoints();
60  const int nz = irz.GetNPoints();
61  SetSize(nx*ny*nz);
62 
63  for (int iz = 0; iz < nz; ++iz)
64  {
65  IntegrationPoint &ipz = irz.IntPoint(iz);
66  for (int iy = 0; iy < ny; ++iy)
67  {
68  IntegrationPoint &ipy = iry.IntPoint(iy);
69  for (int ix = 0; ix < nx; ++ix)
70  {
71  IntegrationPoint &ipx = irx.IntPoint(ix);
72  IntegrationPoint &ip = IntPoint(iz*nx*ny + iy*nx + ix);
73 
74  ip.x = ipx.x;
75  ip.y = ipy.x;
76  ip.z = ipz.x;
77  ip.weight = ipx.weight*ipy.weight*ipz.weight;
78  }
79  }
80  }
81 
82  SetPointIndices();
83 }
84 
85 const Array<double> &IntegrationRule::GetWeights() const
86 {
87  if (weights.Size() != GetNPoints())
88  {
89  weights.SetSize(GetNPoints());
90  for (int i = 0; i < GetNPoints(); i++)
91  {
92  weights[i] = IntPoint(i).weight;
93  }
94  }
95  return weights;
96 }
97 
98 void IntegrationRule::SetPointIndices()
99 {
100  for (int i = 0; i < Size(); i++)
101  {
102  IntPoint(i).index = i;
103  }
104 }
105 
106 void IntegrationRule::GrundmannMollerSimplexRule(int s, int n)
107 {
108  // for pow on older compilers
109  using std::pow;
110  const int d = 2*s + 1;
111  Vector fact(d + n + 1);
112  Array<int> beta(n), sums(n);
113 
114  fact(0) = 1.;
115  for (int i = 1; i < fact.Size(); i++)
116  {
117  fact(i) = fact(i - 1)*i;
118  }
119 
120  // number of points is \binom{n + s + 1}{n + 1}
121  int np = 1, f = 1;
122  for (int i = 0; i <= n; i++)
123  {
124  np *= (s + i + 1), f *= (i + 1);
125  }
126  np /= f;
127  SetSize(np);
128 
129  int pt = 0;
130  for (int i = 0; i <= s; i++)
131  {
132  double weight;
133 
134  weight = pow(2., -2*s)*pow(static_cast<double>(d + n - 2*i),
135  d)/fact(i)/fact(d + n - i);
136  if (i%2)
137  {
138  weight = -weight;
139  }
140 
141  // loop over all beta : beta_0 + ... + beta_{n-1} <= s - i
142  int k = s - i;
143  beta = 0;
144  sums = 0;
145  while (true)
146  {
147  IntegrationPoint &ip = IntPoint(pt++);
148  ip.weight = weight;
149  ip.x = double(2*beta[0] + 1)/(d + n - 2*i);
150  ip.y = double(2*beta[1] + 1)/(d + n - 2*i);
151  if (n == 3)
152  {
153  ip.z = double(2*beta[2] + 1)/(d + n - 2*i);
154  }
155 
156  int j = 0;
157  while (sums[j] == k)
158  {
159  beta[j++] = 0;
160  if (j == n)
161  {
162  goto done_beta;
163  }
164  }
165  beta[j]++;
166  sums[j]++;
167  for (j--; j >= 0; j--)
168  {
169  sums[j] = sums[j+1];
170  }
171  }
172  done_beta:
173  ;
174  }
175 }
176 
177 
178 #ifdef MFEM_USE_MPFR
179 
180 // Class for computing hi-precision (HP) quadrature in 1D
181 class HP_Quadrature1D
182 {
183 protected:
184  mpfr_t pi, z, pp, p1, p2, p3, dz, w, rtol;
185 
186 public:
187  static const mpfr_rnd_t rnd = GMP_RNDN;
188  static const int default_prec = 128;
189 
190  // prec = MPFR precision in bits
191  HP_Quadrature1D(const int prec = default_prec)
192  {
193  mpfr_inits2(prec, pi, z, pp, p1, p2, p3, dz, w, rtol, (mpfr_ptr) 0);
194  mpfr_const_pi(pi, rnd);
195  mpfr_set_si_2exp(rtol, 1, -32, rnd); // 2^(-32) < 2.33e-10
196  }
197 
198  // set rtol = 2^exponent
199  // this is a tolerance for the last correction of x_i in Newton's algorithm;
200  // this gives roughly rtol^2 accuracy for the final x_i.
201  void SetRelTol(const int exponent = -32)
202  {
203  mpfr_set_si_2exp(rtol, 1, exponent, rnd);
204  }
205 
206  // n - number of quadrature points
207  // k - index of the point to compute, 0 <= k < n
208  // see also: QuadratureFunctions1D::GaussLegendre
209  void ComputeGaussLegendrePoint(const int n, const int k)
210  {
211  MFEM_ASSERT(n > 0 && 0 <= k && k < n, "invalid n = " << n
212  << " and/or k = " << k);
213 
214  int i = (k < (n+1)/2) ? k+1 : n-k;
215 
216  // Initial guess for the x-coordinate:
217  // set z = cos(pi * (i - 0.25) / (n + 0.5)) =
218  // = sin(pi * ((n+1-2*i) / (2*n+1)))
219  mpfr_set_si(z, n+1-2*i, rnd);
220  mpfr_div_si(z, z, 2*n+1, rnd);
221  mpfr_mul(z, z, pi, rnd);
222  mpfr_sin(z, z, rnd);
223 
224  bool done = false;
225  while (1)
226  {
227  mpfr_set_si(p2, 1, rnd);
228  mpfr_set(p1, z, rnd);
229  for (int j = 2; j <= n; j++)
230  {
231  mpfr_set(p3, p2, rnd);
232  mpfr_set(p2, p1, rnd);
233  // p1 = ((2 * j - 1) * z * p2 - (j - 1) * p3) / j;
234  mpfr_mul_si(p1, z, 2*j-1, rnd);
235  mpfr_mul_si(p3, p3, j-1, rnd);
236  mpfr_fms(p1, p1, p2, p3, rnd);
237  mpfr_div_si(p1, p1, j, rnd);
238  }
239  // p1 is Legendre polynomial
240 
241  // derivative of the Legendre polynomial:
242  // pp = n * (z*p1-p2) / (z*z - 1);
243  mpfr_fms(pp, z, p1, p2, rnd);
244  mpfr_mul_si(pp, pp, n, rnd);
245  mpfr_sqr(p2, z, rnd);
246  mpfr_sub_si(p2, p2, 1, rnd);
247  mpfr_div(pp, pp, p2, rnd);
248 
249  if (done) { break; }
250 
251  // set delta_z: dz = p1/pp;
252  mpfr_div(dz, p1, pp, rnd);
253  // compute absolute tolerance: atol = rtol*(1-z)
254  mpfr_t &atol = w;
255  mpfr_si_sub(atol, 1, z, rnd);
256  mpfr_mul(atol, atol, rtol, rnd);
257  if (mpfr_cmpabs(dz, atol) <= 0)
258  {
259  done = true;
260  // continue the computation: get pp at the new point, then exit
261  }
262  // update z = z - dz
263  mpfr_sub(z, z, dz, rnd);
264  }
265 
266  // map z to (0,1): z = (1 - z)/2
267  mpfr_si_sub(z, 1, z, rnd);
268  mpfr_div_2si(z, z, 1, rnd);
269 
270  // weight: w = 1/(4*z*(1 - z)*pp*pp)
271  mpfr_sqr(w, pp, rnd);
272  mpfr_mul_2si(w, w, 2, rnd);
273  mpfr_mul(w, w, z, rnd);
274  mpfr_si_sub(p1, 1, z, rnd); // p1 = 1-z
275  mpfr_mul(w, w, p1, rnd);
276  mpfr_si_div(w, 1, w, rnd);
277 
278  if (k >= (n+1)/2) { mpfr_swap(z, p1); }
279  }
280 
281  // n - number of quadrature points
282  // k - index of the point to compute, 0 <= k < n
283  // see also: QuadratureFunctions1D::GaussLobatto
284  void ComputeGaussLobattoPoint(const int n, const int k)
285  {
286  MFEM_ASSERT(n > 1 && 0 <= k && k < n, "invalid n = " << n
287  << " and/or k = " << k);
288 
289  int i = (k < (n+1)/2) ? k : n-1-k;
290 
291  if (i == 0)
292  {
293  mpfr_set_si(z, 0, rnd);
294  mpfr_set_si(p1, 1, rnd);
295  mpfr_set_si(w, n*(n-1), rnd);
296  mpfr_si_div(w, 1, w, rnd); // weight = 1/(n*(n-1))
297  return;
298  }
299  // initial guess is the corresponding Chebyshev point, z:
300  // z = -cos(pi * i/(n-1)) = sin(pi * (2*i-n+1)/(2*n-2))
301  mpfr_set_si(z, 2*i-n+1, rnd);
302  mpfr_div_si(z, z, 2*(n-1), rnd);
303  mpfr_mul(z, pi, z, rnd);
304  mpfr_sin(z, z, rnd);
305  bool done = false;
306  for (int iter = 0 ; true ; ++iter)
307  {
308  // build Legendre polynomials, up to P_{n}(z)
309  mpfr_set_si(p1, 1, rnd);
310  mpfr_set(p2, z, rnd);
311 
312  for (int l = 1 ; l < (n-1) ; ++l)
313  {
314  // P_{l+1}(x) = [ (2*l+1)*x*P_l(x) - l*P_{l-1}(x) ]/(l+1)
315  mpfr_mul_si(p1, p1, l, rnd);
316  mpfr_mul_si(p3, z, 2*l+1, rnd);
317  mpfr_fms(p3, p3, p2, p1, rnd);
318  mpfr_div_si(p3, p3, l+1, rnd);
319 
320  mpfr_set(p1, p2, rnd);
321  mpfr_set(p2, p3, rnd);
322  }
323  if (done) { break; }
324  // compute dz = resid/deriv = (z*p2 - p1) / (n*p2);
325  mpfr_fms(dz, z, p2, p1, rnd);
326  mpfr_mul_si(p3, p2, n, rnd);
327  mpfr_div(dz, dz, p3, rnd);
328  // update: z = z - dz
329  mpfr_sub(z, z, dz, rnd);
330  // compute absolute tolerance: atol = rtol*(1 + z)
331  mpfr_t &atol = w;
332  mpfr_add_si(atol, z, 1, rnd);
333  mpfr_mul(atol, atol, rtol, rnd);
334  // check for convergence
335  if (mpfr_cmpabs(dz, atol) <= 0)
336  {
337  done = true;
338  // continue the computation: get p2 at the new point, then exit
339  }
340  // If the iteration does not converge fast, something is wrong.
341  MFEM_VERIFY(iter < 8, "n = " << n << ", i = " << i
342  << ", dz = " << mpfr_get_d(dz, rnd));
343  }
344  // Map to the interval [0,1] and scale the weights
345  mpfr_add_si(z, z, 1, rnd);
346  mpfr_div_2si(z, z, 1, rnd);
347  // set the symmetric point
348  mpfr_si_sub(p1, 1, z, rnd);
349  // w = 1/[ n*(n-1)*[P_{n-1}(z)]^2 ]
350  mpfr_sqr(w, p2, rnd);
351  mpfr_mul_si(w, w, n*(n-1), rnd);
352  mpfr_si_div(w, 1, w, rnd);
353 
354  if (k >= (n+1)/2) { mpfr_swap(z, p1); }
355  }
356 
357  double GetPoint() const { return mpfr_get_d(z, rnd); }
358  double GetSymmPoint() const { return mpfr_get_d(p1, rnd); }
359  double GetWeight() const { return mpfr_get_d(w, rnd); }
360 
361  const mpfr_t &GetHPPoint() const { return z; }
362  const mpfr_t &GetHPSymmPoint() const { return p1; }
363  const mpfr_t &GetHPWeight() const { return w; }
364 
365  ~HP_Quadrature1D()
366  {
367  mpfr_clears(pi, z, pp, p1, p2, p3, dz, w, rtol, (mpfr_ptr) 0);
368  mpfr_free_cache();
369  }
370 };
371 
372 #endif // MFEM_USE_MPFR
373 
374 
375 void QuadratureFunctions1D::GaussLegendre(const int np, IntegrationRule* ir)
376 {
377  ir->SetSize(np);
378 
379  switch (np)
380  {
381  case 1:
382  ir->IntPoint(0).Set1w(0.5, 1.0);
383  return;
384  case 2:
385  ir->IntPoint(0).Set1w(0.21132486540518711775, 0.5);
386  ir->IntPoint(1).Set1w(0.78867513459481288225, 0.5);
387  return;
388  case 3:
389  ir->IntPoint(0).Set1w(0.11270166537925831148, 5./18.);
390  ir->IntPoint(1).Set1w(0.5, 4./9.);
391  ir->IntPoint(2).Set1w(0.88729833462074168852, 5./18.);
392  return;
393  }
394 
395  const int n = np;
396  const int m = (n+1)/2;
397 
398 #ifndef MFEM_USE_MPFR
399 
400  for (int i = 1; i <= m; i++)
401  {
402  double z = cos(M_PI * (i - 0.25) / (n + 0.5));
403  double pp, p1, dz, xi = 0.;
404  bool done = false;
405  while (1)
406  {
407  double p2 = 1;
408  p1 = z;
409  for (int j = 2; j <= n; j++)
410  {
411  double p3 = p2;
412  p2 = p1;
413  p1 = ((2 * j - 1) * z * p2 - (j - 1) * p3) / j;
414  }
415  // p1 is Legendre polynomial
416 
417  pp = n * (z*p1-p2) / (z*z - 1);
418  if (done) { break; }
419 
420  dz = p1/pp;
421  if (fabs(dz) < 1e-16)
422  {
423  done = true;
424  // map the new point (z-dz) to (0,1):
425  xi = ((1 - z) + dz)/2; // (1 - (z - dz))/2 has bad round-off
426  // continue the computation: get pp at the new point, then exit
427  }
428  // update: z = z - dz
429  z -= dz;
430  }
431 
432  ir->IntPoint(i-1).x = xi;
433  ir->IntPoint(n-i).x = 1 - xi;
434  ir->IntPoint(i-1).weight =
435  ir->IntPoint(n-i).weight = 1./(4*xi*(1 - xi)*pp*pp);
436  }
437 
438 #else // MFEM_USE_MPFR is defined
439 
440  HP_Quadrature1D hp_quad;
441  for (int i = 1; i <= m; i++)
442  {
443  hp_quad.ComputeGaussLegendrePoint(n, i-1);
444 
445  ir->IntPoint(i-1).x = hp_quad.GetPoint();
446  ir->IntPoint(n-i).x = hp_quad.GetSymmPoint();
447  ir->IntPoint(i-1).weight = ir->IntPoint(n-i).weight = hp_quad.GetWeight();
448  }
449 
450 #endif // MFEM_USE_MPFR
451 
452 }
453 
454 void QuadratureFunctions1D::GaussLobatto(const int np, IntegrationRule* ir)
455 {
456  /* An np point Gauss-Lobatto quadrature has (np - 2) free abscissa the other
457  (2) abscissa are the interval endpoints.
458 
459  The interior x_i are the zeros of P'_{np-1}(x). The weights of the
460  interior points on the interval [-1,1] are:
461 
462  w_i = 2/(np*(np-1)*[P_{np-1}(x_i)]^2)
463 
464  The end point weights (on [-1,1]) are: w_{end} = 2/(np*(np-1)).
465 
466  The interior abscissa are found via a nonlinear solve, the initial guess
467  for each point is the corresponding Chebyshev point.
468 
469  After we find all points on the interval [-1,1], we will map and scale the
470  points and weights to the MFEM natural interval [0,1].
471 
472  References:
473  [1] E. E. Lewis and W. F. Millier, "Computational Methods of Neutron
474  Transport", Appendix A
475  [2] the QUADRULE software by John Burkardt,
476  https://people.sc.fsu.edu/~jburkardt/cpp_src/quadrule/quadrule.cpp
477  */
478 
479  ir->SetSize(np);
480  if ( np == 1 )
481  {
482  ir->IntPoint(0).Set1w(0.5, 1.0);
483  }
484  else
485  {
486 
487 #ifndef MFEM_USE_MPFR
488 
489  // endpoints and respective weights
490  ir->IntPoint(0).x = 0.0;
491  ir->IntPoint(np-1).x = 1.0;
492  ir->IntPoint(0).weight = ir->IntPoint(np-1).weight = 1.0/(np*(np-1));
493 
494  // interior points and weights
495  // use symmetry and compute just half of the points
496  for (int i = 1 ; i <= (np-1)/2 ; ++i)
497  {
498  // initial guess is the corresponding Chebyshev point, x_i:
499  // x_i = -cos(\pi * (i / (np-1)))
500  double x_i = std::sin(M_PI * ((double)(i)/(np-1) - 0.5));
501  double z_i = 0., p_l;
502  bool done = false;
503  for (int iter = 0 ; true ; ++iter)
504  {
505  // build Legendre polynomials, up to P_{np}(x_i)
506  double p_lm1 = 1.0;
507  p_l = x_i;
508 
509  for (int l = 1 ; l < (np-1) ; ++l)
510  {
511  // The Legendre polynomials can be built by recursion:
512  // x * P_l(x) = 1/(2*l+1)*[ (l+1)*P_{l+1}(x) + l*P_{l-1} ], i.e.
513  // P_{l+1}(x) = [ (2*l+1)*x*P_l(x) - l*P_{l-1} ]/(l+1)
514  double p_lp1 = ( (2*l + 1)*x_i*p_l - l*p_lm1)/(l + 1);
515 
516  p_lm1 = p_l;
517  p_l = p_lp1;
518  }
519  if (done) { break; }
520  // after this loop, p_l holds P_{np-1}(x_i)
521  // resid = (x^2-1)*P'_{np-1}(x_i)
522  // but use the recurrence relationship
523  // (x^2 -1)P'_l(x) = l*[ x*P_l(x) - P_{l-1}(x) ]
524  // thus, resid = (np-1) * (x_i*p_l - p_lm1)
525 
526  // The derivative of the residual is:
527  // \frac{d}{d x} \left[ (x^2 -1)P'_l(x) ] \right] =
528  // l * (l+1) * P_l(x), with l = np-1,
529  // therefore, deriv = np * (np-1) * p_l;
530 
531  // compute dx = resid/deriv
532  double dx = (x_i*p_l - p_lm1) / (np*p_l);
533  if (std::abs(dx) < 1e-16)
534  {
535  done = true;
536  // Map the point to the interval [0,1]
537  z_i = ((1.0 + x_i) - dx)/2;
538  // continue the computation: get p_l at the new point, then exit
539  }
540  // If the iteration does not converge fast, something is wrong.
541  MFEM_VERIFY(iter < 8, "np = " << np << ", i = " << i
542  << ", dx = " << dx);
543  // update x_i:
544  x_i -= dx;
545  }
546  // Map to the interval [0,1] and scale the weights
547  IntegrationPoint &ip = ir->IntPoint(i);
548  ip.x = z_i;
549  // w_i = (2/[ n*(n-1)*[P_{n-1}(x_i)]^2 ]) / 2
550  ip.weight = (double)(1.0 / (np*(np-1)*p_l*p_l));
551 
552  // set the symmetric point
553  IntegrationPoint &symm_ip = ir->IntPoint(np-1-i);
554  symm_ip.x = 1.0 - z_i;
555  symm_ip.weight = ip.weight;
556  }
557 
558 #else // MFEM_USE_MPFR is defined
559 
560  HP_Quadrature1D hp_quad;
561  // use symmetry and compute just half of the points
562  for (int i = 0 ; i <= (np-1)/2 ; ++i)
563  {
564  hp_quad.ComputeGaussLobattoPoint(np, i);
565  ir->IntPoint(i).x = hp_quad.GetPoint();
566  ir->IntPoint(np-1-i).x = hp_quad.GetSymmPoint();
567  ir->IntPoint(i).weight =
568  ir->IntPoint(np-1-i).weight = hp_quad.GetWeight();
569  }
570 
571 #endif // MFEM_USE_MPFR
572 
573  }
574 }
575 
576 void QuadratureFunctions1D::OpenUniform(const int np, IntegrationRule* ir)
577 {
578  ir->SetSize(np);
579 
580  // The Newton-Cotes quadrature is based on weights that integrate exactly the
581  // interpolatory polynomial through the equally spaced quadrature points.
582  for (int i = 0; i < np ; ++i)
583  {
584  ir->IntPoint(i).x = double(i+1) / double(np + 1);
585  }
586 
587  CalculateUniformWeights(ir, Quadrature1D::OpenUniform);
588 }
589 
590 void QuadratureFunctions1D::ClosedUniform(const int np,
591  IntegrationRule* ir)
592 {
593  ir->SetSize(np);
594  if ( np == 1 ) // allow this case as "closed"
595  {
596  ir->IntPoint(0).Set1w(0.5, 1.0);
597  return;
598  }
599 
600  for (int i = 0; i < np ; ++i)
601  {
602  ir->IntPoint(i).x = double(i) / (np-1);
603  }
604 
605  CalculateUniformWeights(ir, Quadrature1D::ClosedUniform);
606 }
607 
608 void QuadratureFunctions1D::OpenHalfUniform(const int np, IntegrationRule* ir)
609 {
610  ir->SetSize(np);
611 
612  // Open half points: the centers of np uniform intervals
613  for (int i = 0; i < np ; ++i)
614  {
615  ir->IntPoint(i).x = double(2*i+1) / (2*np);
616  }
617 
618  CalculateUniformWeights(ir, Quadrature1D::OpenHalfUniform);
619 }
620 
621 void QuadratureFunctions1D::ClosedGL(const int np, IntegrationRule* ir)
622 {
623  ir->SetSize(np);
624  ir->IntPoint(0).x = 0.0;
625  ir->IntPoint(np-1).x = 1.0;
626 
627  if ( np > 2 )
628  {
629  IntegrationRule gl_ir;
630  GaussLegendre(np-1, &gl_ir);
631 
632  for (int i = 1; i < np-1; ++i)
633  {
634  ir->IntPoint(i).x = (gl_ir.IntPoint(i-1).x + gl_ir.IntPoint(i).x)/2;
635  }
636  }
637 
638  CalculateUniformWeights(ir, Quadrature1D::ClosedGL);
639 }
640 
641 void QuadratureFunctions1D::GivePolyPoints(const int np, double *pts,
642  const int type)
643 {
644  IntegrationRule ir(np);
645 
646  switch (type)
647  {
648  case Quadrature1D::GaussLegendre:
649  {
650  GaussLegendre(np,&ir);
651  break;
652  }
653  case Quadrature1D::GaussLobatto:
654  {
655  GaussLobatto(np, &ir);
656  break;
657  }
658  case Quadrature1D::OpenUniform:
659  {
660  OpenUniform(np,&ir);
661  break;
662  }
663  case Quadrature1D::ClosedUniform:
664  {
665  ClosedUniform(np,&ir);
666  break;
667  }
668  case Quadrature1D::OpenHalfUniform:
669  {
670  OpenHalfUniform(np, &ir);
671  break;
672  }
673  case Quadrature1D::ClosedGL:
674  {
675  ClosedGL(np, &ir);
676  break;
677  }
678  default:
679  {
680  MFEM_ABORT("Asking for an unknown type of 1D Quadrature points, "
681  "type = " << type);
682  }
683  }
684 
685  for (int i = 0 ; i < np ; ++i)
686  {
687  pts[i] = ir.IntPoint(i).x;
688  }
689 }
690 
691 void QuadratureFunctions1D::CalculateUniformWeights(IntegrationRule *ir,
692  const int type)
693 {
694  /* The Lagrange polynomials are:
695  p_i = \prod_{j \neq i} {\frac{x - x_j }{x_i - x_j}}
696 
697  The weight associated with each abscissa is the integral of p_i over
698  [0,1]. To calculate the integral of p_i, we use a Gauss-Legendre
699  quadrature rule. This approach does not suffer from bad round-off/
700  cancellation errors for large number of points.
701  */
702  const int n = ir->Size();
703  switch (n)
704  {
705  case 1:
706  ir->IntPoint(0).weight = 1.;
707  return;
708  case 2:
709  ir->IntPoint(0).weight = .5;
710  ir->IntPoint(1).weight = .5;
711  return;
712  }
713 
714 #ifndef MFEM_USE_MPFR
715 
716  // This algorithm should work for any set of points, not just uniform
717  const IntegrationRule &glob_ir = IntRules.Get(Geometry::SEGMENT, n-1);
718  const int m = glob_ir.GetNPoints();
719  Vector xv(n);
720  for (int j = 0; j < n; j++)
721  {
722  xv(j) = ir->IntPoint(j).x;
723  }
724  Poly_1D::Basis basis(n-1, xv.GetData()); // nodal basis, with nodes at 'xv'
725  Vector w(n);
726  // Integrate all nodal basis functions using 'glob_ir':
727  w = 0.0;
728  for (int i = 0; i < m; i++)
729  {
730  const IntegrationPoint &ip = glob_ir.IntPoint(i);
731  basis.Eval(ip.x, xv);
732  w.Add(ip.weight, xv); // w += ip.weight * xv
733  }
734  for (int j = 0; j < n; j++)
735  {
736  ir->IntPoint(j).weight = w(j);
737  }
738 
739 #else // MFEM_USE_MPFR is defined
740 
741  static const mpfr_rnd_t rnd = HP_Quadrature1D::rnd;
742  HP_Quadrature1D hp_quad;
743  mpfr_t l, lk, w0, wi, tmp, *weights;
744  mpfr_inits2(hp_quad.default_prec, l, lk, w0, wi, tmp, (mpfr_ptr) 0);
745  weights = new mpfr_t[n];
746  for (int i = 0; i < n; i++)
747  {
748  mpfr_init2(weights[i], hp_quad.default_prec);
749  mpfr_set_si(weights[i], 0, rnd);
750  }
751  hp_quad.SetRelTol(-48); // rtol = 2^(-48) ~ 3.5e-15
752  const int p = n-1;
753  const int m = p/2+1; // number of points for Gauss-Legendre quadrature
754  int hinv = 0, ihoffset = 0; // x_i = (i+ihoffset/2)/hinv
755  switch (type)
756  {
757  case Quadrature1D::ClosedUniform:
758  // x_i = i/p, i=0,...,p
759  hinv = p;
760  ihoffset = 0;
761  break;
762  case Quadrature1D::OpenUniform:
763  // x_i = (i+1)/(p+2), i=0,...,p
764  hinv = p+2;
765  ihoffset = 2;
766  break;
767  case Quadrature1D::OpenHalfUniform:
768  // x_i = (i+1/2)/(p+1), i=0,...,p
769  hinv = p+1;
770  ihoffset = 1;
771  break;
772  default:
773  MFEM_ABORT("invalid Quadrature1D type: " << type);
774  }
775  // set w0 = (-1)^p*(p!)/(hinv^p)
776  mpfr_fac_ui(w0, p, rnd);
777  mpfr_ui_pow_ui(tmp, hinv, p, rnd);
778  mpfr_div(w0, w0, tmp, rnd);
779  if (p%2) { mpfr_neg(w0, w0, rnd); }
780 
781  for (int j = 0; j < m; j++)
782  {
783  hp_quad.ComputeGaussLegendrePoint(m, j);
784 
785  // Compute l = \prod_{i=0}^p (x-x_i) and lk = l/(x-x_k), where
786  // x = hp_quad.GetHPPoint(), x_i = (i+ihoffset/2)/hinv, and x_k is the
787  // node closest to x, i.e. k = min(max(round(x*hinv-ihoffset/2),0),p)
788  mpfr_mul_si(tmp, hp_quad.GetHPPoint(), hinv, rnd);
789  mpfr_sub_d(tmp, tmp, 0.5*ihoffset, rnd);
790  mpfr_round(tmp, tmp);
791  int k = min(max((int)mpfr_get_si(tmp, rnd), 0), p);
792  mpfr_set_si(lk, 1, rnd);
793  for (int i = 0; i <= p; i++)
794  {
795  mpfr_set_si(tmp, 2*i+ihoffset, rnd);
796  mpfr_div_si(tmp, tmp, 2*hinv, rnd);
797  mpfr_sub(tmp, hp_quad.GetHPPoint(), tmp, rnd);
798  if (i != k)
799  {
800  mpfr_mul(lk, lk, tmp, rnd);
801  }
802  else
803  {
804  mpfr_set(l, tmp, rnd);
805  }
806  }
807  mpfr_mul(l, l, lk, rnd);
808  mpfr_set(wi, w0, rnd);
809  for (int i = 0; true; i++)
810  {
811  if (i != k)
812  {
813  // tmp = l/(wi*(x - x_i))
814  mpfr_set_si(tmp, 2*i+ihoffset, rnd);
815  mpfr_div_si(tmp, tmp, 2*hinv, rnd);
816  mpfr_sub(tmp, hp_quad.GetHPPoint(), tmp, rnd);
817  mpfr_mul(tmp, tmp, wi, rnd);
818  mpfr_div(tmp, l, tmp, rnd);
819  }
820  else
821  {
822  // tmp = lk/wi
823  mpfr_div(tmp, lk, wi, rnd);
824  }
825  // weights[i] += hp_quad.weight*tmp
826  mpfr_mul(tmp, tmp, hp_quad.GetHPWeight(), rnd);
827  mpfr_add(weights[i], weights[i], tmp, rnd);
828 
829  if (i == p) { break; }
830 
831  // update wi *= (i+1)/(i-p)
832  mpfr_mul_si(wi, wi, i+1, rnd);
833  mpfr_div_si(wi, wi, i-p, rnd);
834  }
835  }
836  for (int i = 0; i < n; i++)
837  {
838  ir->IntPoint(i).weight = mpfr_get_d(weights[i], rnd);
839  mpfr_clear(weights[i]);
840  }
841  delete [] weights;
842  mpfr_clears(l, lk, w0, wi, tmp, (mpfr_ptr) 0);
843 
844 #endif // MFEM_USE_MPFR
845 
846 }
847 
848 
849 int Quadrature1D::CheckClosed(int type)
850 {
851  switch (type)
852  {
853  case GaussLobatto:
854  case ClosedUniform:
855  return type;
856  default:
857  return Invalid;
858  }
859 }
860 
861 int Quadrature1D::CheckOpen(int type)
862 {
863  switch (type)
864  {
865  case GaussLegendre:
866  case GaussLobatto:
867  case OpenUniform:
868  case ClosedUniform:
869  case OpenHalfUniform:
870  return type; // all types can work as open
871  default:
872  return Invalid;
873  }
874 }
875 
876 
877 IntegrationRules IntRules(0, Quadrature1D::GaussLegendre);
878 
879 IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre);
880 
881 IntegrationRules::IntegrationRules(int Ref, int _type):
882  quad_type(_type)
883 {
884  refined = Ref;
885 
886  if (refined < 0) { own_rules = 0; return; }
887 
888  own_rules = 1;
889 
890  const MemoryType h_mt = MemoryType::HOST;
891  PointIntRules.SetSize(2, h_mt);
892  PointIntRules = NULL;
893 
894  SegmentIntRules.SetSize(32, h_mt);
895  SegmentIntRules = NULL;
896 
897  // TriangleIntegrationRule() assumes that this size is >= 26
898  TriangleIntRules.SetSize(32, h_mt);
899  TriangleIntRules = NULL;
900 
901  SquareIntRules.SetSize(32, h_mt);
902  SquareIntRules = NULL;
903 
904  // TetrahedronIntegrationRule() assumes that this size is >= 10
905  TetrahedronIntRules.SetSize(32, h_mt);
906  TetrahedronIntRules = NULL;
907 
908  PrismIntRules.SetSize(32, h_mt);
909  PrismIntRules = NULL;
910 
911  CubeIntRules.SetSize(32, h_mt);
912  CubeIntRules = NULL;
913 }
914 
915 const IntegrationRule &IntegrationRules::Get(int GeomType, int Order)
916 {
917  Array<IntegrationRule *> *ir_array;
918 
919  switch (GeomType)
920  {
921  case Geometry::POINT: ir_array = &PointIntRules; Order = 0; break;
922  case Geometry::SEGMENT: ir_array = &SegmentIntRules; break;
923  case Geometry::TRIANGLE: ir_array = &TriangleIntRules; break;
924  case Geometry::SQUARE: ir_array = &SquareIntRules; break;
925  case Geometry::TETRAHEDRON: ir_array = &TetrahedronIntRules; break;
926  case Geometry::CUBE: ir_array = &CubeIntRules; break;
927  case Geometry::PRISM: ir_array = &PrismIntRules; break;
928  default:
929  mfem_error("IntegrationRules::Get(...) : Unknown geometry type!");
930  ir_array = NULL;
931  }
932 
933  if (Order < 0)
934  {
935  Order = 0;
936  }
937 
938  if (!HaveIntRule(*ir_array, Order))
939  {
940 #ifdef MFEM_USE_LEGACY_OPENMP
941  #pragma omp critical
942 #endif
943  {
944  if (!HaveIntRule(*ir_array, Order))
945  {
946  IntegrationRule *ir = GenerateIntegrationRule(GeomType, Order);
947  int RealOrder = Order;
948  while (RealOrder+1 < ir_array->Size() &&
949  /* */ (*ir_array)[RealOrder+1] == ir)
950  {
951  RealOrder++;
952  }
953  ir->SetOrder(RealOrder);
954  }
955  }
956  }
957 
958  return *(*ir_array)[Order];
959 }
960 
961 void IntegrationRules::Set(int GeomType, int Order, IntegrationRule &IntRule)
962 {
963  Array<IntegrationRule *> *ir_array;
964 
965  switch (GeomType)
966  {
967  case Geometry::POINT: ir_array = &PointIntRules; break;
968  case Geometry::SEGMENT: ir_array = &SegmentIntRules; break;
969  case Geometry::TRIANGLE: ir_array = &TriangleIntRules; break;
970  case Geometry::SQUARE: ir_array = &SquareIntRules; break;
971  case Geometry::TETRAHEDRON: ir_array = &TetrahedronIntRules; break;
972  case Geometry::CUBE: ir_array = &CubeIntRules; break;
973  case Geometry::PRISM: ir_array = &PrismIntRules; break;
974  default:
975  mfem_error("IntegrationRules::Set(...) : Unknown geometry type!");
976  ir_array = NULL;
977  }
978 
979  if (HaveIntRule(*ir_array, Order))
980  {
981  MFEM_ABORT("Overwriting set rules is not supported!");
982  }
983 
984  AllocIntRule(*ir_array, Order);
985 
986  (*ir_array)[Order] = &IntRule;
987 }
988 
989 void IntegrationRules::DeleteIntRuleArray(Array<IntegrationRule *> &ir_array)
990 {
991  int i;
992  IntegrationRule *ir = NULL;
993 
994  // Many of the intrules have multiple contiguous copies in the ir_array
995  // so we have to be careful to not delete them twice.
996  for (i = 0; i < ir_array.Size(); i++)
997  {
998  if (ir_array[i] != NULL && ir_array[i] != ir)
999  {
1000  ir = ir_array[i];
1001  delete ir;
1002  }
1003  }
1004 }
1005 
1007 {
1008  if (!own_rules) { return; }
1009 
1010  DeleteIntRuleArray(PointIntRules);
1011  DeleteIntRuleArray(SegmentIntRules);
1012  DeleteIntRuleArray(TriangleIntRules);
1013  DeleteIntRuleArray(SquareIntRules);
1014  DeleteIntRuleArray(TetrahedronIntRules);
1015  DeleteIntRuleArray(CubeIntRules);
1016  DeleteIntRuleArray(PrismIntRules);
1017 }
1018 
1019 
1020 IntegrationRule *IntegrationRules::GenerateIntegrationRule(int GeomType,
1021  int Order)
1022 {
1023  switch (GeomType)
1024  {
1025  case Geometry::POINT:
1026  return PointIntegrationRule(Order);
1027  case Geometry::SEGMENT:
1028  return SegmentIntegrationRule(Order);
1029  case Geometry::TRIANGLE:
1030  return TriangleIntegrationRule(Order);
1031  case Geometry::SQUARE:
1032  return SquareIntegrationRule(Order);
1033  case Geometry::TETRAHEDRON:
1034  return TetrahedronIntegrationRule(Order);
1035  case Geometry::CUBE:
1036  return CubeIntegrationRule(Order);
1037  case Geometry::PRISM:
1038  return PrismIntegrationRule(Order);
1039  default:
1040  mfem_error("IntegrationRules::Set(...) : Unknown geometry type!");
1041  return NULL;
1042  }
1043 }
1044 
1045 
1046 // Integration rules for a point
1047 IntegrationRule *IntegrationRules::PointIntegrationRule(int Order)
1048 {
1049  if (Order > 1)
1050  {
1051  mfem_error("Point Integration Rule of Order > 1 not defined");
1052  return NULL;
1053  }
1054 
1055  IntegrationRule *ir = new IntegrationRule(1);
1056  ir->IntPoint(0).x = .0;
1057  ir->IntPoint(0).weight = 1.;
1058 
1059  PointIntRules[1] = PointIntRules[0] = ir;
1060 
1061  return ir;
1062 }
1063 
1064 // Integration rules for line segment [0,1]
1065 IntegrationRule *IntegrationRules::SegmentIntegrationRule(int Order)
1066 {
1067  int RealOrder = GetSegmentRealOrder(Order); // RealOrder >= Order
1068  // Order is one of {RealOrder-1,RealOrder}
1069  AllocIntRule(SegmentIntRules, RealOrder);
1070 
1071  IntegrationRule tmp, *ir;
1072  ir = refined ? &tmp : new IntegrationRule;
1073 
1074  int n = 0;
1075  // n is the number of points to achieve the exact integral of a
1076  // degree Order polynomial
1077  switch (quad_type)
1078  {
1080  {
1081  // Gauss-Legendre is exact for 2*n-1
1082  n = Order/2 + 1;
1083  quad_func.GaussLegendre(n, ir);
1084  break;
1085  }
1087  {
1088  // Gauss-Lobatto is exact for 2*n-3
1089  n = Order/2 + 2;
1090  quad_func.GaussLobatto(n, ir);
1091  break;
1092  }
1094  {
1095  // Open Newton Cotes is exact for n-(n+1)%2 = n-1+n%2
1096  n = Order | 1; // n is always odd
1097  quad_func.OpenUniform(n, ir);
1098  break;
1099  }
1101  {
1102  // Closed Newton Cotes is exact for n-(n+1)%2 = n-1+n%2
1103  n = Order | 1; // n is always odd
1104  quad_func.ClosedUniform(n, ir);
1105  break;
1106  }
1108  {
1109  // Open half Newton Cotes is exact for n-(n+1)%2 = n-1+n%2
1110  n = Order | 1; // n is always odd
1111  quad_func.OpenHalfUniform(n, ir);
1112  break;
1113  }
1114  default:
1115  {
1116  MFEM_ABORT("unknown Quadrature1D type: " << quad_type);
1117  }
1118  }
1119  if (refined)
1120  {
1121  // Effectively passing memory management to SegmentIntegrationRules
1122  ir = new IntegrationRule(2*n);
1123  for (int j = 0; j < n; j++)
1124  {
1125  ir->IntPoint(j).x = tmp.IntPoint(j).x/2.0;
1126  ir->IntPoint(j).weight = tmp.IntPoint(j).weight/2.0;
1127  ir->IntPoint(j+n).x = 0.5 + tmp.IntPoint(j).x/2.0;
1128  ir->IntPoint(j+n).weight = tmp.IntPoint(j).weight/2.0;
1129  }
1130  }
1131  SegmentIntRules[RealOrder-1] = SegmentIntRules[RealOrder] = ir;
1132  return ir;
1133 }
1134 
1135 // Integration rules for reference triangle {[0,0],[1,0],[0,1]}
1136 IntegrationRule *IntegrationRules::TriangleIntegrationRule(int Order)
1137 {
1138  IntegrationRule *ir = NULL;
1139  // Note: Set TriangleIntRules[*] to ir only *after* ir is fully constructed.
1140  // This is needed in multithreaded environment.
1141 
1142  // assuming that orders <= 25 are pre-allocated
1143  switch (Order)
1144  {
1145  case 0: // 1 point - 0 degree
1146  case 1:
1147  ir = new IntegrationRule(1);
1148  ir->AddTriMidPoint(0, 0.5);
1149  TriangleIntRules[0] = TriangleIntRules[1] = ir;
1150  return ir;
1151 
1152  case 2: // 3 point - 2 degree
1153  ir = new IntegrationRule(3);
1154  ir->AddTriPoints3(0, 1./6., 1./6.);
1155  TriangleIntRules[2] = ir;
1156  // interior points
1157  return ir;
1158 
1159  case 3: // 4 point - 3 degree (has one negative weight)
1160  ir = new IntegrationRule(4);
1161  ir->AddTriMidPoint(0, -0.28125); // -9./32.
1162  ir->AddTriPoints3(1, 0.2, 25./96.);
1163  TriangleIntRules[3] = ir;
1164  return ir;
1165 
1166  case 4: // 6 point - 4 degree
1167  ir = new IntegrationRule(6);
1168  ir->AddTriPoints3(0, 0.091576213509770743460, 0.054975871827660933819);
1169  ir->AddTriPoints3(3, 0.44594849091596488632, 0.11169079483900573285);
1170  TriangleIntRules[4] = ir;
1171  return ir;
1172 
1173  case 5: // 7 point - 5 degree
1174  ir = new IntegrationRule(7);
1175  ir->AddTriMidPoint(0, 0.1125);
1176  ir->AddTriPoints3(1, 0.10128650732345633880, 0.062969590272413576298);
1177  ir->AddTriPoints3(4, 0.47014206410511508977, 0.066197076394253090369);
1178  TriangleIntRules[5] = ir;
1179  return ir;
1180 
1181  case 6: // 12 point - 6 degree
1182  ir = new IntegrationRule(12);
1183  ir->AddTriPoints3(0, 0.063089014491502228340, 0.025422453185103408460);
1184  ir->AddTriPoints3(3, 0.24928674517091042129, 0.058393137863189683013);
1185  ir->AddTriPoints6(6, 0.053145049844816947353, 0.31035245103378440542,
1186  0.041425537809186787597);
1187  TriangleIntRules[6] = ir;
1188  return ir;
1189 
1190  case 7: // 12 point - degree 7
1191  ir = new IntegrationRule(12);
1192  ir->AddTriPoints3R(0, 0.062382265094402118174, 0.067517867073916085443,
1193  0.026517028157436251429);
1194  ir->AddTriPoints3R(3, 0.055225456656926611737, 0.32150249385198182267,
1195  0.043881408714446055037);
1196  // slightly better with explicit 3rd area coordinate
1197  ir->AddTriPoints3R(6, 0.034324302945097146470, 0.66094919618673565761,
1198  0.30472650086816719592, 0.028775042784981585738);
1199  ir->AddTriPoints3R(9, 0.51584233435359177926, 0.27771616697639178257,
1200  0.20644149867001643817, 0.067493187009802774463);
1201  TriangleIntRules[7] = ir;
1202  return ir;
1203 
1204  case 8: // 16 point - 8 degree
1205  ir = new IntegrationRule(16);
1206  ir->AddTriMidPoint(0, 0.0721578038388935841255455552445323);
1207  ir->AddTriPoints3(1, 0.170569307751760206622293501491464,
1208  0.0516086852673591251408957751460645);
1209  ir->AddTriPoints3(4, 0.0505472283170309754584235505965989,
1210  0.0162292488115990401554629641708902);
1211  ir->AddTriPoints3(7, 0.459292588292723156028815514494169,
1212  0.0475458171336423123969480521942921);
1213  ir->AddTriPoints6(10, 0.008394777409957605337213834539296,
1214  0.263112829634638113421785786284643,
1215  0.0136151570872174971324223450369544);
1216  TriangleIntRules[8] = ir;
1217  return ir;
1218 
1219  case 9: // 19 point - 9 degree
1220  ir = new IntegrationRule(19);
1221  ir->AddTriMidPoint(0, 0.0485678981413994169096209912536443);
1222  ir->AddTriPoints3b(1, 0.020634961602524744433,
1223  0.0156673501135695352684274156436046);
1224  ir->AddTriPoints3b(4, 0.12582081701412672546,
1225  0.0389137705023871396583696781497019);
1226  ir->AddTriPoints3(7, 0.188203535619032730240961280467335,
1227  0.0398238694636051265164458871320226);
1228  ir->AddTriPoints3(10, 0.0447295133944527098651065899662763,
1229  0.0127888378293490156308393992794999);
1230  ir->AddTriPoints6(13, 0.0368384120547362836348175987833851,
1231  0.2219629891607656956751025276931919,
1232  0.0216417696886446886446886446886446);
1233  TriangleIntRules[9] = ir;
1234  return ir;
1235 
1236  case 10: // 25 point - 10 degree
1237  ir = new IntegrationRule(25);
1238  ir->AddTriMidPoint(0, 0.0454089951913767900476432975500142);
1239  ir->AddTriPoints3b(1, 0.028844733232685245264984935583748,
1240  0.0183629788782333523585030359456832);
1241  ir->AddTriPoints3(4, 0.109481575485037054795458631340522,
1242  0.0226605297177639673913028223692986);
1243  ir->AddTriPoints6(7, 0.141707219414879954756683250476361,
1244  0.307939838764120950165155022930631,
1245  0.0363789584227100543021575883096803);
1246  ir->AddTriPoints6(13, 0.025003534762686386073988481007746,
1247  0.246672560639902693917276465411176,
1248  0.0141636212655287424183685307910495);
1249  ir->AddTriPoints6(19, 0.0095408154002994575801528096228873,
1250  0.0668032510122002657735402127620247,
1251  4.71083348186641172996373548344341E-03);
1252  TriangleIntRules[10] = ir;
1253  return ir;
1254 
1255  case 11: // 28 point -- 11 degree
1256  ir = new IntegrationRule(28);
1257  ir->AddTriPoints6(0, 0.0,
1258  0.141129718717363295960826061941652,
1259  3.68119189165027713212944752369032E-03);
1260  ir->AddTriMidPoint(6, 0.0439886505811161193990465846607278);
1261  ir->AddTriPoints3(7, 0.0259891409282873952600324854988407,
1262  4.37215577686801152475821439991262E-03);
1263  ir->AddTriPoints3(10, 0.0942875026479224956305697762754049,
1264  0.0190407859969674687575121697178070);
1265  ir->AddTriPoints3b(13, 0.010726449965572372516734795387128,
1266  9.42772402806564602923839129555767E-03);
1267  ir->AddTriPoints3(16, 0.207343382614511333452934024112966,
1268  0.0360798487723697630620149942932315);
1269  ir->AddTriPoints3b(19, 0.122184388599015809877869236727746,
1270  0.0346645693527679499208828254519072);
1271  ir->AddTriPoints6(22, 0.0448416775891304433090523914688007,
1272  0.2772206675282791551488214673424523,
1273  0.0205281577146442833208261574536469);
1274  TriangleIntRules[11] = ir;
1275  return ir;
1276 
1277  case 12: // 33 point - 12 degree
1278  ir = new IntegrationRule(33);
1279  ir->AddTriPoints3b(0, 2.35652204523900E-02, 1.28655332202275E-02);
1280  ir->AddTriPoints3b(3, 1.20551215411079E-01, 2.18462722690190E-02);
1281  ir->AddTriPoints3(6, 2.71210385012116E-01, 3.14291121089425E-02);
1282  ir->AddTriPoints3(9, 1.27576145541586E-01, 1.73980564653545E-02);
1283  ir->AddTriPoints3(12, 2.13173504532100E-02, 3.08313052577950E-03);
1284  ir->AddTriPoints6(15, 1.15343494534698E-01, 2.75713269685514E-01,
1285  2.01857788831905E-02);
1286  ir->AddTriPoints6(21, 2.28383322222570E-02, 2.81325580989940E-01,
1287  1.11783866011515E-02);
1288  ir->AddTriPoints6(27, 2.57340505483300E-02, 1.16251915907597E-01,
1289  8.65811555432950E-03);
1290  TriangleIntRules[12] = ir;
1291  return ir;
1292 
1293  case 13: // 37 point - 13 degree
1294  ir = new IntegrationRule(37);
1295  ir->AddTriPoints3b(0, 0.0,
1296  2.67845189554543044455908674650066E-03);
1297  ir->AddTriMidPoint(3, 0.0293480398063595158995969648597808);
1298  ir->AddTriPoints3(4, 0.0246071886432302181878499494124643,
1299  3.92538414805004016372590903990464E-03);
1300  ir->AddTriPoints3b(7, 0.159382493797610632566158925635800,
1301  0.0253344765879434817105476355306468);
1302  ir->AddTriPoints3(10, 0.227900255506160619646298948153592,
1303  0.0250401630452545330803738542916538);
1304  ir->AddTriPoints3(13, 0.116213058883517905247155321839271,
1305  0.0158235572961491595176634480481793);
1306  ir->AddTriPoints3b(16, 0.046794039901841694097491569577008,
1307  0.0157462815379843978450278590138683);
1308  ir->AddTriPoints6(19, 0.0227978945382486125477207592747430,
1309  0.1254265183163409177176192369310890,
1310  7.90126610763037567956187298486575E-03);
1311  ir->AddTriPoints6(25, 0.0162757709910885409437036075960413,
1312  0.2909269114422506044621801030055257,
1313  7.99081889046420266145965132482933E-03);
1314  ir->AddTriPoints6(31, 0.0897330604516053590796290561145196,
1315  0.2723110556841851025078181617634414,
1316  0.0182757511120486476280967518782978);
1317  TriangleIntRules[13] = ir;
1318  return ir;
1319 
1320  case 14: // 42 point - 14 degree
1321  ir = new IntegrationRule(42);
1322  ir->AddTriPoints3b(0, 2.20721792756430E-02, 1.09417906847145E-02);
1323  ir->AddTriPoints3b(3, 1.64710561319092E-01, 1.63941767720625E-02);
1324  ir->AddTriPoints3(6, 2.73477528308839E-01, 2.58870522536460E-02);
1325  ir->AddTriPoints3(9, 1.77205532412543E-01, 2.10812943684965E-02);
1326  ir->AddTriPoints3(12, 6.17998830908730E-02, 7.21684983488850E-03);
1327  ir->AddTriPoints3(15, 1.93909612487010E-02, 2.46170180120000E-03);
1328  ir->AddTriPoints6(18, 5.71247574036480E-02, 1.72266687821356E-01,
1329  1.23328766062820E-02);
1330  ir->AddTriPoints6(24, 9.29162493569720E-02, 3.36861459796345E-01,
1331  1.92857553935305E-02);
1332  ir->AddTriPoints6(30, 1.46469500556540E-02, 2.98372882136258E-01,
1333  7.21815405676700E-03);
1334  ir->AddTriPoints6(36, 1.26833093287200E-03, 1.18974497696957E-01,
1335  2.50511441925050E-03);
1336  TriangleIntRules[14] = ir;
1337  return ir;
1338 
1339  case 15: // 54 point - 15 degree
1340  ir = new IntegrationRule(54);
1341  ir->AddTriPoints3b(0, 0.0834384072617499333, 0.016330909424402645);
1342  ir->AddTriPoints3b(3, 0.192779070841738867, 0.01370640901568218);
1343  ir->AddTriPoints3(6, 0.293197167913025367, 0.01325501829935165);
1344  ir->AddTriPoints3(9, 0.146467786942772933, 0.014607981068243055);
1345  ir->AddTriPoints3(12, 0.0563628676656034333, 0.005292304033121995);
1346  ir->AddTriPoints3(15, 0.0165751268583703333, 0.0018073215320460175);
1347  ir->AddTriPoints6(18, 0.0099122033092248, 0.239534554154794445,
1348  0.004263874050854718);
1349  ir->AddTriPoints6(24, 0.015803770630228, 0.404878807318339958,
1350  0.006958088258345965);
1351  ir->AddTriPoints6(30, 0.00514360881697066667, 0.0950021131130448885,
1352  0.0021459664703674175);
1353  ir->AddTriPoints6(36, 0.0489223257529888, 0.149753107322273969,
1354  0.008117664640887445);
1355  ir->AddTriPoints6(42, 0.0687687486325192, 0.286919612441334979,
1356  0.012803670460631195);
1357  ir->AddTriPoints6(48, 0.1684044181246992, 0.281835668099084562,
1358  0.016544097765822835);
1359  TriangleIntRules[15] = ir;
1360  return ir;
1361 
1362  case 16: // 61 point - 17 degree (used for 16 as well)
1363  case 17:
1364  ir = new IntegrationRule(61);
1365  ir->AddTriMidPoint(0, 1.67185996454015E-02);
1366  ir->AddTriPoints3b(1, 5.65891888645200E-03, 2.54670772025350E-03);
1367  ir->AddTriPoints3b(4, 3.56473547507510E-02, 7.33543226381900E-03);
1368  ir->AddTriPoints3b(7, 9.95200619584370E-02, 1.21754391768360E-02);
1369  ir->AddTriPoints3b(10, 1.99467521245206E-01, 1.55537754344845E-02);
1370  ir->AddTriPoints3 (13, 2.52141267970953E-01, 1.56285556093100E-02);
1371  ir->AddTriPoints3 (16, 1.62047004658461E-01, 1.24078271698325E-02);
1372  ir->AddTriPoints3 (19, 7.58758822607460E-02, 7.02803653527850E-03);
1373  ir->AddTriPoints3 (22, 1.56547269678220E-02, 1.59733808688950E-03);
1374  ir->AddTriPoints6 (25, 1.01869288269190E-02, 3.34319867363658E-01,
1375  4.05982765949650E-03);
1376  ir->AddTriPoints6 (31, 1.35440871671036E-01, 2.92221537796944E-01,
1377  1.34028711415815E-02);
1378  ir->AddTriPoints6 (37, 5.44239242905830E-02, 3.19574885423190E-01,
1379  9.22999660541100E-03);
1380  ir->AddTriPoints6 (43, 1.28685608336370E-02, 1.90704224192292E-01,
1381  4.23843426716400E-03);
1382  ir->AddTriPoints6 (49, 6.71657824135240E-02, 1.80483211648746E-01,
1383  9.14639838501250E-03);
1384  ir->AddTriPoints6 (55, 1.46631822248280E-02, 8.07113136795640E-02,
1385  3.33281600208250E-03);
1386  TriangleIntRules[16] = TriangleIntRules[17] = ir;
1387  return ir;
1388 
1389  case 18: // 73 point - 19 degree (used for 18 as well)
1390  case 19:
1391  ir = new IntegrationRule(73);
1392  ir->AddTriMidPoint(0, 0.0164531656944595);
1393  ir->AddTriPoints3b(1, 0.020780025853987, 0.005165365945636);
1394  ir->AddTriPoints3b(4, 0.090926214604215, 0.011193623631508);
1395  ir->AddTriPoints3b(7, 0.197166638701138, 0.015133062934734);
1396  ir->AddTriPoints3 (10, 0.255551654403098, 0.015245483901099);
1397  ir->AddTriPoints3 (13, 0.17707794215213, 0.0120796063708205);
1398  ir->AddTriPoints3 (16, 0.110061053227952, 0.0080254017934005);
1399  ir->AddTriPoints3 (19, 0.05552862425184, 0.004042290130892);
1400  ir->AddTriPoints3 (22, 0.012621863777229, 0.0010396810137425);
1401  ir->AddTriPoints6 (25, 0.003611417848412, 0.395754787356943,
1402  0.0019424384524905);
1403  ir->AddTriPoints6 (31, 0.13446675453078, 0.307929983880436,
1404  0.012787080306011);
1405  ir->AddTriPoints6 (37, 0.014446025776115, 0.26456694840652,
1406  0.004440451786669);
1407  ir->AddTriPoints6 (43, 0.046933578838178, 0.358539352205951,
1408  0.0080622733808655);
1409  ir->AddTriPoints6 (49, 0.002861120350567, 0.157807405968595,
1410  0.0012459709087455);
1411  ir->AddTriPoints6 (55, 0.075050596975911, 0.223861424097916,
1412  0.0091214200594755);
1413  ir->AddTriPoints6 (61, 0.03464707481676, 0.142421601113383,
1414  0.0051292818680995);
1415  ir->AddTriPoints6 (67, 0.065494628082938, 0.010161119296278,
1416  0.001899964427651);
1417  TriangleIntRules[18] = TriangleIntRules[19] = ir;
1418  return ir;
1419 
1420  case 20: // 85 point - 20 degree
1421  ir = new IntegrationRule(85);
1422  ir->AddTriMidPoint(0, 0.01380521349884976);
1423  ir->AddTriPoints3b(1, 0.001500649324429, 0.00088951477366337);
1424  ir->AddTriPoints3b(4, 0.0941397519389508667, 0.010056199056980585);
1425  ir->AddTriPoints3b(7, 0.2044721240895264, 0.013408923629665785);
1426  ir->AddTriPoints3(10, 0.264500202532787333, 0.012261566900751005);
1427  ir->AddTriPoints3(13, 0.211018964092076767, 0.008197289205347695);
1428  ir->AddTriPoints3(16, 0.107735607171271333, 0.0073979536993248);
1429  ir->AddTriPoints3(19, 0.0390690878378026667, 0.0022896411388521255);
1430  ir->AddTriPoints3(22, 0.0111743797293296333, 0.0008259132577881085);
1431  ir->AddTriPoints6(25, 0.00534961818733726667, 0.0635496659083522206,
1432  0.001174585454287792);
1433  ir->AddTriPoints6(31, 0.00795481706619893333, 0.157106918940706982,
1434  0.0022329628770908965);
1435  ir->AddTriPoints6(37, 0.0104223982812638, 0.395642114364374018,
1436  0.003049783403953986);
1437  ir->AddTriPoints6(43, 0.0109644147961233333, 0.273167570712910522,
1438  0.0034455406635941015);
1439  ir->AddTriPoints6(49, 0.0385667120854623333, 0.101785382485017108,
1440  0.0039987375362390815);
1441  ir->AddTriPoints6(55, 0.0355805078172182, 0.446658549176413815,
1442  0.003693067142668012);
1443  ir->AddTriPoints6(61, 0.0496708163627641333, 0.199010794149503095,
1444  0.00639966593932413);
1445  ir->AddTriPoints6(67, 0.0585197250843317333, 0.3242611836922827,
1446  0.008629035587848275);
1447  ir->AddTriPoints6(73, 0.121497787004394267, 0.208531363210132855,
1448  0.009336472951467735);
1449  ir->AddTriPoints6(79, 0.140710844943938733, 0.323170566536257485,
1450  0.01140911202919763);
1451  TriangleIntRules[20] = ir;
1452  return ir;
1453 
1454  case 21: // 126 point - 25 degree (used also for degrees from 21 to 24)
1455  case 22:
1456  case 23:
1457  case 24:
1458  case 25:
1459  ir = new IntegrationRule(126);
1460  ir->AddTriPoints3b(0, 0.0279464830731742, 0.0040027909400102085);
1461  ir->AddTriPoints3b(3, 0.131178601327651467, 0.00797353841619525);
1462  ir->AddTriPoints3b(6, 0.220221729512072267, 0.006554570615397765);
1463  ir->AddTriPoints3 (9, 0.298443234019804467, 0.00979150048281781);
1464  ir->AddTriPoints3(12, 0.2340441723373718, 0.008235442720768635);
1465  ir->AddTriPoints3(15, 0.151468334609017567, 0.00427363953704605);
1466  ir->AddTriPoints3(18, 0.112733893545993667, 0.004080942928613246);
1467  ir->AddTriPoints3(21, 0.0777156920915263, 0.0030605732699918895);
1468  ir->AddTriPoints3(24, 0.034893093614297, 0.0014542491324683325);
1469  ir->AddTriPoints3(27, 0.00725818462093236667, 0.00034613762283099815);
1470  ir->AddTriPoints6(30, 0.0012923527044422, 0.227214452153364077,
1471  0.0006241445996386985);
1472  ir->AddTriPoints6(36, 0.0053997012721162, 0.435010554853571706,
1473  0.001702376454401511);
1474  ir->AddTriPoints6(42, 0.006384003033975, 0.320309599272204437,
1475  0.0016798271630320255);
1476  ir->AddTriPoints6(48, 0.00502821150199306667, 0.0917503222800051889,
1477  0.000858078269748377);
1478  ir->AddTriPoints6(54, 0.00682675862178186667, 0.0380108358587243835,
1479  0.000740428158357803);
1480  ir->AddTriPoints6(60, 0.0100161996399295333, 0.157425218485311668,
1481  0.0017556563053643425);
1482  ir->AddTriPoints6(66, 0.02575781317339, 0.239889659778533193,
1483  0.003696775074853242);
1484  ir->AddTriPoints6(72, 0.0302278981199158, 0.361943118126060531,
1485  0.003991543738688279);
1486  ir->AddTriPoints6(78, 0.0305049901071620667, 0.0835519609548285602,
1487  0.0021779813065790205);
1488  ir->AddTriPoints6(84, 0.0459565473625693333, 0.148443220732418205,
1489  0.003682528350708916);
1490  ir->AddTriPoints6(90, 0.0674428005402775333, 0.283739708727534955,
1491  0.005481786423209775);
1492  ir->AddTriPoints6(96, 0.0700450914159106, 0.406899375118787573,
1493  0.00587498087177056);
1494  ir->AddTriPoints6(102, 0.0839115246401166, 0.194113987024892542,
1495  0.005007800356899285);
1496  ir->AddTriPoints6(108, 0.120375535677152667, 0.32413434700070316,
1497  0.00665482039381434);
1498  ir->AddTriPoints6(114, 0.148066899157366667, 0.229277483555980969,
1499  0.00707722325261307);
1500  ir->AddTriPoints6(120, 0.191771865867325067, 0.325618122595983752,
1501  0.007440689780584005);
1502  TriangleIntRules[21] =
1503  TriangleIntRules[22] =
1504  TriangleIntRules[23] =
1505  TriangleIntRules[24] =
1506  TriangleIntRules[25] = ir;
1507  return ir;
1508 
1509  default:
1510  // Grundmann-Moller rules
1511  int i = (Order / 2) * 2 + 1; // Get closest odd # >= Order
1512  AllocIntRule(TriangleIntRules, i);
1513  ir = new IntegrationRule;
1514  ir->GrundmannMollerSimplexRule(i/2,2);
1515  TriangleIntRules[i-1] = TriangleIntRules[i] = ir;
1516  return ir;
1517  }
1518 }
1519 
1520 // Integration rules for unit square
1521 IntegrationRule *IntegrationRules::SquareIntegrationRule(int Order)
1522 {
1523  int RealOrder = GetSegmentRealOrder(Order);
1524  // Order is one of {RealOrder-1,RealOrder}
1525  if (!HaveIntRule(SegmentIntRules, RealOrder))
1526  {
1527  SegmentIntegrationRule(RealOrder);
1528  }
1529  AllocIntRule(SquareIntRules, RealOrder); // RealOrder >= Order
1530  SquareIntRules[RealOrder-1] =
1531  SquareIntRules[RealOrder] =
1532  new IntegrationRule(*SegmentIntRules[RealOrder],
1533  *SegmentIntRules[RealOrder]);
1534  return SquareIntRules[Order];
1535 }
1536 
1537 /** Integration rules for reference tetrahedron
1538  {[0,0,0],[1,0,0],[0,1,0],[0,0,1]} */
1539 IntegrationRule *IntegrationRules::TetrahedronIntegrationRule(int Order)
1540 {
1541  IntegrationRule *ir;
1542  // Note: Set TetrahedronIntRules[*] to ir only *after* ir is fully
1543  // constructed. This is needed in multithreaded environment.
1544 
1545  // assuming that orders <= 9 are pre-allocated
1546  switch (Order)
1547  {
1548  case 0: // 1 point - degree 1
1549  case 1:
1550  ir = new IntegrationRule(1);
1551  ir->AddTetMidPoint(0, 1./6.);
1552  TetrahedronIntRules[0] = TetrahedronIntRules[1] = ir;
1553  return ir;
1554 
1555  case 2: // 4 points - degree 2
1556  ir = new IntegrationRule(4);
1557  // ir->AddTetPoints4(0, 0.13819660112501051518, 1./24.);
1558  ir->AddTetPoints4b(0, 0.58541019662496845446, 1./24.);
1559  TetrahedronIntRules[2] = ir;
1560  return ir;
1561 
1562  case 3: // 5 points - degree 3 (negative weight)
1563  ir = new IntegrationRule(5);
1564  ir->AddTetMidPoint(0, -2./15.);
1565  ir->AddTetPoints4b(1, 0.5, 0.075);
1566  TetrahedronIntRules[3] = ir;
1567  return ir;
1568 
1569  case 4: // 11 points - degree 4 (negative weight)
1570  ir = new IntegrationRule(11);
1571  ir->AddTetPoints4(0, 1./14., 343./45000.);
1572  ir->AddTetMidPoint(4, -74./5625.);
1573  ir->AddTetPoints6(5, 0.10059642383320079500, 28./1125.);
1574  TetrahedronIntRules[4] = ir;
1575  return ir;
1576 
1577  case 5: // 14 points - degree 5
1578  ir = new IntegrationRule(14);
1579  ir->AddTetPoints6(0, 0.045503704125649649492,
1580  7.0910034628469110730E-03);
1581  ir->AddTetPoints4(6, 0.092735250310891226402, 0.012248840519393658257);
1582  ir->AddTetPoints4b(10, 0.067342242210098170608,
1583  0.018781320953002641800);
1584  TetrahedronIntRules[5] = ir;
1585  return ir;
1586 
1587  case 6: // 24 points - degree 6
1588  ir = new IntegrationRule(24);
1589  ir->AddTetPoints4(0, 0.21460287125915202929,
1590  6.6537917096945820166E-03);
1591  ir->AddTetPoints4(4, 0.040673958534611353116,
1592  1.6795351758867738247E-03);
1593  ir->AddTetPoints4b(8, 0.032986329573173468968,
1594  9.2261969239424536825E-03);
1595  ir->AddTetPoints12(12, 0.063661001875017525299, 0.26967233145831580803,
1596  8.0357142857142857143E-03);
1597  TetrahedronIntRules[6] = ir;
1598  return ir;
1599 
1600  case 7: // 31 points - degree 7 (negative weight)
1601  ir = new IntegrationRule(31);
1602  ir->AddTetPoints6(0, 0.0, 9.7001763668430335097E-04);
1603  ir->AddTetMidPoint(6, 0.018264223466108820291);
1604  ir->AddTetPoints4(7, 0.078213192330318064374, 0.010599941524413686916);
1605  ir->AddTetPoints4(11, 0.12184321666390517465,
1606  -0.062517740114331851691);
1607  ir->AddTetPoints4b(15, 2.3825066607381275412E-03,
1608  4.8914252630734993858E-03);
1609  ir->AddTetPoints12(19, 0.1, 0.2, 0.027557319223985890653);
1610  TetrahedronIntRules[7] = ir;
1611  return ir;
1612 
1613  case 8: // 43 points - degree 8 (negative weight)
1614  ir = new IntegrationRule(43);
1615  ir->AddTetPoints4(0, 5.7819505051979972532E-03,
1616  1.6983410909288737984E-04);
1617  ir->AddTetPoints4(4, 0.082103588310546723091,
1618  1.9670333131339009876E-03);
1619  ir->AddTetPoints12(8, 0.036607749553197423679, 0.19048604193463345570,
1620  2.1405191411620925965E-03);
1621  ir->AddTetPoints6(20, 0.050532740018894224426,
1622  4.5796838244672818007E-03);
1623  ir->AddTetPoints12(26, 0.22906653611681113960, 0.035639582788534043717,
1624  5.7044858086819185068E-03);
1625  ir->AddTetPoints4(38, 0.20682993161067320408, 0.014250305822866901248);
1626  ir->AddTetMidPoint(42, -0.020500188658639915841);
1627  TetrahedronIntRules[8] = ir;
1628  return ir;
1629 
1630  case 9: // orders 9 and higher -- Grundmann-Moller rules
1631  ir = new IntegrationRule;
1632  ir->GrundmannMollerSimplexRule(4,3);
1633  TetrahedronIntRules[9] = ir;
1634  return ir;
1635 
1636  default: // Grundmann-Moller rules
1637  int i = (Order / 2) * 2 + 1; // Get closest odd # >= Order
1638  AllocIntRule(TetrahedronIntRules, i);
1639  ir = new IntegrationRule;
1640  ir->GrundmannMollerSimplexRule(i/2,3);
1641  TetrahedronIntRules[i-1] = TetrahedronIntRules[i] = ir;
1642  return ir;
1643  }
1644 }
1645 
1646 // Integration rules for reference prism
1647 IntegrationRule *IntegrationRules::PrismIntegrationRule(int Order)
1648 {
1649  const IntegrationRule & irt = Get(Geometry::TRIANGLE, Order);
1650  const IntegrationRule & irs = Get(Geometry::SEGMENT, Order);
1651  int nt = irt.GetNPoints();
1652  int ns = irs.GetNPoints();
1653  AllocIntRule(PrismIntRules, Order);
1654  PrismIntRules[Order] = new IntegrationRule(nt * ns);
1655 
1656  for (int ks=0; ks<ns; ks++)
1657  {
1658  const IntegrationPoint & ips = irs.IntPoint(ks);
1659  for (int kt=0; kt<nt; kt++)
1660  {
1661  int kp = ks * nt + kt;
1662  const IntegrationPoint & ipt = irt.IntPoint(kt);
1663  IntegrationPoint & ipp = PrismIntRules[Order]->IntPoint(kp);
1664  ipp.x = ipt.x;
1665  ipp.y = ipt.y;
1666  ipp.z = ips.x;
1667  ipp.weight = ipt.weight * ips.weight;
1668  }
1669  }
1670  return PrismIntRules[Order];
1671 }
1672 
1673 // Integration rules for reference cube
1674 IntegrationRule *IntegrationRules::CubeIntegrationRule(int Order)
1675 {
1676  int RealOrder = GetSegmentRealOrder(Order);
1677  if (!HaveIntRule(SegmentIntRules, RealOrder))
1678  {
1679  SegmentIntegrationRule(RealOrder);
1680  }
1681  AllocIntRule(CubeIntRules, RealOrder);
1682  CubeIntRules[RealOrder-1] =
1683  CubeIntRules[RealOrder] =
1684  new IntegrationRule(*SegmentIntRules[RealOrder],
1685  *SegmentIntRules[RealOrder],
1686  *SegmentIntRules[RealOrder]);
1687  return CubeIntRules[Order];
1688 }
1689 
1690 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:915
void ClosedUniform(const int np, IntegrationRule *ir)
Definition: intrules.cpp:590
~IntegrationRules()
Destroys an IntegrationRules object.
Definition: intrules.cpp:1006
aka closed Newton-Cotes
Definition: intrules.hpp:296
Container class for integration rules.
Definition: intrules.hpp:309
void OpenUniform(const int np, IntegrationRule *ir)
Definition: intrules.cpp:576
aka open Newton-Cotes
Definition: intrules.hpp:295
void GaussLobatto(const int np, IntegrationRule *ir)
Definition: intrules.cpp:454
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
void OpenHalfUniform(const int np, IntegrationRule *ir)
Definition: intrules.cpp:608
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:28
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:654
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.
Definition: intrules.hpp:242
Class for integration point with weight.
Definition: intrules.hpp:25
Host memory; using new[] and delete[].
aka &quot;open half&quot; Newton-Cotes
Definition: intrules.hpp:297
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition: intrules.hpp:381
void pts(int iphi, int t, double x[])
void Set(int GeomType, int Order, IntegrationRule &IntRule)
Definition: intrules.cpp:961
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:378
void Set1w(const double x1, const double w)
Definition: intrules.hpp:84
void GaussLegendre(const int np, IntegrationRule *ir)
Definition: intrules.cpp:375
double f(const Vector &p)