integral.cc
Go to the documentation of this file.
1 //LIC// ====================================================================
2 //LIC// This file forms part of oomph-lib, the object-oriented,
3 //LIC// multi-physics finite-element library, available
4 //LIC// at http://www.oomph-lib.org.
5 //LIC//
6 //LIC// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 //Static data for the Gaussian integration rules
31 
32 //oomph-lib header
33 #include "integral.h"
34 
35 namespace oomph
36 {
37 
38 double ReflectedIntegral::knot(const unsigned &i, const unsigned &j) const
39 {
40  const unsigned Npts_half = Integral_pt->nweight();
41  if (i < Npts_half)
42  {
43  if (j != Index_s)
44  {
45  return Integral_pt->knot(i,j);
46  }
47  else
48  {
51  }
52  }
53  else
54  {
55  if (j != Index_s)
56  {
57  return Integral_pt->knot(i-Npts_half,j);
58  }
59  else
60  {
61  return fraction_to_local(
63  + local_to_fraction(Integral_pt->knot(i-Npts_half,j))
64  * (1.0-S_reflect_fraction));
65  }
66  }
67 }
68 
69 double ReflectedIntegral::local_to_fraction(const double& s) const
70 {
71  return (s - S_min)/(S_max - S_min);
72 }
73 
74 double ReflectedIntegral::fraction_to_local(const double& fraction) const
75 {
76  return fraction*(S_max - S_min) + S_min;
77 }
78 
79 double ReflectedIntegral::weight(const unsigned &i) const
80 {
81  if (i < Integral_pt->nweight())
82  {
84  }
85  else
86  {
87  return Integral_pt->weight(i - Integral_pt->nweight())
88  *(1.0-S_reflect_fraction);
89  }
90 }
91 
92 //Need to define the static data here
93 
94 //============================================================
95 // Define the positions and weights of the 1D Gauss points
96 //============================================================
97 
98 //------------------------------------------------------------
99 //N=2
100 //------------------------------------------------------------
101 const double Gauss<1,2>::Knot[2][1]=
102 {{-0.577350269189626},{0.577350269189626}};
103 
104 const double Gauss<1,2>::Weight[2]={1.0,1.0};
105 
106 //------------------------------------------------------------
107 //N=3
108 //------------------------------------------------------------
109 const double Gauss<1,3>::Knot[3][1]=
110 {{-0.774596669241483},{0.0},{0.774596669241483}};
111 
112 const double Gauss<1,3>::Weight[3]={(5.0/9.0),(8.0/9.0),(5.0/9.0)};
113 
114 //------------------------------------------------------------
115 //N=4
116 //------------------------------------------------------------
117 const double Gauss<1,4>::Knot[4][1]=
118 {{-0.861136311594053},{-0.339981043584856}, {0.339981043584856},
119  {0.861136311594053}};
120 
121 const double Gauss<1,4>::Weight[4]={0.347854845137448, 0.652145154862546, 0.652145154862546, 0.347854845137448};
122 
123 
124 
125 //============================================================
126 // Define the positions and weights of the 2D Gauss points
127 //============================================================
128 
129 //------------------------------------------------------------
130 //N=2
131 //------------------------------------------------------------
132 const double Gauss<2,2>::Knot[4][2]={
133  {-0.577350269189626,-0.577350269189626},
134  {-0.577350269189626,0.577350269189626},
135  {0.577350269189626,-0.577350269189626},
136  {0.577350269189626,0.577350269189626}};
137 const double Gauss<2,2>::Weight[4] = {1.0,1.0,1.0,1.0};
138 
139 //------------------------------------------------------------
140 //N=3
141 //------------------------------------------------------------
142 const double Gauss<2,3>::Knot[9][2]= {
143  {-0.774596669241483,-0.774596669241483},
144  {-0.774596669241483,0.0},
145  {-0.774596669241483,0.774596662941483},
146  {0.0,-0.774596669241483},
147  {0.0,0.0},
148  {0.0,0.774596662941483},
149  {0.774596662941483,-0.774596669241483},
150  {0.774596662941483,0.0},
151  {0.774596662941483,0.774596662941483}};
152 const double Gauss<2,3>::Weight[9] = {(25.0/81.0),(40.0/81.0),(25.0/81.0),
153  (40.0/81.0),(64.0/81.0),(40.0/81.0),(25.0/81.0),(40.0/81.0),(25.0/81.0)};
154 
155 //------------------------------------------------------------
156 //N=4
157 //------------------------------------------------------------
158 const double Gauss<2,4>::Knot[16][2]= {
159  {-0.861136311594053,-0.861136311594053},
160  {-0.339981043584856,-0.861136311594053},
161  { 0.339981043584856,-0.861136311594053},
162  { 0.861136311594053,-0.861136311594053},
163 
164  {-0.861136311594053,-0.339981043584856},
165  {-0.339981043584856,-0.339981043584856},
166  { 0.339981043584856,-0.339981043584856},
167  { 0.861136311594053,-0.339981043584856},
168 
169  {-0.861136311594053, 0.339981043584856},
170  {-0.339981043584856, 0.339981043584856},
171  { 0.339981043584856, 0.339981043584856},
172  { 0.861136311594053, 0.339981043584856},
173 
174  {-0.861136311594053, 0.861136311594053},
175  {-0.339981043584856, 0.861136311594053},
176  { 0.339981043584856, 0.861136311594053},
177  { 0.861136311594053, 0.861136311594053}};
178 
179 // Quick sanity check: they sum to 4 :)
180 const double Gauss<2,4>::Weight[16] =
181 {0.1210029932855979 , 0.2268518518518480 , 0.2268518518518480, 0.1210029932855979,
182  0.2268518518518480 , 0.4252933030106941 , 0.4252933030106941, 0.2268518518518480,
183  0.2268518518518480 , 0.4252933030106941 , 0.4252933030106941, 0.2268518518518480,
184  0.1210029932855979 , 0.2268518518518480 , 0.2268518518518480, 0.1210029932855979};
185 
186 
187 
188 //============================================================
189 // Define the positions and weights of the 3D Gauss points
190 // (produced with utilities/gauss_weights.cc)
191 //============================================================
192 
193 //------------------------------------------------------------
194 //N=2
195 //------------------------------------------------------------
196 const double Gauss<3,2>::Knot[8][3]=
197 {{-0.57735026918963,-0.57735026918963,-0.57735026918963},{-0.57735026918963,-0.57735026918963,0.57735026918963},{-0.57735026918963,0.57735026918963,-0.57735026918963},{-0.57735026918963,0.57735026918963,0.57735026918963},{0.57735026918963,-0.57735026918963,-0.57735026918963},{0.57735026918963,-0.57735026918963,0.57735026918963},{0.57735026918963,0.57735026918963,-0.57735026918963},{0.57735026918963,0.57735026918963,0.57735026918963}};
198 const double Gauss<3,2>::Weight[8] = {1,1,1,1,1,1,1,1};
199 
200 //------------------------------------------------------------
201 //N=3
202 //------------------------------------------------------------
203 const double Gauss<3,3>::Knot[27][3]=
204 {{-0.77459666924148,-0.77459666924148,-0.77459666924148},{-0.77459666924148,-0.77459666924148,0},{-0.77459666924148,-0.77459666924148,0.77459666924148},{-0.77459666924148,0,-0.77459666924148},{-0.77459666924148,0,0},{-0.77459666924148,0,0.77459666924148},{-0.77459666924148,0.77459666924148,-0.77459666924148},{-0.77459666924148,0.77459666924148,0},{-0.77459666924148,0.77459666924148,0.77459666924148},{0,-0.77459666924148,-0.77459666924148},{0,-0.77459666924148,0},{0,-0.77459666924148,0.77459666924148},{0,0,-0.77459666924148},{0,0,0},{0,0,0.77459666924148},{0,0.77459666924148,-0.77459666924148},{0,0.77459666924148,0},{0,0.77459666924148,0.77459666924148},{0.77459666924148,-0.77459666924148,-0.77459666924148},{0.77459666924148,-0.77459666924148,0},{0.77459666924148,-0.77459666924148,0.77459666924148},{0.77459666924148,0,-0.77459666924148},{0.77459666924148,0,0},{0.77459666924148,0,0.77459666924148},{0.77459666924148,0.77459666924148,-0.77459666924148},{0.77459666924148,0.77459666924148,0},{0.77459666924148,0.77459666924148,0.77459666924148}};
205 
206 
207 const double Gauss<3,3>::Weight[27] = {0.17146776406035,0.27434842249657,0.17146776406035,0.27434842249657,0.43895747599451,0.27434842249657,0.17146776406035,0.27434842249657,0.17146776406035,0.27434842249657,0.43895747599451,0.27434842249657,0.43895747599451,0.70233196159122,0.43895747599451,0.27434842249657,0.43895747599451,0.27434842249657,0.17146776406035,0.27434842249657,0.17146776406035,0.27434842249657,0.43895747599451,0.27434842249657,0.17146776406035,0.27434842249657,0.17146776406035};
208 
209 
210 
211 //------------------------------------------------------------
212 //N=4
213 //------------------------------------------------------------
214 const double Gauss<3,4>::Knot[64][3]=
215 {{-0.86113631159405,-0.86113631159405,-0.86113631159405},{-0.86113631159405,-0.86113631159405,-0.33998104358486},{-0.86113631159405,-0.86113631159405,0.33998104358486},{-0.86113631159405,-0.86113631159405,0.86113631159405},{-0.86113631159405,-0.33998104358486,-0.86113631159405},{-0.86113631159405,-0.33998104358486,-0.33998104358486},{-0.86113631159405,-0.33998104358486,0.33998104358486},{-0.86113631159405,-0.33998104358486,0.86113631159405},{-0.86113631159405,0.33998104358486,-0.86113631159405},{-0.86113631159405,0.33998104358486,-0.33998104358486},{-0.86113631159405,0.33998104358486,0.33998104358486},{-0.86113631159405,0.33998104358486,0.86113631159405},{-0.86113631159405,0.86113631159405,-0.86113631159405},{-0.86113631159405,0.86113631159405,-0.33998104358486},{-0.86113631159405,0.86113631159405,0.33998104358486},{-0.86113631159405,0.86113631159405,0.86113631159405},{-0.33998104358486,-0.86113631159405,-0.86113631159405},{-0.33998104358486,-0.86113631159405,-0.33998104358486},{-0.33998104358486,-0.86113631159405,0.33998104358486},{-0.33998104358486,-0.86113631159405,0.86113631159405},{-0.33998104358486,-0.33998104358486,-0.86113631159405},{-0.33998104358486,-0.33998104358486,-0.33998104358486},{-0.33998104358486,-0.33998104358486,0.33998104358486},{-0.33998104358486,-0.33998104358486,0.86113631159405},{-0.33998104358486,0.33998104358486,-0.86113631159405},{-0.33998104358486,0.33998104358486,-0.33998104358486},{-0.33998104358486,0.33998104358486,0.33998104358486},{-0.33998104358486,0.33998104358486,0.86113631159405},{-0.33998104358486,0.86113631159405,-0.86113631159405},{-0.33998104358486,0.86113631159405,-0.33998104358486},{-0.33998104358486,0.86113631159405,0.33998104358486},{-0.33998104358486,0.86113631159405,0.86113631159405},{0.33998104358486,-0.86113631159405,-0.86113631159405},{0.33998104358486,-0.86113631159405,-0.33998104358486},{0.33998104358486,-0.86113631159405,0.33998104358486},{0.33998104358486,-0.86113631159405,0.86113631159405},{0.33998104358486,-0.33998104358486,-0.86113631159405},{0.33998104358486,-0.33998104358486,-0.33998104358486},{0.33998104358486,-0.33998104358486,0.33998104358486},{0.33998104358486,-0.33998104358486,0.86113631159405},{0.33998104358486,0.33998104358486,-0.86113631159405},{0.33998104358486,0.33998104358486,-0.33998104358486},{0.33998104358486,0.33998104358486,0.33998104358486},{0.33998104358486,0.33998104358486,0.86113631159405},{0.33998104358486,0.86113631159405,-0.86113631159405},{0.33998104358486,0.86113631159405,-0.33998104358486},{0.33998104358486,0.86113631159405,0.33998104358486},{0.33998104358486,0.86113631159405,0.86113631159405},{0.86113631159405,-0.86113631159405,-0.86113631159405},{0.86113631159405,-0.86113631159405,-0.33998104358486},{0.86113631159405,-0.86113631159405,0.33998104358486},{0.86113631159405,-0.86113631159405,0.86113631159405},{0.86113631159405,-0.33998104358486,-0.86113631159405},{0.86113631159405,-0.33998104358486,-0.33998104358486},{0.86113631159405,-0.33998104358486,0.33998104358486},{0.86113631159405,-0.33998104358486,0.86113631159405},{0.86113631159405,0.33998104358486,-0.86113631159405},{0.86113631159405,0.33998104358486,-0.33998104358486},{0.86113631159405,0.33998104358486,0.33998104358486},{0.86113631159405,0.33998104358486,0.86113631159405},{0.86113631159405,0.86113631159405,-0.86113631159405},{0.86113631159405,0.86113631159405,-0.33998104358486},{0.86113631159405,0.86113631159405,0.33998104358486},{0.86113631159405,0.86113631159405,0.86113631159405}};
216 
217 
218 
219 const double Gauss<3,4>::Weight[64] =
220 {0.042091477490529,0.078911515795068,0.078911515795068,0.042091477490529,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.042091477490529,0.078911515795068,0.078911515795068,0.042091477490529,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.14794033605678,0.27735296695391,0.27735296695391,0.14794033605678,0.14794033605678,0.27735296695391,0.27735296695391,0.14794033605678,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.14794033605678,0.27735296695391,0.27735296695391,0.14794033605678,0.14794033605678,0.27735296695391,0.27735296695391,0.14794033605678,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.042091477490529,0.078911515795068,0.078911515795068,0.042091477490529,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.078911515795068,0.14794033605678,0.14794033605678,0.078911515795068,0.042091477490529,0.078911515795068,0.078911515795068,0.042091477490529};
221 
222 //1D Triangles, which are just the same as 1D Gauss elements, but scaled
223 //so their coordinate runs from 0 to 1, rather than -1 to 1
224 
225 //------------------------------------------------------------
226 //N=2
227 //------------------------------------------------------------
228 const double TGauss<1,2>::Knot[2][1]=
229 {{0.5*(-0.577350269189626+1.0)},{0.5*(0.577350269189626+1.0)}};
230 
231 const double TGauss<1,2>::Weight[2]={0.5,0.5};
232 
233 //------------------------------------------------------------
234 //N=3
235 //------------------------------------------------------------
236 const double TGauss<1,3>::Knot[3][1]=
237 {{0.5*(-0.774596669241483+1.0)},{0.5},{0.5*(0.774596669241483+1.0)}};
238 
239 const double TGauss<1,3>::Weight[3]={(5.0/18.0),(8.0/18.0),(5.0/18.0)};
240 
241 //------------------------------------------------------------
242 //N=4
243 //------------------------------------------------------------
244 const double TGauss<1,4>::Knot[4][1]=
245 {{0.5*(-0.861136311594053+1.0)},{0.5*(-0.339981043584856+1.0)},
246  {0.5*(0.339981043584856+1.0)}, {0.5*(0.861136311594053+1.0)}};
247 
248 const double TGauss<1,4>::Weight[4]={0.5*0.347854845137448,
249 0.5*0.652145154862546, 0.5*0.652145154862546, 0.5*0.347854845137448};
250 
251 //------------------------------------------------------------
252 //N=5
253 //------------------------------------------------------------
254 const double TGauss<1,5>::Knot[5][1]= {};
255 
256 const double TGauss<1,5>::Weight[5]= {};
257 
258 //// Define the positions and weights of the 2D Gauss points for triangles
259 //
260 //------------------------------------------------------------
261 // "Full integration" weights for linear triangles
262 // accurate up to second order (Bathe p 467)
263 //------------------------------------------------------------
264 const double TGauss<2,2>::Knot[3][2]={
265  {0.1666666666667,0.1666666666667},
266  {0.6666666666667,0.1666666666667},
267  {0.1666666666667,0.6666666666667}};
268 
269 const double TGauss<2,2>::Weight[3] = {
270 0.1666666666667,
271 0.1666666666667,
272 0.1666666666667};
273 
274 //------------------------------------------------------------
275 // "Full integration" weights for quadratic triangles
276 // accurate up to fifth order (Bathe p 467)
277 //------------------------------------------------------------
278 const double TGauss<2,3>::Knot[7][2]={
279  {0.1012865073235,0.1012865073235},
280  {0.7974269853531,0.1012865073235},
281  {0.1012865073235,0.7974269853531},
282  {0.4701420641051,0.0597158717898},
283  {0.4701420641051,0.4701420641051},
284  {0.0597158717898,0.4701420641051},
285  {0.3333333333333,0.3333333333333}};
286 
287 const double TGauss<2,3>::Weight[7] = {
288 0.5*0.1259391805448,
289 0.5*0.1259391805448,
290 0.5*0.1259391805448,
291 0.5*0.1323941527885,
292 0.5*0.1323941527885,
293 0.5*0.1323941527885,
294 0.5*0.225};
295 
296 
297 
298 //------------------------------------------------------------
299 //"Full integration" weights for cubic triangles
300 // accurate up to seventh order (Bathe p 467)
301 //------------------------------------------------------------
302 const double TGauss<2,4>::Knot[13][2]={
303  {0.0651301029022,0.0651301029022},
304  {0.8697397941956,0.0651301029022},
305  {0.0651301029022,0.8697397941956},
306  {0.3128654960049,0.0486903154253},
307  {0.6384441885698,0.3128654960049},
308  {0.0486903154253,0.6384441885698},
309  {0.6384441885698,0.0486903154253},
310  {0.3128654960049,0.6384441885698},
311  {0.0486903154253,0.3128654960049},
312  {0.2603459660790,0.2603459660790},
313  {0.4793080678419,0.2603459660790},
314  {0.2603459660790,0.4793080678419},
315  {0.3333333333333,0.3333333333333}};
316 
317 
318 const double TGauss<2,4>::Weight[13] = {
319 0.5*0.0533472356088,
320 0.5*0.0533472356088,
321 0.5*0.0533472356088,
322 0.5*0.0771137608903,
323 0.5*0.0771137608903,
324 0.5*0.0771137608903,
325 0.5*0.0771137608903,
326 0.5*0.0771137608903,
327 0.5*0.0771137608903,
328 0.5*0.1756152574332,
329 0.5*0.1756152574332,
330 0.5*0.1756152574332,
331 0.5*-0.1495700444677};
332 
333 //------------------------------------------------------------
334 //"Full integration" weights for 2D triangles
335 // accurate up to order 11
336 // http://people.sc.fsu.edu/~jburkardt/datasets/quadrature_rules_tri/quadrature_rules_tri.html
337 // [TOMS706_37, order 37, degree of precision 13, a rule from ACM TOMS algorithm #706.]
338 //------------------------------------------------------------
339  const double TGauss<2,13>::Knot[37][2]={
340  {0.333333333333333333333333333333, 0.333333333333333333333333333333},
341  {0.950275662924105565450352089520, 0.024862168537947217274823955239},
342  {0.024862168537947217274823955239, 0.950275662924105565450352089520},
343  {0.024862168537947217274823955239, 0.024862168537947217274823955239},
344  {0.171614914923835347556304795551, 0.414192542538082326221847602214},
345  {0.414192542538082326221847602214, 0.171614914923835347556304795551},
346  {0.414192542538082326221847602214, 0.414192542538082326221847602214},
347  {0.539412243677190440263092985511, 0.230293878161404779868453507244},
348  {0.230293878161404779868453507244, 0.539412243677190440263092985511},
349  {0.230293878161404779868453507244, 0.230293878161404779868453507244},
350  {0.772160036676532561750285570113, 0.113919981661733719124857214943},
351  {0.113919981661733719124857214943, 0.772160036676532561750285570113},
352  {0.113919981661733719124857214943, 0.113919981661733719124857214943},
353  {0.009085399949835353883572964740, 0.495457300025082323058213517632},
354  {0.495457300025082323058213517632, 0.009085399949835353883572964740},
355  {0.495457300025082323058213517632, 0.495457300025082323058213517632},
356  {0.062277290305886993497083640527, 0.468861354847056503251458179727},
357  {0.468861354847056503251458179727, 0.062277290305886993497083640527},
358  {0.468861354847056503251458179727, 0.468861354847056503251458179727},
359  {0.022076289653624405142446876931, 0.851306504174348550389457672223},
360  {0.022076289653624405142446876931, 0.126617206172027096933163647918},
361  {0.851306504174348550389457672223, 0.022076289653624405142446876931},
362  {0.851306504174348550389457672223, 0.126617206172027096933163647918},
363  {0.126617206172027096933163647918, 0.022076289653624405142446876931},
364  {0.126617206172027096933163647918, 0.851306504174348550389457672223},
365  {0.018620522802520968955913511549, 0.689441970728591295496647976487},
366  {0.018620522802520968955913511549, 0.291937506468887771754472382212},
367  {0.689441970728591295496647976487, 0.018620522802520968955913511549},
368  {0.689441970728591295496647976487, 0.291937506468887771754472382212},
369  {0.291937506468887771754472382212, 0.018620522802520968955913511549},
370  {0.291937506468887771754472382212, 0.689441970728591295496647976487},
371  {0.096506481292159228736516560903, 0.635867859433872768286976979827},
372  {0.096506481292159228736516560903, 0.267625659273967961282458816185},
373  {0.635867859433872768286976979827, 0.096506481292159228736516560903},
374  {0.635867859433872768286976979827, 0.267625659273967961282458816185},
375  {0.267625659273967961282458816185, 0.096506481292159228736516560903},
376  {0.267625659273967961282458816185, 0.635867859433872768286976979827}
377  };
378 
379  // Modified by DR 2017 to include factor 1/2
380  const double TGauss<2,13>::Weight[37] =
381  {
382  0.5*0.051739766065744133555179145422,
383  0.5*0.008007799555564801597804123460,
384  0.5*0.008007799555564801597804123460,
385  0.5*0.008007799555564801597804123460,
386  0.5*0.046868898981821644823226732071,
387  0.5*0.046868898981821644823226732071,
388  0.5*0.046868898981821644823226732071,
389  0.5*0.046590940183976487960361770070,
390  0.5*0.046590940183976487960361770070,
391  0.5*0.046590940183976487960361770070,
392  0.5*0.031016943313796381407646220131,
393  0.5*0.031016943313796381407646220131,
394  0.5*0.031016943313796381407646220131,
395  0.5*0.010791612736631273623178240136,
396  0.5*0.010791612736631273623178240136,
397  0.5*0.010791612736631273623178240136,
398  0.5*0.032195534242431618819414482205,
399  0.5*0.032195534242431618819414482205,
400  0.5*0.032195534242431618819414482205,
401  0.5*0.015445834210701583817692900053,
402  0.5*0.015445834210701583817692900053,
403  0.5*0.015445834210701583817692900053,
404  0.5*0.015445834210701583817692900053,
405  0.5*0.015445834210701583817692900053,
406  0.5*0.015445834210701583817692900053,
407  0.5*0.017822989923178661888748319485,
408  0.5*0.017822989923178661888748319485,
409  0.5*0.017822989923178661888748319485,
410  0.5*0.017822989923178661888748319485,
411  0.5*0.017822989923178661888748319485,
412  0.5*0.017822989923178661888748319485,
413  0.5*0.037038683681384627918546472190,
414  0.5*0.037038683681384627918546472190,
415  0.5*0.037038683681384627918546472190,
416  0.5*0.037038683681384627918546472190,
417  0.5*0.037038683681384627918546472190,
418  0.5*0.037038683681384627918546472190
419  };
420 
421 //------------------------------------------------------------
422 //"Full integration" weights for 2D triangles
423 // accurate up to order 16
424 // http://people.sc.fsu.edu/~jburkardt/datasets/quadrature_rules_tri/
425 // quadrature_rules_tri.html
426 // This uses the order 16 Dunavant rule, computed using software available from
427 // https://people.sc.fsu.edu/~jburkardt/cpp_src/triangle_dunavant_rule/triangle_dunavant_rule.html
428 // This integegration scheme is accurate up to order 16 and uses 52 points
429 //------------------------------------------------------------
430  const double TGauss<2,16>::Knot[52][2]={
431  {0.333333333333333, 0.333333333333333 },
432  {0.005238916103123, 0.497380541948438 },
433  {0.497380541948438, 0.497380541948438 },
434  {0.497380541948438, 0.005238916103123 },
435  {0.173061122901295, 0.413469438549352 },
436  {0.413469438549352, 0.413469438549352 },
437  {0.413469438549352, 0.173061122901295 },
438  {0.059082801866017, 0.470458599066991 },
439  {0.470458599066991, 0.470458599066991 },
440  {0.470458599066991, 0.059082801866017 },
441  {0.518892500060958, 0.240553749969521 },
442  {0.240553749969521, 0.240553749969521 },
443  {0.240553749969521, 0.518892500060958 },
444  {0.704068411554854, 0.147965794222573 },
445  {0.147965794222573, 0.147965794222573 },
446  {0.147965794222573, 0.704068411554854 },
447  {0.849069624685052, 0.07546518765747399 },
448  {0.07546518765747399,0.07546518765747399},
449  {0.07546518765747399,0.849069624685052 },
450  {0.96680719475395,0.016596402623025 },
451  {0.016596402623025, 0.016596402623025 },
452  {0.016596402623025, 0.96680719475395 },
453  {0.103575692245252, 0.296555596579887 },
454  {0.296555596579887, 0.599868711174861 },
455  {0.599868711174861, 0.103575692245252 },
456  {0.296555596579887, 0.103575692245252 },
457  {0.599868711174861, 0.296555596579887 },
458  {0.103575692245252, 0.599868711174861 },
459  {0.020083411655416, 0.337723063403079 },
460  {0.337723063403079, 0.6421935249415049 },
461  {0.6421935249415049,0.020083411655416 },
462  {0.337723063403079, 0.020083411655416 },
463  {0.6421935249415049, 0.337723063403079 },
464  {0.020083411655416, 0.6421935249415049 },
465  {-0.004341002614139, 0.204748281642812 },
466  {0.204748281642812, 0.7995927209713271 },
467  {0.7995927209713271, -0.004341002614139 },
468  {0.204748281642812, -0.004341002614139 },
469  {0.7995927209713271, 0.204748281642812 },
470  {-0.004341002614139, 0.7995927209713271 },
471  {0.04194178646801, 0.189358492130623 },
472  {0.189358492130623, 0.768699721401368 },
473  {0.768699721401368, 0.04194178646801 },
474  {0.189358492130623, 0.04194178646801 },
475  {0.768699721401368, 0.189358492130623 },
476  {0.04194178646801, 0.768699721401368 },
477  {0.014317320230681, 0.08528361568265699 },
478  {0.08528361568265699, 0.900399064086661 },
479  {0.900399064086661, 0.014317320230681 },
480  {0.08528361568265699, 0.014317320230681 },
481  {0.900399064086661, 0.08528361568265699 },
482  {0.014317320230681, 0.900399064086661 }
483  };
484  // Modified by DR 2017 to include factor 1/2
485  const double TGauss<2,16>::Weight[52] =
486  {
487  2* 0.046875697427642,
488  2* 0.006405878578585,
489  2* 0.006405878578585,
490  2* 0.006405878578585,
491  2* 0.041710296739387,
492  2* 0.041710296739387,
493  2* 0.041710296739387,
494  2* 0.026891484250064,
495  2* 0.026891484250064,
496  2* 0.026891484250064,
497  2* 0.04213252276165 ,
498  2* 0.04213252276165 ,
499  2* 0.04213252276165 ,
500  2* 0.030000266842773,
501  2* 0.030000266842773,
502  2* 0.030000266842773,
503  2* 0.014200098925024,
504  2* 0.014200098925024,
505  2* 0.014200098925024,
506  2* 0.003582462351273,
507  2* 0.003582462351273,
508  2* 0.003582462351273,
509  2* 0.032773147460627,
510  2* 0.032773147460627,
511  2* 0.032773147460627,
512  2* 0.032773147460627,
513  2* 0.032773147460627,
514  2* 0.032773147460627,
515  2* 0.015298306248441,
516  2* 0.015298306248441,
517  2* 0.015298306248441,
518  2* 0.015298306248441,
519  2* 0.015298306248441,
520  2* 0.015298306248441,
521  2* 0.002386244192839,
522  2* 0.002386244192839,
523  2* 0.002386244192839,
524  2* 0.002386244192839,
525  2* 0.002386244192839,
526  2* 0.002386244192839,
527  2* 0.019084792755899,
528  2* 0.019084792755899,
529  2* 0.019084792755899,
530  2* 0.019084792755899,
531  2* 0.019084792755899,
532  2* 0.019084792755899,
533  2* 0.006850054546542,
534  2* 0.006850054546542,
535  2* 0.006850054546542,
536  2* 0.006850054546542,
537  2* 0.006850054546542,
538  2* 0.006850054546542
539  };
540 //------------------------------------------------------------
541 //"Full integration" weights for 2D triangles
542 // accurate up to order 15
543 // http://people.sc.fsu.edu/~jburkardt/datasets/quadrature_rules_tri/quadrature_rules_tri.html
544 // [GAUSS8X8, order 64, degree of precision 15, (essentially a product of two 8 point 1D Gauss-Legendre rules).]
545 //------------------------------------------------------------
546  const double TGauss<2,5>::Knot[64][2]={
547  {0.9553660447100000, 0.8862103848242247e-3},
548  {0.9553660447100000, 0.4537789678039195e-2},
549  {0.9553660447100000, 0.1058868260117431e-1},
550  {0.9553660447100000, 0.1822327082910602e-1},
551  {0.9553660447100000, 0.2641068446089399e-1},
552  {0.9553660447100000, 0.3404527268882569e-1},
553  {0.9553660447100000, 0.4009616561196080e-1},
554  {0.9553660447100000, 0.4374774490517578e-1},
555  {0.8556337429600001, 0.2866402391985981e-2},
556  {0.8556337429600001, 0.1467724979327651e-1},
557  {0.8556337429600001, 0.3424855503358430e-1},
558  {0.8556337429600001, 0.5894224214571626e-1},
559  {0.8556337429600001, 0.8542401489428375e-1},
560  {0.8556337429600001, 0.1101177020064157},
561  {0.8556337429600001, 0.1296890072467235},
562  {0.8556337429600001, 0.1414998546480140},
563  {0.7131752428600000, 0.5694926133044352e-2},
564  {0.7131752428600000, 0.2916054411712861e-1},
565  {0.7131752428600000, 0.6804452564827500e-1},
566  {0.7131752428600000, 0.1171055801775613},
567  {0.7131752428600000, 0.1697191769624387},
568  {0.7131752428600000, 0.2187802314917250},
569  {0.7131752428600000, 0.2576642130228714},
570  {0.7131752428600000, 0.2811298310069557},
571  {0.5451866848000000, 0.9030351006711630e-2},
572  {0.5451866848000000, 0.4623939674940125e-1},
573  {0.5451866848000000, 0.1078970888004545},
574  {0.5451866848000000, 0.1856923986620134},
575  {0.5451866848000000, 0.2691209165379867},
576  {0.5451866848000000, 0.3469162263995455},
577  {0.5451866848000000, 0.4085739184505988},
578  {0.5451866848000000, 0.4457829641932884},
579  {0.3719321645800000, 0.1247033193690498e-1},
580  {0.3719321645800000, 0.6385362269957356e-1},
581  {0.3719321645800000, 0.1489989161403976},
582  {0.3719321645800000, 0.2564292182833579},
583  {0.3719321645800000, 0.3716386171366422},
584  {0.3719321645800000, 0.4790689192796024},
585  {0.3719321645800000, 0.5642142127204264},
586  {0.3719321645800000, 0.6155975034830951},
587  {0.2143084794000000, 0.1559996151584746e-1},
588  {0.2143084794000000, 0.7987871227492103e-1},
589  {0.2143084794000000, 0.1863925811641285},
590  {0.2143084794000000, 0.3207842387034378},
591  {0.2143084794000000, 0.4649072818965623},
592  {0.2143084794000000, 0.5992989394358715},
593  {0.2143084794000000, 0.7058128083250790},
594  {0.2143084794000000, 0.7700915590841526},
595  {0.9132360790000005e-1, 0.1804183496379599e-1},
596  {0.9132360790000005e-1, 0.9238218584838476e-1},
597  {0.9132360790000005e-1, 0.2155687489628060},
598  {0.9132360790000005e-1, 0.3709968314854498},
599  {0.9132360790000005e-1, 0.5376795606145502},
600  {0.9132360790000005e-1, 0.6931076431371940},
601  {0.9132360790000005e-1, 0.8162942062516152},
602  {0.9132360790000005e-1, 0.8906345571362040},
603  {0.1777991514999999e-1, 0.1950205026019779e-1},
604  {0.1777991514999999e-1, 0.9985913490381848e-1},
605  {0.1777991514999999e-1, 0.2330157982952792},
606  {0.1777991514999999e-1, 0.4010234473667467},
607  {0.1777991514999999e-1, 0.5811966374832533},
608  {0.1777991514999999e-1, 0.7492042865547208},
609  {0.1777991514999999e-1, 0.8823609499461815},
610  {0.1777991514999999e-1, 0.9627180345898023}
611  };
612  // Modified by DR 2017 to include factor 1/2
613  const double TGauss<2,5>::Weight[64] =
614  {
615  0.5*0.3335674062677772e-3,
616  0.5*0.7327880811491046e-3,
617  0.5*0.1033723454167925e-2,
618  0.5*0.1195112498415193e-2,
619  0.5*0.1195112498415193e-2,
620  0.5*0.1033723454167925e-2,
621  0.5*0.7327880811491046e-3,
622  0.5*0.3335674062677772e-3,
623  0.5*0.1806210919443461e-2,
624  0.5*0.3967923151181667e-2,
625  0.5*0.5597437146194232e-2,
626  0.5*0.6471331443180639e-2,
627  0.5*0.6471331443180639e-2,
628  0.5*0.5597437146194232e-2,
629  0.5*0.3967923151181667e-2,
630  0.5*0.1806210919443461e-2,
631  0.5*0.4599755803015752e-2,
632  0.5*0.1010484287526739e-1,
633  0.5*0.1425461651131868e-1,
634  0.5*0.1648010431039818e-1,
635  0.5*0.1648010431039818e-1,
636  0.5*0.1425461651131868e-1,
637  0.5*0.1010484287526739e-1,
638  0.5*0.4599755803015752e-2,
639  0.5*0.8017259531156730e-2,
640  0.5*0.1761248886287915e-1,
641  0.5*0.2484544071087993e-1,
642  0.5*0.2872441038508419e-1,
643  0.5*0.2872441038508419e-1,
644  0.5*0.2484544071087993e-1,
645  0.5*0.1761248886287915e-1,
646  0.5*0.8017259531156730e-2,
647  0.5*0.1073501897357062e-1,
648  0.5*0.2358292149331603e-1,
649  0.5*0.3326776143412911e-1,
650  0.5*0.3846165753898425e-1,
651  0.5*0.3846165753898425e-1,
652  0.5*0.3326776143412911e-1,
653  0.5*0.2358292149331603e-1,
654  0.5*0.1073501897357062e-1,
655  0.5*0.1138879740452669e-1,
656  0.5*0.2501915606814251e-1,
657  0.5*0.3529381699354388e-1,
658  0.5*0.4080402900378691e-1,
659  0.5*0.4080402900378691e-1,
660  0.5*0.3529381699354388e-1,
661  0.5*0.2501915606814251e-1,
662  0.5*0.1138879740452669e-1,
663  0.5*0.9223845391285393e-2,
664  0.5*0.2026314273544469e-1,
665  0.5*0.2858464328177232e-1,
666  0.5*0.3304739223149761e-1,
667  0.5*0.3304739223149761e-1,
668  0.5*0.2858464328177232e-1,
669  0.5*0.2026314273544469e-1,
670  0.5*0.9223845391285393e-2,
671  0.5*0.4509812715921713e-2,
672  0.5*0.9907253959306707e-2,
673  0.5*0.1397588340693756e-1,
674  0.5*0.1615785427783403e-1,
675  0.5*0.1615785427783403e-1,
676  0.5*0.1397588340693756e-1,
677  0.5*0.9907253959306707e-2,
678  0.5*0.4509812715921713e-2
679  };
680 
681 //------------------------------------------------------------
682 //"Full integration" weights for 2D triangles
683 // accurate up to points 19, degree of precision 8, a rule from ACM TOMS
684 // algorithm #584.
685 // http://people.sc.fsu.edu/~jburkardt/datasets/quadrature_rules_tri/quadrature_rules_tri.html
686 // TOMS589_19 a 19 point Gauss rule accurate to order 8
687 //------------------------------------------------------------
688  const double TGauss<2,9>::Knot[19][2]={
689  {0.3333333333333333, 0.3333333333333333},
690  {0.7974269853530872, 0.1012865073234563},
691  {0.1012865073234563, 0.7974269853530872},
692  {0.1012865073234563, 0.1012865073234563},
693  {0.0597158717897698, 0.4701420641051151},
694  {0.4701420641051151, 0.0597158717897698},
695  {0.4701420641051151, 0.4701420641051151},
696  {0.5357953464498992, 0.2321023267750504},
697  {0.2321023267750504, 0.5357953464498992},
698  {0.2321023267750504, 0.2321023267750504},
699  {0.9410382782311209, 0.0294808608844396},
700  {0.0294808608844396, 0.9410382782311209},
701  {0.0294808608844396, 0.0294808608844396},
702  {0.7384168123405100, 0.2321023267750504},
703  {0.7384168123405100, 0.0294808608844396},
704  {0.2321023267750504, 0.7384168123405100},
705  {0.2321023267750504, 0.0294808608844396},
706  {0.0294808608844396, 0.7384168123405100},
707  {0.0294808608844396, 0.2321023267750504}
708  };
709 
710  const double TGauss<2,9>::Weight[19]={
711  2*0.0378610912003147,
712  2*0.0376204254131829,
713  2*0.0376204254131829,
714  2*0.0376204254131829,
715  2*0.0783573522441174,
716  2*0.0783573522441174,
717  2*0.0783573522441174,
718  2*0.1162714796569659,
719  2*0.1162714796569659,
720  2*0.1162714796569659,
721  2*0.0134442673751655,
722  2*0.0134442673751655,
723  2*0.0134442673751655,
724  2*0.0375097224552317,
725  2*0.0375097224552317,
726  2*0.0375097224552317,
727  2*0.0375097224552317,
728  2*0.0375097224552317,
729  2*0.0375097224552317
730  };
731 
732 //------------------------------------------------------------
733 //// Define the positions and weights of the 3D Gauss points for tets
734 
735 //------------------------------------------------------------
736 // "Full integration" weights for linear tets
737 // accurate up to second order (e.g. from German Zienkiwicz p 200)
738 //------------------------------------------------------------
739 const double TGauss<3,2>::Knot[4][3]={
740  {0.138196601125011,0.138196601125011,0.585410196624969},
741  {0.138196601125011,0.585410196624969,0.138196601125011},
742  {0.585410196624969,0.138196601125011,0.138196601125011},
743  {0.138196601125011,0.138196601125011,0.138196601125011}};
744 
745 
746 const double TGauss<3,2>::Weight[4] = {
747 0.0416666666667,
748 0.0416666666667,
749 0.0416666666667,
750 0.0416666666667};
751 
752 
753 
754 
755 
756 //------------------------------------------------------------
757 //"Full integration" weights for quadratic tets
758 // accurate up to fifth order
759 // The numbers are from Keast CMAME 55 pp339-348 (1986)
760 //------------------------------------------------------------
761 const double TGauss<3,3>::Knot[11][3]={
762  {0.25,0.25,0.25},
763  {0.785714285714286,0.071428571428571,0.071428571428571},
764  {0.071428571428571,0.071428571428571,0.071428571428571},
765  {0.071428571428571,0.785714285714286,0.071428571428571},
766  {0.071428571428571,0.071428571428571,0.785714285714286},
767  {0.399403576166799,0.399403576166799,0.100596423833201},
768  {0.399403576166799,0.100596423833201,0.399403576166799},
769  {0.100596423833201,0.399403576166799,0.399403576166799},
770  {0.399403576166799,0.100596423833201,0.100596423833201},
771  {0.100596423833201,0.399403576166799,0.100596423833201},
772  {0.100596423833201,0.100596423833201,0.399403576166799}};
773 
774 
775 const double TGauss<3,3>::Weight[11] ={
776 -0.01315555555556,
777  0.00762222222222,
778  0.00762222222222,
779  0.00762222222222,
780  0.00762222222222,
781  0.02488888888889,
782  0.02488888888889,
783  0.02488888888889,
784  0.02488888888889,
785  0.02488888888889,
786  0.02488888888889
787 };
788 
789 //------------------------------------------------------------
790 //"Full integration" weights for quartic tets
791 // accurate up to eighth order
792 // The numbers are from Keast CMAME 55 pp339-348 (1986)
793 //------------------------------------------------------------
794 const double TGauss<3,5>::Knot[45][3]={
795  {2.50000000000000000e-01, 2.50000000000000000e-01, 2.50000000000000000e-01},
796  {1.27470936566639015e-01, 1.27470936566639015e-01, 1.27470936566639015e-01},
797  {1.27470936566639015e-01, 1.27470936566639015e-01, 6.17587190300082967e-01},
798  {1.27470936566639015e-01, 6.17587190300082967e-01, 1.27470936566639015e-01},
799  {6.17587190300082967e-01, 1.27470936566639015e-01, 1.27470936566639015e-01},
800  {3.20788303926322960e-02, 3.20788303926322960e-02, 3.20788303926322960e-02},
801  {3.20788303926322960e-02, 3.20788303926322960e-02, 9.03763508822103123e-01},
802  {3.20788303926322960e-02, 9.03763508822103123e-01, 3.20788303926322960e-02},
803  {9.03763508822103123e-01, 3.20788303926322960e-02, 3.20788303926322960e-02},
804  {4.97770956432810185e-02, 4.97770956432810185e-02, 4.50222904356718978e-01},
805  {4.97770956432810185e-02, 4.50222904356718978e-01, 4.50222904356718978e-01},
806  {4.50222904356718978e-01, 4.50222904356718978e-01, 4.97770956432810185e-02},
807  {4.50222904356718978e-01, 4.97770956432810185e-02, 4.97770956432810185e-02},
808  {4.97770956432810185e-02, 4.50222904356718978e-01, 4.97770956432810185e-02},
809  {4.50222904356718978e-01, 4.97770956432810185e-02, 4.50222904356718978e-01},
810  {1.83730447398549945e-01, 1.83730447398549945e-01, 3.16269552601450060e-01},
811  {1.83730447398549945e-01, 3.16269552601450060e-01, 3.16269552601450060e-01},
812  {3.16269552601450060e-01, 3.16269552601450060e-01, 1.83730447398549945e-01},
813  {3.16269552601450060e-01, 1.83730447398549945e-01, 1.83730447398549945e-01},
814  {1.83730447398549945e-01, 3.16269552601450060e-01, 1.83730447398549945e-01},
815  {3.16269552601450060e-01, 1.83730447398549945e-01, 3.16269552601450060e-01},
816  {2.31901089397150906e-01, 2.31901089397150906e-01, 2.29177878448171174e-02},
817  {2.31901089397150906e-01, 2.29177878448171174e-02, 5.13280033360881072e-01},
818  {2.29177878448171174e-02, 5.13280033360881072e-01, 2.31901089397150906e-01},
819  {5.13280033360881072e-01, 2.31901089397150906e-01, 2.31901089397150906e-01},
820  {2.31901089397150906e-01, 5.13280033360881072e-01, 2.31901089397150906e-01},
821  {5.13280033360881072e-01, 2.31901089397150906e-01, 2.29177878448171174e-02},
822  {2.31901089397150906e-01, 2.29177878448171174e-02, 2.31901089397150906e-01},
823  {2.31901089397150906e-01, 5.13280033360881072e-01, 2.29177878448171174e-02},
824  {2.29177878448171174e-02, 2.31901089397150906e-01, 5.13280033360881072e-01},
825  {5.13280033360881072e-01, 2.29177878448171174e-02, 2.31901089397150906e-01},
826  {2.29177878448171174e-02, 2.31901089397150906e-01, 2.31901089397150906e-01},
827  {2.31901089397150906e-01, 2.31901089397150906e-01, 5.13280033360881072e-01},
828  {3.79700484718286102e-02, 3.79700484718286102e-02, 7.30313427807538396e-01},
829  {3.79700484718286102e-02, 7.30313427807538396e-01, 1.93746475248804382e-01},
830  {7.30313427807538396e-01, 1.93746475248804382e-01, 3.79700484718286102e-02},
831  {1.93746475248804382e-01, 3.79700484718286102e-02, 3.79700484718286102e-02},
832  {3.79700484718286102e-02, 1.93746475248804382e-01, 3.79700484718286102e-02},
833  {1.93746475248804382e-01, 3.79700484718286102e-02, 7.30313427807538396e-01},
834  {3.79700484718286102e-02, 7.30313427807538396e-01, 3.79700484718286102e-02},
835  {3.79700484718286102e-02, 1.93746475248804382e-01, 7.30313427807538396e-01},
836  {7.30313427807538396e-01, 3.79700484718286102e-02, 1.93746475248804382e-01},
837  {1.93746475248804382e-01, 7.30313427807538396e-01, 3.79700484718286102e-02},
838  {7.30313427807538396e-01, 3.79700484718286102e-02, 3.79700484718286102e-02},
839  {3.79700484718286102e-02, 3.79700484718286102e-02, 1.93746475248804382e-01}
840 };
841 
842 const double TGauss<3,5>::Weight[45] ={
843 -3.93270066412926145e-02, 4.08131605934270525e-03,
844  4.08131605934270525e-03, 4.08131605934270525e-03,
845  4.08131605934270525e-03, 6.58086773304341943e-04,
846  6.58086773304341943e-04, 6.58086773304341943e-04,
847  6.58086773304341943e-04, 4.38425882512284693e-03,
848  4.38425882512284693e-03, 4.38425882512284693e-03,
849  4.38425882512284693e-03, 4.38425882512284693e-03,
850  4.38425882512284693e-03, 1.38300638425098166e-02,
851  1.38300638425098166e-02, 1.38300638425098166e-02,
852  1.38300638425098166e-02, 1.38300638425098166e-02,
853  1.38300638425098166e-02, 4.24043742468372453e-03,
854  4.24043742468372453e-03, 4.24043742468372453e-03,
855  4.24043742468372453e-03, 4.24043742468372453e-03,
856  4.24043742468372453e-03, 4.24043742468372453e-03,
857  4.24043742468372453e-03, 4.24043742468372453e-03,
858  4.24043742468372453e-03, 4.24043742468372453e-03,
859  4.24043742468372453e-03, 2.23873973961420164e-03,
860  2.23873973961420164e-03, 2.23873973961420164e-03,
861  2.23873973961420164e-03, 2.23873973961420164e-03,
862  2.23873973961420164e-03, 2.23873973961420164e-03,
863  2.23873973961420164e-03, 2.23873973961420164e-03,
864  2.23873973961420164e-03, 2.23873973961420164e-03,
865  2.23873973961420164e-03};
866 
867 }
868 
869 
870 
871 
double weight(const unsigned &i) const
Return weight of i-th integration point.
Definition: integral.cc:79
double fraction_to_local(const double &fraction) const
Helper function to convert from fraction [0,1] to local coordinate (s in [-1,1] for quads) ...
Definition: integral.cc:74
cstr elem_len * i
Definition: cfortran.h:607
const unsigned Index_s
Index of local coord (s) which scheme will be reflected along.
Definition: integral.h:145
const double S_min
Minimum value of local coord (s) over whole domain.
Definition: integral.h:151
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
double local_to_fraction(const double &local_coord) const
Helper function to convert from local coordinate (s in [-1,1] for quads) to fraction [0...
Definition: integral.cc:69
double knot(const unsigned &i, const unsigned &j) const
Return local coordinate s[j] of i-th integration point.
Definition: integral.cc:38
static char t char * s
Definition: cfortran.h:572
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
const Integral * Integral_pt
Pointer to integral which we will reflect.
Definition: integral.h:142
const double S_reflect_fraction
Position where scheme will be reflected as of fraction of domain [0,1].
Definition: integral.h:157
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
const double S_max
Maximum value of local coord (s) over whole domain.
Definition: integral.h:154
unsigned nweight() const
Number of integration points of the scheme, two times the base scheme.
Definition: integral.h:123