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