MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
intrules_cut.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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#ifndef MFEM_CUTINTRULES
13#define MFEM_CUTINTRULES
14
15#include "../config/config.hpp"
16#include "../general/array.hpp"
17#include "intrules.hpp"
18#include "eltrans.hpp"
19#include "coefficient.hpp"
20
21namespace mfem
22{
23/**
24 @brief Abstract class for construction of IntegrationRules in cut elements.
25
26 Interface for construction of cut-surface and cut-volume IntegrationRules. The
27 cut is specified by the zero level set of a given Coefficient.
28*/
30{
31protected:
32 /// Order of the IntegrationRule.
33 int Order;
34 /// The zero level set of this Coefficient defines the cut surface.
36 /// Space order for the LS projection.
38
39 /** @brief Constructor to set up the generated cut IntegrationRules.
40
41 @param [in] order Order of the constructed IntegrationRule.
42 @param [in] lvlset Coefficient whose zero level set specifies the cut.
43 @param [in] lsO Polynomial degree for projecting the level-set
44 Coefficient to a GridFunction, which is used to
45 compute gradients and normals. */
46 CutIntegrationRules(int order, Coefficient& lvlset, int lsO = 2)
47 : Order(order), LvlSet(&lvlset), lsOrder(lsO)
48 { MFEM_VERIFY(order > 0 && lsO > 0, "Invalid input") }
49
50public:
51
52 /// Change the order of the constructed IntegrationRule.
53 virtual void SetOrder(int order);
54
55 /// Change the Coefficient whose zero level set specifies the cut.
56 virtual void SetLevelSetCoefficient(Coefficient &ls) { LvlSet = &ls; }
57
58 /// Change the polynomial degree for projecting the level set Coefficient
59 /// to a GridFunction, which is used to compute local gradients and normals.
60 virtual void SetLevelSetProjectionOrder(int order);
61
62 /**
63 @brief Construct a cut-surface IntegrationRule.
64
65 Construct an IntegrationRule to integrate on the surface given by the
66 already specified level set function, for the element given by @a Tr.
67
68 @param [in] Tr Specifies the IntegrationRule's associated mesh element.
69 @param [out] result IntegrationRule on the cut-surface
70 */
72 IntegrationRule& result) = 0;
73
74 /**
75 @brief Construct a cut-volume IntegrationRule.
76
77 Construct an IntegrationRule to integrate in the subdomain given by the
78 positive values of the already specified level set function, for the element
79 given by @a Tr.
80
81 @param [in] Tr Specifies the IntegrationRule's associated mesh element.
82 @param [out] result IntegrationRule for the cut-volume
83 @param [in] sir Corresponding IntegrationRule for the surface, which can
84 be used to avoid computations.
85 */
87 IntegrationRule& result,
88 const IntegrationRule* sir = NULL) = 0;
89
90 /**
91 @brief Compute transformation quadrature weights for surface integration.
92
93 Compute the transformation weights for integration over the cut-surface in
94 reference space.
95
96 @param [in] Tr Specifies the IntegrationRule's associated element.
97 @param [in] sir IntegrationRule defining the IntegrationPoints
98 @param [out] weights Vector containing the transformation weights.
99 */
101 const IntegrationRule &sir,
102 Vector &weights) = 0;
103
104 /// @brief Destructor of CutIntegrationRules
106};
107
108#ifdef MFEM_USE_LAPACK
109
110/**
111 @brief Class for subdomain IntegrationRules by means of moment-fitting
112
113 Class for subdomain (surface and subdomain) IntegrationRules by means of
114 moment-fitting. The class provides different functions to construct the
115 IntegrationRules. The construction is done as described in Mueller et al. in
116 "Highly accurate surface and volume integration on implicit domains by means of
117 moment-fitting" (2013, see
118 https://onlinelibrary.wiley.com/doi/full/10.1002/nme.4569).
119*/
121{
122protected:
123 /// @brief Space Dimension of the element
124 int dim;
125 /// @brief Number of divergence-free basis functions for surface integration
127 /// @brief Number of basis functions for volume integration
129 /// @brief IntegrationRule representing the reused IntegrationPoints
131 /// @brief SVD of the matrix for volumetric IntegrationRules
133 /// @brief Array of face integration points
135 /// @brief Column-wise Matrix for the face quadrature weights
137 /// @brief Indicates the already computed face IntegrationRules
139
140 /**
141 @brief Initialization for surface IntegrationRule
142
143 Initialize the members for computation of surface IntegrationRule. This is
144 called when the first IntegrationRule is computed or when Order or level-set
145 change.
146
147 @param [in] order Order of the IntegrationRule
148 @param [in] levelset level-set function defining the implicit interface
149 @param [in] lsO polynomial degree for approximation of level-set function
150 @param [in] Tr ElemenTransformation to initalize the members with
151 */
152 void InitSurface(int order, Coefficient& levelset, int lsO,
154
155 /**
156 @brief Initialization for volume IntegrationRule
157
158 Initialize the members for computation of surface IntegrationRule. This is
159 called when the first IntegrationRule is computed or when Order or level-set
160 change.
161
162 @param [in] order Order of the IntegrationRule
163 @param [in] levelset level-set function defining the implicit interface
164 @param [in] lsO polynomial degree for approximation of level-set function
165 @param [in] Tr ElemenTransformation to initalize the members with
166 */
167 void InitVolume(int order, Coefficient& levelset, int lsO,
169
170 /// @brief Initialize the MomentFittingIntRules
171 void Init(int order, Coefficient& levelset, int lsO)
172 { Order = order; LvlSet = &levelset; lsOrder = lsO; FaceWeightsComp = 0.;}
173
174 /// @brief Clear stored data of the MomentFittingIntRules
175 void Clear();
176
177 /**
178 @brief Compute the IntegrationRules on the faces
179
180 Compute the IntegrationRules on the (cut) faces of an element. These will
181 be saved and reused if possible.
182
183 @param [in] Tr ElementTransformation of the element the IntegrationRules on
184 the faces are computed
185 */
187
188 /**
189 @brief Compute 1D quadrature weights
190
191 Compute the quadrature weights for the 1D surface quadrature rule.
192
193 @param [in] Tr ElementTransformation of the current element
194 */
196
197 /**
198 @brief Compute the 1D quadrature weights
199
200 Compute the 1D quadrature weights for the volumetric subdomain quadrature
201 rule.
202
203 @param [in] Tr ElementTransformation of the current element
204 @param [in] sir corresponding IntegrationRule on surface
205 */
207 const IntegrationRule* sir);
208
209 /**
210 @brief Compute 2D quadrature weights
211
212 Compute the quadrature weights for the 2D surface quadrature rule by means
213 of moment-fitting. To construct the quadrature rule, special integrals are
214 reduced to integrals over the edges of the subcell where the level-set is
215 positive.
216
217 @param [in] Tr ElementTransformation of the current element
218 */
220
221 /**
222 @brief Compute the 2D quadrature weights
223
224 Compute the 2D quadrature weights for the volumetric subdomain quadrature
225 rule by means of moment-fitting. To construct the quadrature rule, special
226 integrals are reduced to integrals over the boundary of the subcell where
227 the level-set is positive.
228
229 @param [in] Tr ElementTransformation of the current element
230 @param [in] sir corresponding IntegrationRule on surface
231 */
233 const IntegrationRule* sir);
234
235 /**
236 @brief Compute 2D quadrature weights
237
238 Compute the quadrature weights for the 3D surface quadrature rule by means
239 of moment-fitting. To construct the quadrature rule, special integrals are
240 reduced to integrals over the edges of the subcell where the level-set is
241 positive.
242
243 @param [in] Tr ElementTransformation of the current element
244 */
246
247
248 /**
249 @brief Compute the 3D quadrature weights
250
251 Compute the 3D quadrature weights for the volumetric subdomain quadrature
252 rule by means of moment-fitting. To construct the quadrature rule, special
253 integrals are reduced to integrals over the boundary of the subcell where
254 the level-set is positive.
255
256 @param [in] Tr ElementTransformation of the current element
257 @param [in] sir corresponding IntegrationRule on surface
258 */
260 const IntegrationRule* sir);
261
262 /// @brief A divergence free basis on the element [-1,1]^2
263 void DivFreeBasis2D(const IntegrationPoint& ip, DenseMatrix& shape);
264 /// @brief A orthogonalized divergence free basis on the element [-1,1]^2
265 void OrthoBasis2D(const IntegrationPoint& ip, DenseMatrix& shape);
266 /// @brief A orthogonalized divergence free basis on the element [-1,1]^3
267 void OrthoBasis3D(const IntegrationPoint& ip, DenseMatrix& shape);
268 /// @brief A step of the modified Gram-Schmidt algorithm
269 void mGSStep(DenseMatrix& shape, DenseTensor& shapeMFN, int step);
270
271 /// @brief Monomial basis on the element [-1,1]^2
272 void Basis2D(const IntegrationPoint& ip, Vector& shape);
273 /// @brief Antiderivatives of the monomial basis on the element [-1,1]^2
274 void BasisAD2D(const IntegrationPoint& ip, DenseMatrix& shape);
275 /// @brief Monomial basis on the element [-1,1]^3
276 void Basis3D(const IntegrationPoint& ip, Vector& shape);
277 /// @brief Antiderivatives of the monomial basis on the element [-1,1]^3
278 void BasisAD3D(const IntegrationPoint& ip, DenseMatrix& shape);
279
280public:
281
282 /** @brief Constructor to set up the generated cut IntegrationRules.
283
284 @param [in] order Order of the constructed IntegrationRule.
285 @param [in] lvlset Coefficient whose zero level set specifies the cut.
286 @param [in] lsO Polynomial degree for projecting the level-set
287 Coefficient to a GridFunction, which is used to
288 compute gradients and normals. */
290 : CutIntegrationRules(order, lvlset, lsO),
291 dim(-1), nBasis(-1), nBasisVolume(-1), VolumeSVD(nullptr)
293
294 /// Change the order of the constructed IntegrationRule.
295 void SetOrder(int order) override;
296
299
300 /**
301 @brief Construct a cut-surface IntegrationRule.
302
303 Construct an IntegrationRule to integrate on the surface given by the
304 already specified level set function, for the element given by @a Tr.
305
306 @param [in] Tr Specifies the IntegrationRule's associated mesh element.
307 @param [out] result IntegrationRule on the cut-surface
308 */
310 IntegrationRule& result) override;
311
312 /**
313 @brief Construct a cut-volume IntegrationRule.
314
315 Construct an IntegrationRule to integrate in the subdomain given by the
316 positive values of the already specified level set function, for the element
317 given by @a Tr.
318
319 @param [in] Tr Specifies the IntegrationRule's associated mesh element.
320 @param [out] result IntegrationRule for the cut-volume
321 @param [in] sir Corresponding IntegrationRule for the surface, which can
322 be used to avoid computations.
323 */
325 IntegrationRule& result,
326 const IntegrationRule* sir = nullptr) override;
327
328 /**
329 @brief Compute transformation quadrature weights for surface integration.
330
331 Compute the transformation weights for integration over the cut-surface in
332 reference space.
333
334 @param [in] Tr Specifies the IntegrationRule's associated element.
335 @param [in] sir IntegrationRule defining the IntegrationPoints
336 @param [out] weights Vector containing the transformation weights.
337 */
339 const IntegrationRule &sir,
340 Vector &weights) override;
341
342 /// @brief Destructor of MomentFittingIntRules
343 ~MomentFittingIntRules() override { delete VolumeSVD; }
344};
345
346#endif //MFEM_USE_LAPACK
347
348namespace DivFreeBasis
349{
350
351/// @brief 3 dimensional divergence free basis functions on [-1,1]^3
352inline void GetDivFree3DBasis(const Vector& X, DenseMatrix& shape, int Order)
353{
354 int nBasis;
355 if (Order == 0) { nBasis = 3; }
356 else if (Order == 1) { nBasis = 11; }
357 else if (Order == 2) { nBasis = 26; }
358 else if (Order == 3) { nBasis = 50; }
359 else if (Order == 4) { nBasis = 85; }
360 else if (Order == 5) { nBasis = 133; }
361 else if (Order == 6) { nBasis = 196; }
362 else { nBasis = 276; Order = 7; }
363
364 shape.SetSize(nBasis, 3);
365 shape = 0.;
366
367 real_t xi = X(0);
368 real_t eta = X(1);
369 real_t nu = X(2);
370
371 switch (Order)
372 {
373 case 7:
374 shape(196,0) = -2.995357736356377*nu + 26.958219627207395*pow(nu,3.)
375 - 59.30808317985627*pow(nu,5.)
376 + 36.714527682768164*pow(nu,7.);
377 shape(197,1) = -2.995357736356377*nu + 26.958219627207395*pow(nu,3.)
378 - 59.30808317985627*pow(nu,5.)
379 + 36.714527682768164*pow(nu,7.);
380 shape(198,0) = -0.689981317681863*eta
381 + 14.489607671319124*eta*pow(nu,2.)
382 - 43.46882301395737*eta*pow(nu,4.)
383 + 31.87713687690207*eta*pow(nu,6.);
384 shape(199,1) = -0.689981317681863*eta
385 + 14.489607671319124*eta*pow(nu,2.)
386 - 43.46882301395737*eta*pow(nu,4.)
387 + 31.87713687690207*eta*pow(nu,6.);
388 shape(199,2) = 0.689981317681863*nu - 4.829869223773041*pow(nu,3.)
389 + 8.693764602791473*pow(nu,5.)
390 - 4.553876696700296*pow(nu,7.);
391 shape(200,0) = -2.4581457378987928*nu + 11.471346776861033*pow(nu,3.)
392 - 10.324212099174929*pow(nu,5.)
393 + 7.374437213696378*pow(eta,2.)*nu
394 - 34.4140403305831*pow(eta,2.)*pow(nu,3.)
395 + 30.972636297524787*pow(eta,2.)*pow(nu,5.);
396 shape(201,1) = -2.4581457378987928*nu + 11.471346776861033*pow(nu,3.)
397 - 10.324212099174929*pow(nu,5.)
398 + 7.374437213696378*pow(eta,2.)*nu
399 - 34.4140403305831*pow(eta,2.)*pow(nu,3.)
400 + 30.972636297524787*pow(eta,2.)*pow(nu,5.);
401 shape(201,2) = 0.49162914757975856*eta
402 - 7.374437213696378*eta*pow(nu,2.)
403 + 17.20702016529155*eta*pow(nu,4.)
404 - 10.324212099174929*eta*pow(nu,6.);
405 shape(202,0) = -1.5785117100452566*eta
406 + 15.785117100452565*eta*pow(nu,2.)
407 - 18.415969950527995*eta*pow(nu,4.)
408 + 2.6308528500754274*pow(eta,3.)
409 - 26.308528500754274*pow(eta,3.)*pow(nu,2.)
410 + 30.69328325087999*pow(eta,3.)*pow(nu,4.);
411 shape(203,1) = -1.5785117100452566*eta
412 + 15.785117100452565*eta*pow(nu,2.)
413 - 18.415969950527995*eta*pow(nu,4.)
414 + 2.6308528500754274*pow(eta,3.)
415 - 26.308528500754274*pow(eta,3.)*pow(nu,2.)
416 + 30.69328325087999*pow(eta,3.)*pow(nu,4.);
417 shape(203,2) = 1.5785117100452566*nu - 5.261705700150855*pow(nu,3.)
418 + 3.6831939901055986*pow(nu,5.)
419 - 7.8925585502262825*pow(eta,2.)*nu
420 + 26.308528500754274*pow(eta,2.)*pow(nu,3.)
421 - 18.415969950527995*pow(eta,2.)*pow(nu,5.);
422 shape(204,0) = -1.5785117100452566*nu + 2.6308528500754274*pow(nu,3.)
423 + 15.785117100452565*pow(eta,2.)*nu
424 - 26.308528500754274*pow(eta,2.)*pow(nu,3.)
425 - 18.415969950527995*pow(eta,4.)*nu
426 + 30.69328325087999*pow(eta,4.)*pow(nu,3.);
427 shape(205,1) = -1.5785117100452566*nu + 2.6308528500754274*pow(nu,3.)
428 + 15.785117100452565*pow(eta,2.)*nu
429 - 26.308528500754274*pow(eta,2.)*pow(nu,3.)
430 - 18.415969950527995*pow(eta,4.)*nu
431 + 30.69328325087999*pow(eta,4.)*pow(nu,3.);
432 shape(205,2) = 2.6308528500754274*eta
433 - 15.785117100452565*eta*pow(nu,2.)
434 + 13.154264250377137*eta*pow(nu,4.)
435 - 6.138656650175998*pow(eta,3.)
436 + 36.83193990105599*pow(eta,3.)*pow(nu,2.)
437 - 30.69328325087999*pow(eta,3.)*pow(nu,4.);
438 shape(206,0) = -2.4581457378987928*eta
439 + 7.374437213696378*eta*pow(nu,2.)
440 + 11.471346776861033*pow(eta,3.)
441 - 34.4140403305831*pow(eta,3.)*pow(nu,2.)
442 - 10.324212099174929*pow(eta,5.)
443 + 30.972636297524787*pow(eta,5.)*pow(nu,2.);
444 shape(207,1) = -2.4581457378987928*eta
445 + 7.374437213696378*eta*pow(nu,2.)
446 + 11.471346776861033*pow(eta,3.)
447 - 34.4140403305831*pow(eta,3.)*pow(nu,2.)
448 - 10.324212099174929*pow(eta,5.)
449 + 30.972636297524787*pow(eta,5.)*pow(nu,2.);
450 shape(207,2) = 2.4581457378987928*nu - 2.4581457378987928*pow(nu,3.)
451 - 34.4140403305831*pow(eta,2.)*nu
452 + 34.4140403305831*pow(eta,2.)*pow(nu,3.)
453 + 51.621060495874644*pow(eta,4.)*nu
454 - 51.621060495874644*pow(eta,4.)*pow(nu,3.);
455 shape(208,0) = -0.689981317681863*nu
456 + 14.489607671319124*pow(eta,2.)*nu
457 - 43.46882301395737*pow(eta,4.)*nu
458 + 31.87713687690207*pow(eta,6.)*nu;
459 shape(209,1) = -0.689981317681863*nu
460 + 14.489607671319124*pow(eta,2.)*nu
461 - 43.46882301395737*pow(eta,4.)*nu
462 + 31.87713687690207*pow(eta,6.)*nu;
463 shape(209,2) = 4.829869223773041*eta
464 - 14.489607671319124*eta*pow(nu,2.)
465 - 28.97921534263825*pow(eta,3.)
466 + 86.93764602791474*pow(eta,3.)*pow(nu,2.)
467 + 31.87713687690207*pow(eta,5.)
468 - 95.63141063070621*pow(eta,5.)*pow(nu,2.);
469 shape(210,0) = -2.995357736356377*eta + 26.958219627207395*pow(eta,3.)
470 - 59.30808317985627*pow(eta,5.)
471 + 36.714527682768164*pow(eta,7.);
472 shape(211,1) = -2.995357736356377*eta + 26.958219627207395*pow(eta,3.)
473 - 59.30808317985627*pow(eta,5.)
474 + 36.714527682768164*pow(eta,7.);
475 shape(211,2) = 2.995357736356377*nu
476 - 80.87465888162218*pow(eta,2.)*nu
477 + 296.54041589928136*pow(eta,4.)*nu
478 - 257.00169377937715*pow(eta,6.)*nu;
479 shape(212,2) = -2.995357736356377*eta + 26.958219627207395*pow(eta,3.)
480 - 59.30808317985627*pow(eta,5.)
481 + 36.714527682768164*pow(eta,7.);
482 shape(213,0) = -0.689981317681863*xi
483 + 14.489607671319124*xi*pow(nu,2.)
484 - 43.46882301395737*xi*pow(nu,4.)
485 + 31.87713687690207*xi*pow(nu,6.);
486 shape(213,2) = 0.689981317681863*nu - 4.829869223773041*pow(nu,3.)
487 + 8.693764602791473*pow(nu,5.)
488 - 4.553876696700296*pow(nu,7.);
489 shape(214,1) = -0.689981317681863*xi
490 + 14.489607671319124*xi*pow(nu,2.)
491 - 43.46882301395737*xi*pow(nu,4.)
492 + 31.87713687690207*xi*pow(nu,6.);
493 shape(215,0) = 6.595897162251698*xi*eta*nu
494 - 30.780853423841258*xi*eta*pow(nu,3.)
495 + 27.70276808145713*xi*eta*pow(nu,5.);
496 shape(215,2) = 0.21986323874172325*eta
497 - 3.297948581125849*eta*pow(nu,2.)
498 + 7.695213355960314*eta*pow(nu,4.)
499 - 4.617128013576188*eta*pow(nu,6.);
500 shape(216,1) = 6.595897162251698*xi*eta*nu
501 - 30.780853423841258*xi*eta*pow(nu,3.)
502 + 27.70276808145713*xi*eta*pow(nu,5.);
503 shape(216,2) = 0.21986323874172325*xi
504 - 3.297948581125849*xi*pow(nu,2.)
505 + 7.695213355960314*xi*pow(nu,4.)
506 - 4.617128013576188*xi*pow(nu,6.);
507 shape(217,0) = -0.7702348464916399*xi
508 + 7.7023484649163985*xi*pow(nu,2.)
509 - 8.986073209069131*xi*pow(nu,4.)
510 + 2.3107045394749197*xi*pow(eta,2.)
511 - 23.107045394749196*xi*pow(eta,2.)*pow(nu,2.)
512 + 26.958219627207395*xi*pow(eta,2.)*pow(nu,4.);
513 shape(217,2) = 0.7702348464916399*nu - 2.567449488305466*pow(nu,3.)
514 + 1.7972146418138264*pow(nu,5.)
515 - 2.3107045394749197*pow(eta,2.)*nu
516 + 7.7023484649163985*pow(eta,2.)*pow(nu,3.)
517 - 5.391643925441479*pow(eta,2.)*pow(nu,5.);
518 shape(218,1) = -0.7702348464916399*xi
519 + 7.7023484649163985*xi*pow(nu,2.)
520 - 8.986073209069131*xi*pow(nu,4.)
521 + 2.3107045394749197*xi*pow(eta,2.)
522 - 23.107045394749196*xi*pow(eta,2.)*pow(nu,2.)
523 + 26.958219627207395*xi*pow(eta,2.)*pow(nu,4.);
524 shape(218,2) = -4.621409078949839*xi*eta*nu
525 + 15.404696929832797*xi*eta*pow(nu,3.)
526 - 10.783287850882958*xi*eta*pow(nu,5.);
527 shape(219,0) = 9.644865862208764*xi*eta*nu
528 - 16.074776437014606*xi*eta*pow(nu,3.)
529 - 16.074776437014606*xi*pow(eta,3.)*nu
530 + 26.79129406169101*xi*pow(eta,3.)*pow(nu,3.);
531 shape(219,2) = 0.8037388218507303*eta
532 - 4.822432931104382*eta*pow(nu,2.)
533 + 4.0186941092536514*eta*pow(nu,4.)
534 - 1.3395647030845506*pow(eta,3.)
535 + 8.037388218507303*pow(eta,3.)*pow(nu,2.)
536 - 6.697823515422753*pow(eta,3.)*pow(nu,4.);
537 shape(220,1) = 9.644865862208764*xi*eta*nu
538 - 16.074776437014606*xi*eta*pow(nu,3.)
539 - 16.074776437014606*xi*pow(eta,3.)*nu
540 + 26.79129406169101*xi*pow(eta,3.)*pow(nu,3.);
541 shape(220,2) = 0.8037388218507303*xi
542 - 4.822432931104382*xi*pow(nu,2.)
543 + 4.0186941092536514*xi*pow(nu,4.)
544 - 4.0186941092536514*xi*pow(eta,2.)
545 + 24.11216465552191*xi*pow(eta,2.)*pow(nu,2.)
546 - 20.093470546268257*xi*pow(eta,2.)*pow(nu,4.);
547 shape(221,0) = -0.7702348464916399*xi
548 + 2.3107045394749197*xi*pow(nu,2.)
549 + 7.7023484649163985*xi*pow(eta,2.)
550 - 23.107045394749196*xi*pow(eta,2.)*pow(nu,2.)
551 - 8.986073209069131*xi*pow(eta,4.)
552 + 26.958219627207395*xi*pow(eta,4.)*pow(nu,2.);
553 shape(221,2) = 0.7702348464916399*nu - 0.7702348464916399*pow(nu,3.)
554 - 7.7023484649163985*pow(eta,2.)*nu
555 + 7.7023484649163985*pow(eta,2.)*pow(nu,3.)
556 + 8.986073209069131*pow(eta,4.)*nu
557 - 8.986073209069131*pow(eta,4.)*pow(nu,3.);
558 shape(222,1) = -0.7702348464916399*xi
559 + 2.3107045394749197*xi*pow(nu,2.)
560 + 7.7023484649163985*xi*pow(eta,2.)
561 - 23.107045394749196*xi*pow(eta,2.)*pow(nu,2.)
562 - 8.986073209069131*xi*pow(eta,4.)
563 + 26.958219627207395*xi*pow(eta,4.)*pow(nu,2.);
564 shape(222,2) = -15.404696929832797*xi*eta*nu
565 + 15.404696929832797*xi*eta*pow(nu,3.)
566 + 35.944292836276524*xi*pow(eta,3.)*nu
567 - 35.944292836276524*xi*pow(eta,3.)*pow(nu,3.);
568 shape(223,0) = 6.595897162251698*xi*eta*nu
569 - 30.780853423841258*xi*pow(eta,3.)*nu
570 + 27.70276808145713*xi*pow(eta,5.)*nu;
571 shape(223,2) = 1.0993161937086162*eta
572 - 3.297948581125849*eta*pow(nu,2.)
573 - 5.130142237306876*pow(eta,3.)
574 + 15.390426711920629*pow(eta,3.)*pow(nu,2.)
575 + 4.617128013576188*pow(eta,5.)
576 - 13.851384040728565*pow(eta,5.)*pow(nu,2.);
577 shape(224,1) = 6.595897162251698*xi*eta*nu
578 - 30.780853423841258*xi*pow(eta,3.)*nu
579 + 27.70276808145713*xi*pow(eta,5.)*nu;
580 shape(224,2) = 1.0993161937086162*xi
581 - 3.297948581125849*xi*pow(nu,2.)
582 - 15.390426711920629*xi*pow(eta,2.)
583 + 46.17128013576188*xi*pow(eta,2.)*pow(nu,2.)
584 + 23.08564006788094*xi*pow(eta,4.)
585 - 69.25692020364282*xi*pow(eta,4.)*pow(nu,2.);
586 shape(225,0) = -0.689981317681863*xi
587 + 14.489607671319124*xi*pow(eta,2.)
588 - 43.46882301395737*xi*pow(eta,4.)
589 + 31.87713687690207*xi*pow(eta,6.);
590 shape(225,2) = 0.689981317681863*nu
591 - 14.489607671319124*pow(eta,2.)*nu
592 + 43.46882301395737*pow(eta,4.)*nu
593 - 31.87713687690207*pow(eta,6.)*nu;
594 shape(226,1) = -0.689981317681863*xi
595 + 14.489607671319124*xi*pow(eta,2.)
596 - 43.46882301395737*xi*pow(eta,4.)
597 + 31.87713687690207*xi*pow(eta,6.);
598 shape(226,2) = -28.97921534263825*xi*eta*nu
599 + 173.87529205582948*xi*pow(eta,3.)*nu
600 - 191.26282126141243*xi*pow(eta,5.)*nu;
601 shape(227,2) = -0.689981317681863*xi
602 + 14.489607671319124*xi*pow(eta,2.)
603 - 43.46882301395737*xi*pow(eta,4.)
604 + 31.87713687690207*xi*pow(eta,6.);
605 shape(228,0) = -2.4581457378987928*nu + 11.471346776861033*pow(nu,3.)
606 - 10.324212099174929*pow(nu,5.)
607 + 7.374437213696378*pow(xi,2.)*nu
608 - 34.4140403305831*pow(xi,2.)*pow(nu,3.)
609 + 30.972636297524787*pow(xi,2.)*pow(nu,5.);
610 shape(228,2) = 0.49162914757975856*xi
611 - 7.374437213696378*xi*pow(nu,2.)
612 + 17.20702016529155*xi*pow(nu,4.)
613 - 10.324212099174929*xi*pow(nu,6.);
614 shape(229,1) = 2.4581457378987928*nu + 11.471346776861033*pow(nu,3.)
615 - 10.324212099174929*pow(nu,5.)
616 + 7.374437213696378*pow(xi,2.)*nu
617 - 34.4140403305831*pow(xi,2.)*pow(nu,3.)
618 + 30.972636297524787*pow(xi,2.)*pow(nu,5.);
619 shape(230,0) = -0.7702348464916399*eta
620 + 7.7023484649163985*eta*pow(nu,2.)
621 - 8.986073209069131*eta*pow(nu,4.)
622 + 2.3107045394749197*pow(xi,2.)*eta
623 - 23.107045394749196*pow(xi,2.)*eta*pow(nu,2.)
624 + 26.958219627207395*pow(xi,2.)*eta*pow(nu,4.);
625 shape(230,2) = -4.621409078949839*xi*eta*nu
626 + 15.404696929832797*xi*eta*pow(nu,3.)
627 - 10.783287850882958*xi*eta*pow(nu,5.);
628 shape(231,1) = -0.7702348464916399*eta
629 + 7.7023484649163985*eta*pow(nu,2.)
630 - 8.986073209069131*eta*pow(nu,4.)
631 + 2.3107045394749197*pow(xi,2.)*eta
632 - 23.107045394749196*pow(xi,2.)*eta*pow(nu,2.)
633 + 26.958219627207395*pow(xi,2.)*eta*pow(nu,4.);
634 shape(231,2) = 0.7702348464916399*nu - 2.567449488305466*pow(nu,3.)
635 + 1.7972146418138264*pow(nu,5.)
636 - 2.3107045394749197*pow(xi,2.)*nu
637 + 7.7023484649163985*pow(xi,2.)*pow(nu,3.)
638 - 5.391643925441479*pow(xi,2.)*pow(nu,5.);
639 shape(232,0) = -1.753901900050285*nu + 2.9231698334171416*pow(nu,3.)
640 + 5.261705700150855*pow(eta,2.)*nu
641 - 8.769509500251425*pow(eta,2.)*pow(nu,3.)
642 + 5.261705700150855*pow(xi,2.)*nu
643 - 8.769509500251425*pow(xi,2.)*pow(nu,3.)
644 - 15.785117100452565*pow(xi,2.)*pow(eta,2.)*nu
645 + 26.308528500754274*pow(xi,2.)*pow(eta,2.)*pow(nu,3.);
646 shape(232,2) = 0.8769509500251426*xi
647 - 5.261705700150855*xi*pow(nu,2.)
648 + 4.384754750125713*xi*pow(nu,4.)
649 - 2.6308528500754274*xi*pow(eta,2.)
650 + 15.785117100452565*xi*pow(eta,2.)*pow(nu,2.)
651 - 13.154264250377137*xi*pow(eta,2.)*pow(nu,4.);
652 shape(233,1) = 1.753901900050285*nu + 2.9231698334171416*pow(nu,3.)
653 + 5.261705700150855*pow(eta,2.)*nu
654 - 8.769509500251425*pow(eta,2.)*pow(nu,3.)
655 + 5.261705700150855*pow(xi,2.)*nu
656 - 8.769509500251425*pow(xi,2.)*pow(nu,3.)
657 - 15.785117100452565*pow(xi,2.)*pow(eta,2.)*nu
658 + 26.308528500754274*pow(xi,2.)*pow(eta,2.)*pow(nu,3.);
659 shape(233,2) = 0.8769509500251426*eta
660 - 5.261705700150855*eta*pow(nu,2.)
661 + 4.384754750125713*eta*pow(nu,4.)
662 - 2.6308528500754274*pow(xi,2.)*eta
663 + 15.785117100452565*pow(xi,2.)*eta*pow(nu,2.)
664 - 13.154264250377137*pow(xi,2.)*eta*pow(nu,4.);
665 shape(234,0) = -1.753901900050285*eta
666 + 5.261705700150855*eta*pow(nu,2.)
667 + 2.9231698334171416*pow(eta,3.)
668 - 8.769509500251425*pow(eta,3.)*pow(nu,2.)
669 + 5.261705700150855*pow(xi,2.)*eta
670 - 15.785117100452565*pow(xi,2.)*eta*pow(nu,2.)
671 - 8.769509500251425*pow(xi,2.)*pow(eta,3.)
672 + 26.308528500754274*pow(xi,2.)*pow(eta,3.)*pow(nu,2.);
673 shape(234,2) = -10.52341140030171*xi*eta*nu
674 + 10.52341140030171*xi*eta*pow(nu,3.)
675 + 17.53901900050285*xi*pow(eta,3.)*nu
676 - 17.53901900050285*xi*pow(eta,3.)*pow(nu,3.);
677 shape(235,1) = -1.753901900050285*eta
678 + 5.261705700150855*eta*pow(nu,2.)
679 + 2.9231698334171416*pow(eta,3.)
680 - 8.769509500251425*pow(eta,3.)*pow(nu,2.)
681 + 5.261705700150855*pow(xi,2.)*eta
682 - 15.785117100452565*pow(xi,2.)*eta*pow(nu,2.)
683 - 8.769509500251425*pow(xi,2.)*pow(eta,3.)
684 + 26.308528500754274*pow(xi,2.)*pow(eta,3.)*pow(nu,2.);
685 shape(235,2) = 1.753901900050285*nu - 1.753901900050285*pow(nu,3.)
686 - 8.769509500251425*pow(eta,2.)*nu
687 + 8.769509500251425*pow(eta,2.)*pow(nu,3.)
688 - 5.261705700150855*pow(xi,2.)*nu
689 + 5.261705700150855*pow(xi,2.)*pow(nu,3.)
690 + 26.308528500754274*pow(xi,2.)*pow(eta,2.)*nu
691 - 26.308528500754274*pow(xi,2.)*pow(eta,2.)*pow(nu,3.);
692 shape(236,0) = -0.7702348464916399*nu
693 + 7.7023484649163985*pow(eta,2.)*nu
694 - 8.986073209069131*pow(eta,4.)*nu
695 + 2.3107045394749197*pow(xi,2.)*nu
696 - 23.107045394749196*pow(xi,2.)*pow(eta,2.)*nu
697 + 26.958219627207395*pow(xi,2.)*pow(eta,4.)*nu;
698 shape(236,2) = 0.7702348464916399*xi
699 - 2.3107045394749197*xi*pow(nu,2.)
700 - 7.7023484649163985*xi*pow(eta,2.)
701 + 23.107045394749196*xi*pow(eta,2.)*pow(nu,2.)
702 + 8.986073209069131*xi*pow(eta,4.)
703 - 26.958219627207395*xi*pow(eta,4.)*pow(nu,2.);
704 shape(237,1) = -0.7702348464916399*nu
705 + 7.7023484649163985*pow(eta,2.)*nu
706 - 8.986073209069131*pow(eta,4.)*nu
707 + 2.3107045394749197*pow(xi,2.)*nu
708 - 23.107045394749196*pow(xi,2.)*pow(eta,2.)*nu
709 + 26.958219627207395*pow(xi,2.)*pow(eta,4.)*nu;
710 shape(237,2) = 2.567449488305466*eta
711 - 7.7023484649163985*eta*pow(nu,2.)
712 - 5.990715472712754*pow(eta,3.)
713 + 17.972146418138262*pow(eta,3.)*pow(nu,2.)
714 - 7.7023484649163985*pow(xi,2.)*eta
715 + 23.107045394749196*pow(xi,2.)*eta*pow(nu,2.)
716 + 17.972146418138262*pow(xi,2.)*pow(eta,3.)
717 - 53.91643925441479*pow(xi,2.)*pow(eta,3.)*pow(nu,2.);
718 shape(238,0) = -2.4581457378987928*eta
719 + 11.471346776861033*pow(eta,3.)
720 - 10.324212099174929*pow(eta,5.)
721 + 7.374437213696378*pow(xi,2.)*eta
722 - 34.4140403305831*pow(xi,2.)*pow(eta,3.)
723 + 30.972636297524787*pow(xi,2.)*pow(eta,5.);
724 shape(238,2) = -14.748874427392757*xi*eta*nu
725 + 68.8280806611662*xi*pow(eta,3.)*nu
726 - 61.94527259504957*xi*pow(eta,5.)*nu;
727 shape(239,1) = -2.4581457378987928*eta
728 + 11.471346776861033*pow(eta,3.)
729 - 10.324212099174929*pow(eta,5.)
730 + 7.374437213696378*pow(xi,2.)*eta
731 - 34.4140403305831*pow(xi,2.)*pow(eta,3.)
732 + 30.972636297524787*pow(xi,2.)*pow(eta,5.);
733 shape(239,2) = 2.4581457378987928*nu
734 - 34.4140403305831*pow(eta,2.)*nu
735 + 51.621060495874644*pow(eta,4.)*nu
736 - 7.374437213696378*pow(xi,2.)*nu
737 + 103.24212099174929*pow(xi,2.)*pow(eta,2.)*nu
738 - 154.86318148762393*pow(xi,2.)*pow(eta,4.)*nu;
739 shape(240,2) = -2.4581457378987928*eta
740 + 11.471346776861033*pow(eta,3.)
741 - 10.324212099174929*pow(eta,5.)
742 + 7.374437213696378*pow(xi,2.)*eta
743 - 34.4140403305831*pow(xi,2.)*pow(eta,3.)
744 + 30.972636297524787*pow(xi,2.)*pow(eta,5.);
745 shape(241,0) = -1.5785117100452566*xi
746 + 15.785117100452565*xi*pow(nu,2.)
747 - 18.415969950527995*xi*pow(nu,4.)
748 + 2.6308528500754274*pow(xi,3.)
749 - 26.308528500754274*pow(xi,3.)*pow(nu,2.)
750 + 30.69328325087999*pow(xi,3.)*pow(nu,4.);
751 shape(241,2) = 1.5785117100452566*nu - 5.261705700150855*pow(nu,3.)
752 + 3.6831939901055986*pow(nu,5.)
753 - 7.8925585502262825*pow(xi,2.)*nu
754 + 26.308528500754274*pow(xi,2.)*pow(nu,3.)
755 - 18.415969950527995*pow(xi,2.)*pow(nu,5.);
756 shape(242,1) = -1.5785117100452566*xi
757 + 15.785117100452565*xi*pow(nu,2.)
758 - 18.415969950527995*xi*pow(nu,4.)
759 + 2.6308528500754274*pow(xi,3.)
760 - 26.308528500754274*pow(xi,3.)*pow(nu,2.)
761 + 30.69328325087999*pow(xi,3.)*pow(nu,4.);
762 shape(243,0) = 9.644865862208764*xi*eta*nu
763 - 16.074776437014606*xi*eta*pow(nu,3.)
764 - 16.074776437014606*pow(xi,3.)*eta*nu
765 + 26.79129406169101*pow(xi,3.)*eta*pow(nu,3.);
766 shape(243,2) = 0.8037388218507303*eta
767 - 4.822432931104382*eta*pow(nu,2.)
768 + 4.0186941092536514*eta*pow(nu,4.)
769 - 4.0186941092536514*pow(xi,2.)*eta
770 + 24.11216465552191*pow(xi,2.)*eta*pow(nu,2.)
771 - 20.093470546268257*pow(xi,2.)*eta*pow(nu,4.);
772 shape(244,1) = 9.644865862208764*xi*eta*nu
773 - 16.074776437014606*xi*eta*pow(nu,3.)
774 - 16.074776437014606*pow(xi,3.)*eta*nu
775 + 26.79129406169101*pow(xi,3.)*eta*pow(nu,3.);
776 shape(244,2) = 0.8037388218507303*xi
777 - 4.822432931104382*xi*pow(nu,2.)
778 + 4.0186941092536514*xi*pow(nu,4.)
779 - 1.3395647030845506*pow(xi,3.)
780 + 8.037388218507303*pow(xi,3.)*pow(nu,2.)
781 - 6.697823515422753*pow(xi,3.)*pow(nu,4.);
782 shape(245,0) = -1.753901900050285*xi + 5.261705700150855*xi*pow(nu,2.)
783 + 5.261705700150855*xi*pow(eta,2.)
784 - 15.785117100452565*xi*pow(eta,2.)*pow(nu,2.)
785 + 2.9231698334171416*pow(xi,3.)
786 - 8.769509500251425*pow(xi,3.)*pow(nu,2.)
787 - 8.769509500251425*pow(xi,3.)*pow(eta,2.)
788 + 26.308528500754274*pow(xi,3.)*pow(eta,2.)*pow(nu,2.);
789 shape(245,2) = 1.753901900050285*nu - 1.753901900050285*pow(nu,3.)
790 - 5.261705700150855*pow(eta,2.)*nu
791 + 5.261705700150855*pow(eta,2.)*pow(nu,3.)
792 - 8.769509500251425*pow(xi,2.)*nu
793 + 8.769509500251425*pow(xi,2.)*pow(nu,3.)
794 + 26.308528500754274*pow(xi,2.)*pow(eta,2.)*nu
795 - 26.308528500754274*pow(xi,2.)*pow(eta,2.)*pow(nu,3.);
796 shape(246,1) = -1.753901900050285*xi + 5.261705700150855*xi*pow(nu,2.)
797 + 5.261705700150855*xi*pow(eta,2.)
798 - 15.785117100452565*xi*pow(eta,2.)*pow(nu,2.)
799 + 2.9231698334171416*pow(xi,3.)
800 - 8.769509500251425*pow(xi,3.)*pow(nu,2.)
801 - 8.769509500251425*pow(xi,3.)*pow(eta,2.)
802 + 26.308528500754274*pow(xi,3.)*pow(eta,2.)*pow(nu,2.);
803 shape(246,2) = -10.52341140030171*xi*eta*nu
804 + 10.52341140030171*xi*eta*pow(nu,3.)
805 + 17.53901900050285*pow(xi,3.)*eta*nu
806 - 17.53901900050285*pow(xi,3.)*eta*pow(nu,3.);
807 shape(247,0) = 9.644865862208764*xi*eta*nu
808 - 16.074776437014606*xi*pow(eta,3.)*nu
809 - 16.074776437014606*pow(xi,3.)*eta*nu
810 + 26.79129406169101*pow(xi,3.)*pow(eta,3.)*nu;
811 shape(247,2) = 1.6074776437014606*eta
812 - 4.822432931104382*eta*pow(nu,2.)
813 - 2.679129406169101*pow(eta,3.)
814 + 8.037388218507303*pow(eta,3.)*pow(nu,2.)
815 - 8.037388218507303*pow(xi,2.)*eta
816 + 24.11216465552191*pow(xi,2.)*eta*pow(nu,2.)
817 + 13.395647030845504*pow(xi,2.)*pow(eta,3.)
818 - 40.186941092536514*pow(xi,2.)*pow(eta,3.)*pow(nu,2.);
819 shape(248,1) = 9.644865862208764*xi*eta*nu
820 - 16.074776437014606*xi*pow(eta,3.)*nu
821 - 16.074776437014606*pow(xi,3.)*eta*nu
822 + 26.79129406169101*pow(xi,3.)*pow(eta,3.)*nu;
823 shape(248,2) = 1.6074776437014606*xi
824 - 4.822432931104382*xi*pow(nu,2.)
825 - 8.037388218507303*xi*pow(eta,2.)
826 + 24.11216465552191*xi*pow(eta,2.)*pow(nu,2.)
827 - 2.679129406169101*pow(xi,3.)
828 + 8.037388218507303*pow(xi,3.)*pow(nu,2.)
829 + 13.395647030845504*pow(xi,3.)*pow(eta,2.)
830 - 40.186941092536514*pow(xi,3.)*pow(eta,2.)*pow(nu,2.);
831 shape(249,0) = -1.5785117100452566*xi
832 + 15.785117100452565*xi*pow(eta,2.)
833 - 18.415969950527995*xi*pow(eta,4.)
834 + 2.6308528500754274*pow(xi,3.)
835 - 26.308528500754274*pow(xi,3.)*pow(eta,2.)
836 + 30.69328325087999*pow(xi,3.)*pow(eta,4.);
837 shape(249,2) = 1.5785117100452566*nu
838 - 15.785117100452565*pow(eta,2.)*nu
839 + 18.415969950527995*pow(eta,4.)*nu
840 - 7.8925585502262825*pow(xi,2.)*nu
841 + 78.92558550226282*pow(xi,2.)*pow(eta,2.)*nu
842 - 92.07984975263996*pow(xi,2.)*pow(eta,4.)*nu;
843 shape(250,1) = -1.5785117100452566*xi
844 + 15.785117100452565*xi*pow(eta,2.)
845 - 18.415969950527995*xi*pow(eta,4.)
846 + 2.6308528500754274*pow(xi,3.)
847 - 26.308528500754274*pow(xi,3.)*pow(eta,2.)
848 + 30.69328325087999*pow(xi,3.)*pow(eta,4.);
849 shape(250,2) = -31.57023420090513*xi*eta*nu
850 + 73.66387980211196*xi*pow(eta,3.)*nu
851 + 52.61705700150855*pow(xi,3.)*eta*nu
852 - 122.77313300351994*pow(xi,3.)*pow(eta,3.)*nu;
853 shape(251,2) = -1.5785117100452566*xi
854 + 15.785117100452565*xi*pow(eta,2.)
855 - 18.415969950527995*xi*pow(eta,4.)
856 + 2.6308528500754274*pow(xi,3.)
857 - 26.308528500754274*pow(xi,3.)*pow(eta,2.)
858 + 30.69328325087999*pow(xi,3.)*pow(eta,4.);
859 shape(252,0) = -1.5785117100452566*nu + 2.6308528500754274*pow(nu,3.)
860 + 15.785117100452565*pow(xi,2.)*nu
861 - 26.308528500754274*pow(xi,2.)*pow(nu,3.)
862 - 18.415969950527995*pow(xi,4.)*nu
863 + 30.69328325087999*pow(xi,4.)*pow(nu,3.);
864 shape(252,2) = 2.6308528500754274*xi
865 - 15.785117100452565*xi*pow(nu,2.)
866 + 13.154264250377137*xi*pow(nu,4.)
867 - 6.138656650175998*pow(xi,3.)
868 + 36.83193990105599*pow(xi,3.)*pow(nu,2.)
869 - 30.69328325087999*pow(xi,3.)*pow(nu,4.);
870 shape(253,1) = -1.5785117100452566*nu + 2.6308528500754274*pow(nu,3.)
871 + 15.785117100452565*pow(xi,2.)*nu
872 - 26.308528500754274*pow(xi,2.)*pow(nu,3.)
873 - 18.415969950527995*pow(xi,4.)*nu
874 + 30.69328325087999*pow(xi,4.)*pow(nu,3.);
875 shape(254,0) = -0.7702348464916399*eta
876 + 2.3107045394749197*eta*pow(nu,2.)
877 + 7.7023484649163985*pow(xi,2.)*eta
878 - 23.107045394749196*pow(xi,2.)*eta*pow(nu,2.)
879 - 8.986073209069131*pow(xi,4.)*eta
880 + 26.958219627207395*pow(xi,4.)*eta*pow(nu,2.);
881 shape(254,2) = -15.404696929832797*xi*eta*nu
882 + 15.404696929832797*xi*eta*pow(nu,3.)
883 + 35.944292836276524*pow(xi,3.)*eta*nu
884 - 35.944292836276524*pow(xi,3.)*eta*pow(nu,3.);
885 shape(255,1) = -0.7702348464916399*eta
886 + 2.3107045394749197*eta*pow(nu,2.)
887 + 7.7023484649163985*pow(xi,2.)*eta
888 - 23.107045394749196*pow(xi,2.)*eta*pow(nu,2.)
889 - 8.986073209069131*pow(xi,4.)*eta
890 + 26.958219627207395*pow(xi,4.)*eta*pow(nu,2.);
891 shape(255,2) = 0.7702348464916399*nu - 0.7702348464916399*pow(nu,3.)
892 - 7.7023484649163985*pow(xi,2.)*nu
893 + 7.7023484649163985*pow(xi,2.)*pow(nu,3.)
894 + 8.986073209069131*pow(xi,4.)*nu
895 - 8.986073209069131*pow(xi,4.)*pow(nu,3.);
896 shape(256,0) = -0.7702348464916399*nu
897 + 2.3107045394749197*pow(eta,2.)*nu
898 + 7.7023484649163985*pow(xi,2.)*nu
899 - 23.107045394749196*pow(xi,2.)*pow(eta,2.)*nu
900 - 8.986073209069131*pow(xi,4.)*nu
901 + 26.958219627207395*pow(xi,4.)*pow(eta,2.)*nu;
902 shape(256,2) = 2.567449488305466*xi
903 - 7.7023484649163985*xi*pow(nu,2.)
904 - 7.7023484649163985*xi*pow(eta,2.)
905 + 23.107045394749196*xi*pow(eta,2.)*pow(nu,2.)
906 - 5.990715472712754*pow(xi,3.)
907 + 17.972146418138262*pow(xi,3.)*pow(nu,2.)
908 + 17.972146418138262*pow(xi,3.)*pow(eta,2.)
909 - 53.91643925441479*pow(xi,3.)*pow(eta,2.)*pow(nu,2.);
910 shape(257,1) = -0.7702348464916399*nu
911 + 2.3107045394749197*pow(eta,2.)*nu
912 + 7.7023484649163985*pow(xi,2.)*nu
913 - 23.107045394749196*pow(xi,2.)*pow(eta,2.)*nu
914 - 8.986073209069131*pow(xi,4.)*nu
915 + 26.958219627207395*pow(xi,4.)*pow(eta,2.)*nu;
916 shape(257,2) = 0.7702348464916399*eta
917 - 2.3107045394749197*eta*pow(nu,2.)
918 - 7.7023484649163985*pow(xi,2.)*eta
919 + 23.107045394749196*pow(xi,2.)*eta*pow(nu,2.)
920 + 8.986073209069131*pow(xi,4.)*eta
921 - 26.958219627207395*pow(xi,4.)*eta*pow(nu,2.);
922 shape(258,0) = -1.5785117100452566*eta
923 + 2.6308528500754274*pow(eta,3.)
924 + 15.785117100452565*pow(xi,2.)*eta
925 - 26.308528500754274*pow(xi,2.)*pow(eta,3.)
926 - 18.415969950527995*pow(xi,4.)*eta
927 + 30.69328325087999*pow(xi,4.)*pow(eta,3.);
928 shape(258,2) = -31.57023420090513*xi*eta*nu
929 + 52.61705700150855*xi*pow(eta,3.)*nu
930 + 73.66387980211196*pow(xi,3.)*eta*nu
931 - 122.77313300351994*pow(xi,3.)*pow(eta,3.)*nu;
932 shape(259,1) = -1.5785117100452566*eta
933 + 2.6308528500754274*pow(eta,3.)
934 + 15.785117100452565*pow(xi,2.)*eta
935 - 26.308528500754274*pow(xi,2.)*pow(eta,3.)
936 - 18.415969950527995*pow(xi,4.)*eta
937 + 30.69328325087999*pow(xi,4.)*pow(eta,3.);
938 shape(259,2) = 1.5785117100452566*nu
939 - 7.8925585502262825*pow(eta,2.)*nu
940 - 15.785117100452565*pow(xi,2.)*nu
941 + 78.92558550226282*pow(xi,2.)*pow(eta,2.)*nu
942 + 18.415969950527995*pow(xi,4.)*nu
943 - 92.07984975263996*pow(xi,4.)*pow(eta,2.)*nu;
944 shape(260,2) = -1.5785117100452566*eta
945 + 2.6308528500754274*pow(eta,3.)
946 + 15.785117100452565*pow(xi,2.)*eta
947 - 26.308528500754274*pow(xi,2.)*pow(eta,3.)
948 - 18.415969950527995*pow(xi,4.)*eta
949 + 30.69328325087999*pow(xi,4.)*pow(eta,3.);
950 shape(261,0) = -2.4581457378987928*xi
951 + 7.374437213696378*xi*pow(nu,2.)
952 + 11.471346776861033*pow(xi,3.)
953 - 34.4140403305831*pow(xi,3.)*pow(nu,2.)
954 - 10.324212099174929*pow(xi,5.)
955 + 30.972636297524787*pow(xi,5.)*pow(nu,2.);
956 shape(261,2) = 2.4581457378987928*nu - 2.4581457378987928*pow(nu,3.)
957 - 34.4140403305831*pow(xi,2.)*nu
958 + 34.4140403305831*pow(xi,2.)*pow(nu,3.)
959 + 51.621060495874644*pow(xi,4.)*nu
960 - 51.621060495874644*pow(xi,4.)*pow(nu,3.);
961 shape(262,1) = -2.4581457378987928*xi
962 + 7.374437213696378*xi*pow(nu,2.)
963 + 11.471346776861033*pow(xi,3.)
964 - 34.4140403305831*pow(xi,3.)*pow(nu,2.)
965 - 10.324212099174929*pow(xi,5.)
966 + 30.972636297524787*pow(xi,5.)*pow(nu,2.);
967 shape(263,0) = 6.595897162251698*xi*eta*nu
968 - 30.780853423841258*pow(xi,3.)*eta*nu
969 + 27.70276808145713*pow(xi,5.)*eta*nu;
970 shape(263,2) = 1.0993161937086162*eta
971 - 3.297948581125849*eta*pow(nu,2.)
972 - 15.390426711920629*pow(xi,2.)*eta
973 + 46.17128013576188*pow(xi,2.)*eta*pow(nu,2.)
974 + 23.08564006788094*pow(xi,4.)*eta
975 - 69.25692020364282*pow(xi,4.)*eta*pow(nu,2.);
976 shape(264,1) = 6.595897162251698*xi*eta*nu
977 - 30.780853423841258*pow(xi,3.)*eta*nu
978 + 27.70276808145713*pow(xi,5.)*eta*nu;
979 shape(264,2) = 1.0993161937086162*xi
980 - 3.297948581125849*xi*pow(nu,2.)
981 - 5.130142237306876*pow(xi,3.)
982 + 15.390426711920629*pow(xi,3.)*pow(nu,2.)
983 + 4.617128013576188*pow(xi,5.)
984 - 13.851384040728565*pow(xi,5.)*pow(nu,2.);
985 shape(265,0) = -2.4581457378987928*xi
986 + 7.374437213696378*xi*pow(eta,2.)
987 + 11.471346776861033*pow(xi,3.)
988 - 34.4140403305831*pow(xi,3.)*pow(eta,2.)
989 - 10.324212099174929*pow(xi,5.)
990 + 30.972636297524787*pow(xi,5.)*pow(eta,2.);
991 shape(265,2) = 2.4581457378987928*nu
992 - 7.374437213696378*pow(eta,2.)*nu
993 - 34.4140403305831*pow(xi,2.)*nu
994 + 103.24212099174929*pow(xi,2.)*pow(eta,2.)*nu
995 + 51.621060495874644*pow(xi,4.)*nu
996 - 154.86318148762393*pow(xi,4.)*pow(eta,2.)*nu;
997 shape(266,1) = -2.4581457378987928*xi
998 + 7.374437213696378*xi*pow(eta,2.)
999 + 11.471346776861033*pow(xi,3.)
1000 - 34.4140403305831*pow(xi,3.)*pow(eta,2.)
1001 - 10.324212099174929*pow(xi,5.)
1002 + 30.972636297524787*pow(xi,5.)*pow(eta,2.);
1003 shape(266,2) = -14.748874427392757*xi*eta*nu
1004 + 68.8280806611662*pow(xi,3.)*eta*nu
1005 - 61.94527259504957*pow(xi,5.)*eta*nu;
1006 shape(267,2) = -2.4581457378987928*xi
1007 + 7.374437213696378*xi*pow(eta,2.)
1008 + 11.471346776861033*pow(xi,3.)
1009 - 34.4140403305831*pow(xi,3.)*pow(eta,2.)
1010 - 10.324212099174929*pow(xi,5.)
1011 + 30.972636297524787*pow(xi,5.)*pow(eta,2.);
1012 shape(268,0) = -0.689981317681863*nu
1013 + 14.489607671319124*pow(xi,2.)*nu
1014 - 43.46882301395737*pow(xi,4.)*nu
1015 + 31.87713687690207*pow(xi,6.)*nu;
1016 shape(268,2) = 4.829869223773041*xi
1017 - 14.489607671319124*xi*pow(nu,2.)
1018 - 28.97921534263825*pow(xi,3.)
1019 + 86.93764602791474*pow(xi,3.)*pow(nu,2.)
1020 + 31.87713687690207*pow(xi,5.)
1021 - 95.63141063070621*pow(xi,5.)*pow(nu,2.);
1022 shape(269,1) = -0.689981317681863*nu
1023 + 14.489607671319124*pow(xi,2.)*nu
1024 - 43.46882301395737*pow(xi,4.)*nu
1025 + 31.87713687690207*pow(xi,6.)*nu;
1026 shape(270,0) = -0.689981317681863*eta
1027 + 14.489607671319124*pow(xi,2.)*eta
1028 - 43.46882301395737*pow(xi,4.)*eta
1029 + 31.87713687690207*pow(xi,6.)*eta;
1030 shape(270,2) = -28.97921534263825*xi*eta*nu
1031 + 173.87529205582948*pow(xi,3.)*eta*nu
1032 - 191.26282126141243*pow(xi,5.)*eta*nu;;
1033 shape(271,1) = -0.689981317681863*eta
1034 + 14.489607671319124*pow(xi,2.)*eta
1035 - 43.46882301395737*pow(xi,4.)*eta
1036 + 31.87713687690207*pow(xi,6.)*eta;
1037 shape(271,2) = 0.689981317681863*nu
1038 - 14.489607671319124*pow(xi,2.)*nu
1039 + 43.46882301395737*pow(xi,4.)*nu
1040 - 31.87713687690207*pow(xi,6.)*nu;
1041 shape(272,2) = -0.689981317681863*eta
1042 + 14.489607671319124*pow(xi,2.)*eta
1043 - 43.46882301395737*pow(xi,4.)*eta
1044 + 31.87713687690207*pow(xi,6.)*eta;
1045 shape(273,0) = -2.995357736356377*xi + 26.958219627207395*pow(xi,3.)
1046 - 59.30808317985627*pow(xi,5.)
1047 + 36.714527682768164*pow(xi,7.);
1048 shape(273,2) = 2.995357736356377*nu - 80.87465888162218*pow(xi,2.)*nu
1049 + 296.54041589928136*pow(xi,4.)*nu
1050 - 257.00169377937715*pow(xi,6.)*nu;
1051 shape(274,1) = -2.995357736356377*xi + 26.958219627207395*pow(xi,3.)
1052 - 59.30808317985627*pow(xi,5.)
1053 + 36.714527682768164*pow(xi,7.);
1054 shape(275,2) = -2.995357736356377*xi + 26.958219627207395*pow(xi,3.)
1055 - 59.30808317985627*pow(xi,5.)
1056 + 36.714527682768164*pow(xi,7.);
1057 case 6:
1058 shape(133,0) = -0.3983608994994363 + 8.365578889488162*pow(nu,2.)
1059 - 25.096736668464487*pow(nu,4.)
1060 + 18.404273556873957*pow(nu,6.);
1061 shape(134,1) = -0.3983608994994363 + 8.365578889488162*pow(nu,2.)
1062 - 25.096736668464487*pow(nu,4.)
1063 + 18.404273556873957*pow(nu,6.);
1064 shape(135,0) = 3.8081430021731064*eta*nu
1065 - 17.771334010141164*eta*pow(nu,3.)
1066 + 15.994200609127047*eta*pow(nu,5.);
1067 shape(136,1) = 3.8081430021731064*eta*nu
1068 - 17.771334010141164*eta*pow(nu,3.)
1069 + 15.994200609127047*eta*pow(nu,5.);
1070 shape(136,2) = 0.1269381000724369 - 1.9040715010865532*pow(nu,2.)
1071 + 4.442833502535291*pow(nu,4.)
1072 - 2.6657001015211743*pow(nu,6.);
1073 shape(137,0) = -0.44469529596117835 + 4.446952959611783*pow(nu,2.)
1074 - 5.188111786213748*pow(nu,4.)
1075 + 1.334085887883535*pow(eta,2.)
1076 - 13.34085887883535*pow(eta,2.)*pow(nu,2.)
1077 + 15.564335358641243*pow(eta,2.)*pow(nu,4.);
1078 shape(138,1) = -0.44469529596117835 + 4.446952959611783*pow(nu,2.)
1079 - 5.188111786213748*pow(nu,4.)
1080 + 1.334085887883535*pow(eta,2.)
1081 - 13.34085887883535*pow(eta,2.)*pow(nu,2.)
1082 + 15.564335358641243*pow(eta,2.)*pow(nu,4.);
1083 shape(138,2) = -2.66817177576707*eta*nu
1084 + 8.893905919223567*eta*pow(nu,3.)
1085 - 6.225734143456497*eta*pow(nu,5.);
1086 shape(139,0) = 5.568465901844061*eta*nu
1087 - 9.280776503073437*eta*pow(nu,3.)
1088 - 9.280776503073437*pow(eta,3.)*nu
1089 + 15.467960838455728*pow(eta,3.)*pow(nu,3.);
1090 shape(140,1) = 5.568465901844061*eta*nu
1091 - 9.280776503073437*eta*pow(nu,3.)
1092 - 9.280776503073437*pow(eta,3.)*nu
1093 + 15.467960838455728*pow(eta,3.)*pow(nu,3.);
1094 shape(140,2) = 0.4640388251536718 - 2.7842329509220307*pow(nu,2.)
1095 + 2.320194125768359*pow(nu,4.)
1096 - 2.320194125768359*pow(eta,2.)
1097 + 13.921164754610153*pow(eta,2.)*pow(nu,2.)
1098 - 11.600970628841795*pow(eta,2.)*pow(nu,4.);
1099 shape(141,0) = -0.44469529596117835 + 1.334085887883535*pow(nu,2.)
1100 + 4.446952959611783*pow(eta,2.)
1101 - 13.34085887883535*pow(eta,2.)*pow(nu,2.)
1102 - 5.188111786213748*pow(eta,4.)
1103 + 15.564335358641243*pow(eta,4.)*pow(nu,2.);
1104 shape(142,1) = -0.44469529596117835 + 1.334085887883535*pow(nu,2.)
1105 + 4.446952959611783*pow(eta,2.)
1106 - 13.34085887883535*pow(eta,2.)*pow(nu,2.)
1107 - 5.188111786213748*pow(eta,4.)
1108 + 15.564335358641243*pow(eta,4.)*pow(nu,2.);
1109 shape(142,2) = -8.893905919223567*eta*nu
1110 + 8.893905919223567*eta*pow(nu,3.)
1111 + 20.75244714485499*pow(eta,3.)*nu
1112 - 20.75244714485499*pow(eta,3.)*pow(nu,3.);
1113 shape(143,0) = 3.8081430021731064*eta*nu
1114 - 17.771334010141164*pow(eta,3.)*nu
1115 + 15.994200609127047*pow(eta,5.)*nu;
1116 shape(144,1) = 3.8081430021731064*eta*nu
1117 - 17.771334010141164*pow(eta,3.)*nu
1118 + 15.994200609127047*pow(eta,5.)*nu;
1119 shape(144,2) = 0.6346905003621844 - 1.9040715010865532*pow(nu,2.)
1120 - 8.885667005070582*pow(eta,2.)
1121 + 26.657001015211744*pow(eta,2.)*pow(nu,2.)
1122 + 13.328500507605872*pow(eta,4.)
1123 - 39.985501522817614*pow(eta,4.)*pow(nu,2.);
1124 shape(145,0) = -0.3983608994994363 + 8.365578889488162*pow(eta,2.)
1125 - 25.096736668464487*pow(eta,4.)
1126 + 18.404273556873957*pow(eta,6.);
1127 shape(146,1) = -0.3983608994994363 + 8.365578889488162*pow(eta,2.)
1128 - 25.096736668464487*pow(eta,4.)
1129 + 18.404273556873957*pow(eta,6.);
1130 shape(146,2) = -16.731157778976325*eta*nu
1131 + 100.38694667385795*pow(eta,3.)*nu
1132 - 110.42564134124375*pow(eta,5.)*nu;
1133 shape(147,2) = -0.3983608994994363 + 8.365578889488162*pow(eta,2.)
1134 - 25.096736668464487*pow(eta,4.)
1135 + 18.404273556873957*pow(eta,6.);
1136 shape(148,0) = 3.8081430021731064*xi*nu
1137 - 17.771334010141164*xi*pow(nu,3.)
1138 + 15.994200609127047*xi*pow(nu,5.);
1139 shape(148,2) = 0.1269381000724369 - 1.9040715010865532*pow(nu,2.)
1140 + 4.442833502535291*pow(nu,4.)
1141 - 2.6657001015211743*pow(nu,6.);
1142 shape(149,1) = 3.8081430021731064*xi*nu
1143 - 17.771334010141164*xi*pow(nu,3.)
1144 + 15.994200609127047*xi*pow(nu,5.);
1145 shape(150,0) = 1.1932426932522988*xi*eta
1146 - 11.932426932522988*xi*eta*pow(nu,2.)
1147 + 13.921164754610153*xi*eta*pow(nu,4.);
1148 shape(150,2) = -1.1932426932522988*eta*nu
1149 + 3.97747564417433*eta*pow(nu,3.)
1150 - 2.7842329509220307*eta*pow(nu,5.);
1151 shape(151,1) = 1.1932426932522988*xi*eta
1152 - 11.932426932522988*xi*eta*pow(nu,2.)
1153 + 13.921164754610153*xi*eta*pow(nu,4.);
1154 shape(151,2) = -1.1932426932522988*xi*nu
1155 + 3.97747564417433*xi*pow(nu,3.)
1156 - 2.7842329509220307*xi*pow(nu,5.);
1157 shape(152,0) = 2.7171331399105196*xi*nu
1158 - 4.5285552331842*xi*pow(nu,3.)
1159 - 8.15139941973156*xi*pow(eta,2.)*nu
1160 + 13.5856656995526*xi*pow(eta,2.)*pow(nu,3.);
1161 shape(152,2) = 0.22642776165920997 - 1.3585665699552598*pow(nu,2.)
1162 + 1.13213880829605*pow(nu,4.)
1163 - 0.6792832849776299*pow(eta,2.)
1164 + 4.07569970986578*pow(eta,2.)*pow(nu,2.)
1165 - 3.39641642488815*pow(eta,2.)*pow(nu,4.);
1166 shape(153,1) = 2.7171331399105196*xi*nu
1167 - 4.5285552331842*xi*pow(nu,3.)
1168 - 8.15139941973156*xi*pow(eta,2.)*nu
1169 + 13.5856656995526*xi*pow(eta,2.)*pow(nu,3.);
1170 shape(153,2) = -1.3585665699552598*xi*eta
1171 + 8.15139941973156*xi*eta*pow(nu,2.)
1172 - 6.792832849776299*xi*eta*pow(nu,4.);
1173 shape(154,0) = 2.7171331399105196*xi*eta
1174 - 8.15139941973156*xi*eta*pow(nu,2.)
1175 - 4.5285552331842*xi*pow(eta,3.)
1176 + 13.5856656995526*xi*pow(eta,3.)*pow(nu,2.);
1177 shape(154,2) = -2.7171331399105196*eta*nu
1178 + 2.7171331399105196*eta*pow(nu,3.)
1179 + 4.5285552331842*pow(eta,3.)*nu
1180 - 4.5285552331842*pow(eta,3.)*pow(nu,3.);
1181 shape(155,1) = 2.7171331399105196*xi*eta
1182 - 8.15139941973156*xi*eta*pow(nu,2.)
1183 - 4.5285552331842*xi*pow(eta,3.)
1184 + 13.5856656995526*xi*pow(eta,3.)*pow(nu,2.);
1185 shape(155,2) = -2.7171331399105196*xi*nu
1186 + 2.7171331399105196*xi*pow(nu,3.)
1187 + 13.5856656995526*xi*pow(eta,2.)*nu
1188 - 13.5856656995526*xi*pow(eta,2.)*pow(nu,3.);
1189 shape(156,0) = 1.1932426932522988*xi*nu
1190 - 11.932426932522988*xi*pow(eta,2.)*nu
1191 + 13.921164754610153*xi*pow(eta,4.)*nu;
1192 shape(156,2) = 0.1988737822087165 - 0.5966213466261495*pow(nu,2.)
1193 - 1.988737822087165*pow(eta,2.)
1194 + 5.966213466261495*pow(eta,2.)*pow(nu,2.)
1195 + 2.320194125768359*pow(eta,4.)
1196 - 6.9605823773050775*pow(eta,4.)*pow(nu,2.);
1197 shape(157,1) = 1.1932426932522988*xi*nu
1198 - 11.932426932522988*xi*pow(eta,2.)*nu
1199 + 13.921164754610153*xi*pow(eta,4.)*nu;
1200 shape(157,2) = -3.97747564417433*xi*eta
1201 + 11.932426932522988*xi*eta*pow(nu,2.)
1202 + 9.280776503073437*xi*pow(eta,3.)
1203 - 27.84232950922031*xi*pow(eta,3.)*pow(nu,2.);
1204 shape(158,0) = 3.8081430021731064*xi*eta
1205 - 17.771334010141164*xi*pow(eta,3.)
1206 + 15.994200609127047*xi*pow(eta,5.);
1207 shape(158,2) = -3.8081430021731064*eta*nu
1208 + 17.771334010141164*pow(eta,3.)*nu
1209 - 15.994200609127047*pow(eta,5.)*nu;
1210 shape(159,1) = 3.8081430021731064*xi*eta
1211 - 17.771334010141164*xi*pow(eta,3.)
1212 + 15.994200609127047*xi*pow(eta,5.);
1213 shape(159,2) = -3.8081430021731064*xi*nu
1214 + 53.31400203042349*xi*pow(eta,2.)*nu
1215 - 79.97100304563523*xi*pow(eta,4.)*nu;
1216 shape(160,2) = 3.8081430021731064*xi*eta
1217 - 17.771334010141164*xi*pow(eta,3.)
1218 + 15.994200609127047*xi*pow(eta,5.);
1219 shape(161,0) = -0.44469529596117835 + 4.446952959611783*pow(nu,2.)
1220 - 5.188111786213748*pow(nu,4.)
1221 + 1.334085887883535*pow(xi,2.)
1222 - 13.34085887883535*pow(xi,2.)*pow(nu,2.)
1223 + 15.564335358641243*pow(xi,2.)*pow(nu,4.);
1224 shape(161,2) = -2.66817177576707*xi*nu
1225 + 8.893905919223567*xi*pow(nu,3.)
1226 - 6.225734143456497*xi*pow(nu,5.);
1227 shape(162,1) = -0.44469529596117835 + 4.446952959611783*pow(nu,2.)
1228 - 5.188111786213748*pow(nu,4.)
1229 + 1.334085887883535*pow(xi,2.)
1230 - 13.34085887883535*pow(xi,2.)*pow(nu,2.)
1231 + 15.564335358641243*pow(xi,2.)*pow(nu,4.);
1232 shape(163,0) = 2.7171331399105196*eta*nu
1233 - 4.5285552331842*eta*pow(nu,3.)
1234 - 8.15139941973156*pow(xi,2.)*eta*nu
1235 + 13.5856656995526*pow(xi,2.)*eta*pow(nu,3.);
1236 shape(163,2) = -1.3585665699552598*xi*eta
1237 + 8.15139941973156*xi*eta*pow(nu,2.)
1238 - 6.792832849776299*xi*eta*pow(nu,4.);
1239 shape(164,1) = 2.7171331399105196*eta*nu
1240 - 4.5285552331842*eta*pow(nu,3.)
1241 - 8.15139941973156*pow(xi,2.)*eta*nu
1242 + 13.5856656995526*pow(xi,2.)*eta*pow(nu,3.);
1243 shape(164,2) = 0.22642776165920997 - 1.3585665699552598*pow(nu,2.)
1244 + 1.13213880829605*pow(nu,4.)
1245 - 0.6792832849776299*pow(xi,2.)
1246 + 4.07569970986578*pow(xi,2.)*pow(nu,2.)
1247 - 3.39641642488815*pow(xi,2.)*pow(nu,4.);
1248 shape(165,0) = -0.49410588440130926 + 1.4823176532039277*pow(nu,2.)
1249 + 1.4823176532039277*pow(eta,2.)
1250 - 4.446952959611783*pow(eta,2.)*pow(nu,2.)
1251 + 1.4823176532039277*pow(xi,2.)
1252 - 4.446952959611783*pow(xi,2.)*pow(nu,2.)
1253 - 4.446952959611783*pow(xi,2.)*pow(eta,2.)
1254 + 13.34085887883535*pow(xi,2.)*pow(eta,2.)*pow(nu,2.);
1255 shape(165,2) = -2.9646353064078554*xi*nu
1256 + 2.9646353064078554*xi*pow(nu,3.)
1257 + 8.893905919223567*xi*pow(eta,2.)*nu
1258 - 8.893905919223567*xi*pow(eta,2.)*pow(nu,3.);
1259 shape(166,1) = -0.49410588440130926 + 1.4823176532039277*pow(nu,2.)
1260 + 1.4823176532039277*pow(eta,2.)
1261 - 4.446952959611783*pow(eta,2.)*pow(nu,2.)
1262 + 1.4823176532039277*pow(xi,2.)
1263 - 4.446952959611783*pow(xi,2.)*pow(nu,2.)
1264 - 4.446952959611783*pow(xi,2.)*pow(eta,2.)
1265 + 13.34085887883535*pow(xi,2.)*pow(eta,2.)*pow(nu,2.);
1266 shape(166,2) = -2.9646353064078554*eta*nu
1267 + 2.9646353064078554*eta*pow(nu,3.)
1268 + 8.893905919223567*pow(xi,2.)*eta*nu
1269 - 8.893905919223567*pow(xi,2.)*eta*pow(nu,3.);
1270 shape(167,0) = 2.7171331399105196*eta*nu
1271 - 4.5285552331842*pow(eta,3.)*nu
1272 - 8.15139941973156*pow(xi,2.)*eta*nu
1273 + 13.5856656995526*pow(xi,2.)*pow(eta,3.)*nu;
1274 shape(167,2) = -2.7171331399105196*xi*eta
1275 + 8.15139941973156*xi*eta*pow(nu,2.)
1276 + 4.5285552331842*xi*pow(eta,3.)
1277 - 13.5856656995526*xi*pow(eta,3.)*pow(nu,2.);
1278 shape(168,1) = 2.7171331399105196*eta*nu
1279 - 4.5285552331842*pow(eta,3.)*nu
1280 - 8.15139941973156*pow(xi,2.)*eta*nu
1281 + 13.5856656995526*pow(xi,2.)*pow(eta,3.)*nu;
1282 shape(168,2) = 0.45285552331841994 - 1.3585665699552598*pow(nu,2.)
1283 - 2.2642776165921*pow(eta,2.)
1284 + 6.792832849776299*pow(eta,2.)*pow(nu,2.)
1285 - 1.3585665699552598*pow(xi,2.)
1286 + 4.07569970986578*pow(xi,2.)*pow(nu,2.)
1287 + 6.792832849776299*pow(xi,2.)*pow(eta,2.)
1288 - 20.3784985493289*pow(xi,2.)*pow(eta,2.)*pow(nu,2.);
1289 shape(169,0) = -0.44469529596117835 + 4.446952959611783*pow(eta,2.)
1290 - 5.188111786213748*pow(eta,4.)
1291 + 1.334085887883535*pow(xi,2.)
1292 - 13.34085887883535*pow(xi,2.)*pow(eta,2.)
1293 + 15.564335358641243*pow(xi,2.)*pow(eta,4.);
1294 shape(169,2) = -2.66817177576707*xi*nu
1295 + 26.6817177576707*xi*pow(eta,2.)*nu
1296 - 31.128670717282485*xi*pow(eta,4.)*nu;
1297 shape(170,1) = -0.44469529596117835 + 4.446952959611783*pow(eta,2.)
1298 - 5.188111786213748*pow(eta,4.)
1299 + 1.334085887883535*pow(xi,2.)
1300 - 13.34085887883535*pow(xi,2.)*pow(eta,2.)
1301 + 15.564335358641243*pow(xi,2.)*pow(eta,4.);
1302 shape(170,2) = -8.893905919223567*eta*nu
1303 + 20.75244714485499*pow(eta,3.)*nu
1304 + 26.6817177576707*pow(xi,2.)*eta*nu
1305 - 62.25734143456497*pow(xi,2.)*pow(eta,3.)*nu;
1306 shape(171,2) = -0.44469529596117835 + 4.446952959611783*pow(eta,2.)
1307 - 5.188111786213748*pow(eta,4.)
1308 + 1.334085887883535*pow(xi,2.)
1309 - 13.34085887883535*pow(xi,2.)*pow(eta,2.)
1310 + 15.564335358641243*pow(xi,2.)*pow(eta,4.);
1311 shape(172,0) = 5.568465901844061*xi*nu
1312 - 9.280776503073437*xi*pow(nu,3.)
1313 - 9.280776503073437*pow(xi,3.)*nu
1314 + 15.467960838455728*pow(xi,3.)*pow(nu,3.);
1315 shape(172,2) = 0.4640388251536718 - 2.7842329509220307*pow(nu,2.)
1316 + 2.320194125768359*pow(nu,4.)
1317 - 2.320194125768359*pow(xi,2.)
1318 + 13.921164754610153*pow(xi,2.)*pow(nu,2.)
1319 - 11.600970628841795*pow(xi,2.)*pow(nu,4.);
1320 shape(173,1) = 5.568465901844061*xi*nu
1321 - 9.280776503073437*xi*pow(nu,3.)
1322 - 9.280776503073437*pow(xi,3.)*nu
1323 + 15.467960838455728*pow(xi,3.)*pow(nu,3.);
1324 shape(174,0) = 2.7171331399105196*xi*eta
1325 - 8.15139941973156*xi*eta*pow(nu,2.)
1326 - 4.5285552331842*pow(xi,3.)*eta
1327 + 13.5856656995526*pow(xi,3.)*eta*pow(nu,2.);
1328 shape(174,2) = -2.7171331399105196*eta*nu
1329 + 2.7171331399105196*eta*pow(nu,3.)
1330 + 13.5856656995526*pow(xi,2.)*eta*nu
1331 - 13.5856656995526*pow(xi,2.)*eta*pow(nu,3.);
1332 shape(175,1) = 2.7171331399105196*xi*eta
1333 - 8.15139941973156*xi*eta*pow(nu,2.)
1334 - 4.5285552331842*pow(xi,3.)*eta
1335 + 13.5856656995526*pow(xi,3.)*eta*pow(nu,2.);
1336 shape(175,2) = -2.7171331399105196*xi*nu
1337 + 2.7171331399105196*xi*pow(nu,3.)
1338 + 4.5285552331842*pow(xi,3.)*nu
1339 - 4.5285552331842*pow(xi,3.)*pow(nu,3.);
1340 shape(176,0) = 2.7171331399105196*xi*nu
1341 - 8.15139941973156*xi*pow(eta,2.)*nu
1342 - 4.5285552331842*pow(xi,3.)*nu
1343 + 13.5856656995526*pow(xi,3.)*pow(eta,2.)*nu;
1344 shape(176,2) = 0.45285552331841994 - 1.3585665699552598*pow(nu,2.)
1345 - 1.3585665699552598*pow(eta,2.)
1346 + 4.07569970986578*pow(eta,2.)*pow(nu,2.)
1347 - 2.2642776165921*pow(xi,2.)
1348 + 6.792832849776299*pow(xi,2.)*pow(nu,2.)
1349 + 6.792832849776299*pow(xi,2.)*pow(eta,2.)
1350 - 20.3784985493289*pow(xi,2.)*pow(eta,2.)*pow(nu,2.);
1351 shape(177,1) = 2.7171331399105196*xi*nu
1352 - 8.15139941973156*xi*pow(eta,2.)*nu
1353 - 4.5285552331842*pow(xi,3.)*nu
1354 + 13.5856656995526*pow(xi,3.)*pow(eta,2.)*nu;
1355 shape(177,2) = -2.7171331399105196*xi*eta
1356 + 8.15139941973156*xi*eta*pow(nu,2.)
1357 + 4.5285552331842*pow(xi,3.)*eta
1358 - 13.5856656995526*pow(xi,3.)*eta*pow(nu,2.);
1359 shape(178,0) = 5.568465901844061*xi*eta
1360 - 9.280776503073437*xi*pow(eta,3.)
1361 - 9.280776503073437*pow(xi,3.)*eta
1362 + 15.467960838455728*pow(xi,3.)*pow(eta,3.);
1363 shape(178,2) = -5.568465901844061*eta*nu
1364 + 9.280776503073437*pow(eta,3.)*nu
1365 + 27.84232950922031*pow(xi,2.)*eta*nu
1366 - 46.40388251536718*pow(xi,2.)*pow(eta,3.)*nu;
1367 shape(179,1) = 5.568465901844061*xi*eta
1368 - 9.280776503073437*xi*pow(eta,3.)
1369 - 9.280776503073437*pow(xi,3.)*eta
1370 + 15.467960838455728*pow(xi,3.)*pow(eta,3.);
1371 shape(179,2) = -5.568465901844061*xi*nu
1372 + 27.84232950922031*xi*pow(eta,2.)*nu
1373 + 9.280776503073437*pow(xi,3.)*nu
1374 - 46.40388251536718*pow(xi,3.)*pow(eta,2.)*nu;
1375 shape(180,2) = 5.568465901844061*xi*eta
1376 - 9.280776503073437*xi*pow(eta,3.)
1377 - 9.280776503073437*pow(xi,3.)*eta
1378 + 15.467960838455728*pow(xi,3.)*pow(eta,3.);
1379 shape(181,0) = -0.44469529596117835 + 1.334085887883535*pow(nu,2.)
1380 + 4.446952959611783*pow(xi,2.)
1381 - 13.34085887883535*pow(xi,2.)*pow(nu,2.)
1382 - 5.188111786213748*pow(xi,4.)
1383 + 15.564335358641243*pow(xi,4.)*pow(nu,2.);
1384 shape(181,2) = -8.893905919223567*xi*nu
1385 + 8.893905919223567*xi*pow(nu,3.)
1386 + 20.75244714485499*pow(xi,3.)*nu
1387 - 20.75244714485499*pow(xi,3.)*pow(nu,3.);
1388 shape(182,1) = -0.44469529596117835 + 1.334085887883535*pow(nu,2.)
1389 + 4.446952959611783*pow(xi,2.)
1390 - 13.34085887883535*pow(xi,2.)*pow(nu,2.)
1391 - 5.188111786213748*pow(xi,4.)
1392 + 15.564335358641243*pow(xi,4.)*pow(nu,2.);
1393 shape(183,0) = 1.1932426932522988*eta*nu
1394 - 11.932426932522988*pow(xi,2.)*eta*nu
1395 + 13.921164754610153*pow(xi,4.)*eta*nu;
1396 shape(183,2) = -3.97747564417433*xi*eta
1397 + 11.932426932522988*xi*eta*pow(nu,2.)
1398 + 9.280776503073437*pow(xi,3.)*eta
1399 - 27.84232950922031*pow(xi,3.)*eta*pow(nu,2.);
1400 shape(184,1) = 1.1932426932522988*eta*nu
1401 - 11.932426932522988*pow(xi,2.)*eta*nu
1402 + 13.921164754610153*pow(xi,4.)*eta*nu;
1403 shape(184,2) = 0.1988737822087165 - 0.5966213466261495*pow(nu,2.)
1404 - 1.988737822087165*pow(xi,2.)
1405 + 5.966213466261495*pow(xi,2.)*pow(nu,2.)
1406 + 2.320194125768359*pow(xi,4.)
1407 - 6.9605823773050775*pow(xi,4.)*pow(nu,2.);
1408 shape(185,0) = -0.44469529596117835 + 1.334085887883535*pow(eta,2.)
1409 + 4.446952959611783*pow(xi,2.)
1410 - 13.34085887883535*pow(xi,2.)*pow(eta,2.)
1411 - 5.188111786213748*pow(xi,4.)
1412 + 15.564335358641243*pow(xi,4.)*pow(eta,2.);
1413 shape(185,2) = -8.893905919223567*xi*nu
1414 + 26.6817177576707*xi*pow(eta,2.)*nu
1415 + 20.75244714485499*pow(xi,3.)*nu
1416 - 62.25734143456497*pow(xi,3.)*pow(eta,2.)*nu;
1417 shape(186,1) = -0.44469529596117835 + 1.334085887883535*pow(eta,2.)
1418 + 4.446952959611783*pow(xi,2.)
1419 - 13.34085887883535*pow(xi,2.)*pow(eta,2.)
1420 - 5.188111786213748*pow(xi,4.)
1421 + 15.564335358641243*pow(xi,4.)*pow(eta,2.);
1422 shape(186,2) = -2.66817177576707*eta*nu
1423 + 26.6817177576707*pow(xi,2.)*eta*nu
1424 - 31.128670717282485*pow(xi,4.)*eta*nu;
1425 shape(187,2) = -0.44469529596117835 + 1.334085887883535*pow(eta,2.)
1426 + 4.446952959611783*pow(xi,2.)
1427 - 13.34085887883535*pow(xi,2.)*pow(eta,2.)
1428 - 5.188111786213748*pow(xi,4.)
1429 + 15.564335358641243*pow(xi,4.)*pow(eta,2.);
1430 shape(188,0) = 3.8081430021731064*xi*nu
1431 - 17.771334010141164*pow(xi,3.)*nu
1432 + 15.994200609127047*pow(xi,5.)*nu;
1433 shape(188,2) = 0.6346905003621844 - 1.9040715010865532*pow(nu,2.)
1434 - 8.885667005070582*pow(xi,2.)
1435 + 26.657001015211744*pow(xi,2.)*pow(nu,2.)
1436 + 13.328500507605872*pow(xi,4.)
1437 - 39.985501522817614*pow(xi,4.)*pow(nu,2.);
1438 shape(189,1) = 3.8081430021731064*xi*nu
1439 - 17.771334010141164*pow(xi,3.)*nu
1440 + 15.994200609127047*pow(xi,5.)*nu;
1441 shape(190,0) = 3.8081430021731064*xi*eta
1442 - 17.771334010141164*pow(xi,3.)*eta
1443 + 15.994200609127047*pow(xi,5.)*eta;
1444 shape(190,2) = -3.8081430021731064*eta*nu
1445 + 53.31400203042349*pow(xi,2.)*eta*nu
1446 - 79.97100304563523*pow(xi,4.)*eta*nu;
1447 shape(191,1) = 3.8081430021731064*xi*eta
1448 - 17.771334010141164*pow(xi,3.)*eta
1449 + 15.994200609127047*pow(xi,5.)*eta;
1450 shape(191,2) = -3.8081430021731064*xi*nu
1451 + 17.771334010141164*pow(xi,3.)*nu
1452 - 15.994200609127047*pow(xi,5.)*nu;
1453 shape(192,2) = 3.8081430021731064*xi*eta
1454 - 17.771334010141164*pow(xi,3.)*eta
1455 + 15.994200609127047*pow(xi,5.)*eta;
1456 shape(193,0) = -0.3983608994994363 + 8.365578889488162*pow(xi,2.)
1457 - 25.096736668464487*pow(xi,4.)
1458 + 18.404273556873957*pow(xi,6.);
1459 shape(193,2) = -16.731157778976325*xi*nu
1460 + 100.38694667385795*pow(xi,3.)*nu
1461 - 110.42564134124375*pow(xi,5.)*nu;
1462 shape(194,1) = -0.3983608994994363 + 8.365578889488162*pow(xi,2.)
1463 - 25.096736668464487*pow(xi,4.)
1464 + 18.404273556873957*pow(xi,6.);
1465 shape(195,2) = -0.3983608994994363 + 8.365578889488162*pow(xi,2.)
1466 - 25.096736668464487*pow(xi,4.)
1467 + 18.404273556873957*pow(xi,6.);
1468 case 5:
1469 shape( 85,0) = 2.1986323874172324*nu - 10.260284474613751*pow(nu,3.)
1470 + 9.234256027152377*pow(nu,5.);
1471 shape( 86,1) = 2.1986323874172324*nu - 10.260284474613751*pow(nu,3.)
1472 + 9.234256027152377*pow(nu,5.);
1473 shape( 87,0) = 0.6889189901577688*eta
1474 - 6.889189901577688*eta*pow(nu,2.)
1475 + 8.037388218507303*eta*pow(nu,4.);
1476 shape( 88,1) = 0.6889189901577688*eta
1477 - 6.889189901577688*eta*pow(nu,2.)
1478 + 8.037388218507303*eta*pow(nu,4.);
1479 shape( 88,2) = -0.6889189901577688*nu + 2.2963966338592297*pow(nu,3.)
1480 - 1.6074776437014606*pow(nu,5.);
1481 shape( 89,0) = 1.5687375497513918*nu - 2.6145625829189862*pow(nu,3.)
1482 - 4.706212649254175*pow(eta,2.)*nu
1483 + 7.843687748756959*pow(eta,2.)*pow(nu,3.);
1484 shape( 90,1) = 1.5687375497513918*nu - 2.6145625829189862*pow(nu,3.)
1485 - 4.706212649254175*pow(eta,2.)*nu
1486 + 7.843687748756959*pow(eta,2.)*pow(nu,3.);
1487 shape( 90,2) = -0.7843687748756958*eta
1488 + 4.706212649254175*eta*pow(nu,2.)
1489 - 3.921843874378479*eta*pow(nu,4.);
1490 shape( 91,0) = 1.5687375497513918*eta
1491 - 4.706212649254175*eta*pow(nu,2.)
1492 - 2.6145625829189862*pow(eta,3.)
1493 + 7.843687748756959*pow(eta,3.)*pow(nu,2.);
1494 shape( 92,1) = 1.5687375497513918*eta
1495 - 4.706212649254175*eta*pow(nu,2.)
1496 - 2.6145625829189862*pow(eta,3.)
1497 + 7.843687748756959*pow(eta,3.)*pow(nu,2.);
1498 shape( 92,2) = -1.5687375497513918*nu + 1.5687375497513918*pow(nu,3.)
1499 + 7.843687748756959*pow(eta,2.)*nu
1500 - 7.843687748756959*pow(eta,2.)*pow(nu,3.);
1501 shape( 93,0) = 0.6889189901577688*nu
1502 - 6.889189901577688*pow(eta,2.)*nu
1503 + 8.037388218507303*pow(eta,4.)*nu;
1504 shape( 94,1) = 0.6889189901577688*nu
1505 - 6.889189901577688*pow(eta,2.)*nu
1506 + 8.037388218507303*pow(eta,4.)*nu;
1507 shape( 94,2) = -2.2963966338592297*eta
1508 + 6.889189901577688*eta*pow(nu,2.)
1509 + 5.358258812338202*pow(eta,3.)
1510 - 16.074776437014606*pow(eta,3.)*pow(nu,2.);
1511 shape( 95,0) = 2.1986323874172324*eta
1512 - 10.260284474613751*pow(eta,3.)
1513 + 9.234256027152377*pow(eta,5.);
1514 shape( 96,1) = 2.1986323874172324*eta
1515 - 10.260284474613751*pow(eta,3.)
1516 + 9.234256027152377*pow(eta,5.);
1517 shape( 96,2) = -2.1986323874172324*nu
1518 + 30.780853423841258*pow(eta,2.)*nu
1519 - 46.17128013576188*pow(eta,4.)*nu;
1520 shape( 97,2) = 2.1986323874172324*eta
1521 - 10.260284474613751*pow(eta,3.)
1522 + 9.234256027152377*pow(eta,5.);
1523 shape( 98,0) = 0.6889189901577688*xi
1524 - 6.889189901577688*xi*pow(nu,2.)
1525 + 8.037388218507303*xi*pow(nu,4.);
1526 shape( 98,2) = -0.6889189901577688*nu + 2.2963966338592297*pow(nu,3.)
1527 - 1.6074776437014606*pow(nu,5.);
1528 shape( 99,1) = 0.6889189901577688*xi
1529 - 6.889189901577688*xi*pow(nu,2.)
1530 + 8.037388218507303*xi*pow(nu,4.);
1531 shape(100,0) = -4.209364560120684*xi*eta*nu
1532 + 7.0156076002011405*xi*eta*pow(nu,3.);
1533 shape(100,2) = -0.350780380010057*eta
1534 + 2.104682280060342*eta*pow(nu,2.)
1535 - 1.753901900050285*eta*pow(nu,4.);
1536 shape(101,1) = -4.209364560120684*xi*eta*nu
1537 + 7.0156076002011405*xi*eta*pow(nu,3.);
1538 shape(101,2) = -0.350780380010057*xi + 2.104682280060342*xi*pow(nu,2.)
1539 - 1.753901900050285*xi*pow(nu,4.);
1540 shape(102,0) = 0.7654655446197431*xi
1541 - 2.2963966338592297*xi*pow(nu,2.)
1542 - 2.2963966338592297*xi*pow(eta,2.)
1543 + 6.889189901577688*xi*pow(eta,2.)*pow(nu,2.);
1544 shape(102,2) = -0.7654655446197431*nu + 0.7654655446197431*pow(nu,3.)
1545 + 2.2963966338592297*pow(eta,2.)*nu
1546 - 2.2963966338592297*pow(eta,2.)*pow(nu,3.);
1547 shape(103,1) = 0.7654655446197431*xi
1548 - 2.2963966338592297*xi*pow(nu,2.)
1549 - 2.2963966338592297*xi*pow(eta,2.)
1550 + 6.889189901577688*xi*pow(eta,2.)*pow(nu,2.);
1551 shape(103,2) = 4.592793267718459*xi*eta*nu
1552 - 4.592793267718459*xi*eta*pow(nu,3.);
1553 shape(104,0) = -4.209364560120684*xi*eta*nu
1554 + 7.0156076002011405*xi*pow(eta,3.)*nu;
1555 shape(104,2) = -0.701560760020114*eta
1556 + 2.104682280060342*eta*pow(nu,2.)
1557 + 1.1692679333668567*pow(eta,3.)
1558 - 3.50780380010057*pow(eta,3.)*pow(nu,2.);
1559 shape(105,1) = -4.209364560120684*xi*eta*nu
1560 + 7.0156076002011405*xi*pow(eta,3.)*nu;
1561 shape(105,2) = -0.701560760020114*xi + 2.104682280060342*xi*pow(nu,2.)
1562 + 3.50780380010057*xi*pow(eta,2.)
1563 - 10.52341140030171*xi*pow(eta,2.)*pow(nu,2.);
1564 shape(106,0) = 0.6889189901577688*xi
1565 - 6.889189901577688*xi*pow(eta,2.)
1566 + 8.037388218507303*xi*pow(eta,4.);
1567 shape(106,2) = -0.6889189901577688*nu
1568 + 6.889189901577688*pow(eta,2.)*nu
1569 - 8.037388218507303*pow(eta,4.)*nu;
1570 shape(107,1) = 0.6889189901577688*xi
1571 - 6.889189901577688*xi*pow(eta,2.)
1572 + 8.037388218507303*xi*pow(eta,4.);
1573 shape(107,2) = 13.778379803155376*xi*eta*nu
1574 - 32.14955287402921*xi*pow(eta,3.)*nu;
1575 shape(108,2) = 0.6889189901577688*xi
1576 - 6.889189901577688*xi*pow(eta,2.)
1577 + 8.037388218507303*xi*pow(eta,4.);
1578 shape(109,0) = 1.5687375497513918*nu - 2.6145625829189862*pow(nu,3.)
1579 - 4.706212649254175*pow(xi,2.)*nu
1580 + 7.843687748756959*pow(xi,2.)*pow(nu,3.);
1581 shape(109,2) = -0.7843687748756958*xi
1582 + 4.706212649254175*xi*pow(nu,2.)
1583 - 3.921843874378479*xi*pow(nu,4.);
1584 shape(110,1) = 1.5687375497513918*nu - 2.6145625829189862*pow(nu,3.)
1585 - 4.706212649254175*pow(xi,2.)*nu
1586 + 7.843687748756959*pow(xi,2.)*pow(nu,3.);
1587 shape(111,0) = 0.7654655446197431*eta
1588 - 2.2963966338592297*eta*pow(nu,2.)
1589 - 2.2963966338592297*pow(xi,2.)*eta
1590 + 6.889189901577688*pow(xi,2.)*eta*pow(nu,2.);
1591 shape(111,2) = 4.592793267718459*xi*eta*nu
1592 - 4.592793267718459*xi*eta*pow(nu,3.);
1593 shape(112,1) = 0.7654655446197431*eta
1594 - 2.2963966338592297*eta*pow(nu,2.)
1595 - 2.2963966338592297*pow(xi,2.)*eta
1596 + 6.889189901577688*pow(xi,2.)*eta*pow(nu,2.);
1597 shape(112,2) = -0.7654655446197431*nu + 0.7654655446197431*pow(nu,3.)
1598 + 2.2963966338592297*pow(xi,2.)*nu
1599 - 2.2963966338592297*pow(xi,2.)*pow(nu,3.);
1600 shape(113,0) = 0.7654655446197431*nu
1601 - 2.2963966338592297*pow(eta,2.)*nu
1602 - 2.2963966338592297*pow(xi,2.)*nu
1603 + 6.889189901577688*pow(xi,2.)*pow(eta,2.)*nu;
1604 shape(113,2) = -0.7654655446197431*xi
1605 + 2.2963966338592297*xi*pow(nu,2.)
1606 + 2.2963966338592297*xi*pow(eta,2.)
1607 - 6.889189901577688*xi*pow(eta,2.)*pow(nu,2.);
1608 shape(114,1) = 0.7654655446197431*nu
1609 - 2.2963966338592297*pow(eta,2.)*nu
1610 - 2.2963966338592297*pow(xi,2.)*nu
1611 + 6.889189901577688*pow(xi,2.)*pow(eta,2.)*nu;
1612 shape(114,2) = -0.7654655446197431*eta
1613 + 2.2963966338592297*eta*pow(nu,2.)
1614 + 2.2963966338592297*pow(xi,2.)*eta
1615 - 6.889189901577688*pow(xi,2.)*eta*pow(nu,2.);
1616 shape(115,0) = 1.5687375497513918*eta
1617 - 2.6145625829189862*pow(eta,3.)
1618 - 4.706212649254175*pow(xi,2.)*eta
1619 + 7.843687748756959*pow(xi,2.)*pow(eta,3.);
1620 shape(115,2) = 9.41242529850835*xi*eta*nu
1621 - 15.687375497513917*xi*pow(eta,3.)*nu;
1622 shape(116,1) = 1.5687375497513918*eta
1623 - 2.6145625829189862*pow(eta,3.)
1624 - 4.706212649254175*pow(xi,2.)*eta
1625 + 7.843687748756959*pow(xi,2.)*pow(eta,3.);
1626 shape(116,2) = -1.5687375497513918*nu
1627 + 7.843687748756959*pow(eta,2.)*nu
1628 + 4.706212649254175*pow(xi,2.)*nu
1629 - 23.531063246270875*pow(xi,2.)*pow(eta,2.)*nu;
1630 shape(117,2) = 1.5687375497513918*eta
1631 - 2.6145625829189862*pow(eta,3.)
1632 - 4.706212649254175*pow(xi,2.)*eta
1633 + 7.843687748756959*pow(xi,2.)*pow(eta,3.);
1634 shape(118,0) = 1.5687375497513918*xi
1635 - 4.706212649254175*xi*pow(nu,2.)
1636 - 2.6145625829189862*pow(xi,3.)
1637 + 7.843687748756959*pow(xi,3.)*pow(nu,2.);
1638 shape(118,2) = -1.5687375497513918*nu + 1.5687375497513918*pow(nu,3.)
1639 + 7.843687748756959*pow(xi,2.)*nu
1640 - 7.843687748756959*pow(xi,2.)*pow(nu,3.);
1641 shape(119,1) = 1.5687375497513918*xi
1642 - 4.706212649254175*xi*pow(nu,2.)
1643 - 2.6145625829189862*pow(xi,3.)
1644 + 7.843687748756959*pow(xi,3.)*pow(nu,2.);
1645 shape(120,0) = -4.209364560120684*xi*eta*nu
1646 + 7.0156076002011405*pow(xi,3.)*eta*nu;
1647 shape(120,2) = -0.701560760020114*eta
1648 + 2.104682280060342*eta*pow(nu,2.)
1649 + 3.50780380010057*pow(xi,2.)*eta
1650 - 10.52341140030171*pow(xi,2.)*eta*pow(nu,2.);
1651 shape(121,1) = -4.209364560120684*xi*eta*nu
1652 + 7.0156076002011405*pow(xi,3.)*eta*nu;
1653 shape(121,2) = -0.701560760020114*xi + 2.104682280060342*xi*pow(nu,2.)
1654 + 1.1692679333668567*pow(xi,3.)
1655 - 3.50780380010057*pow(xi,3.)*pow(nu,2.);
1656 shape(122,0) = 1.5687375497513918*xi
1657 - 4.706212649254175*xi*pow(eta,2.)
1658 - 2.6145625829189862*pow(xi,3.)
1659 + 7.843687748756959*pow(xi,3.)*pow(eta,2.);
1660 shape(122,2) = -1.5687375497513918*nu
1661 + 4.706212649254175*pow(eta,2.)*nu
1662 + 7.843687748756959*pow(xi,2.)*nu
1663 - 23.531063246270875*pow(xi,2.)*pow(eta,2.)*nu;
1664 shape(123,1) = 1.5687375497513918*xi
1665 - 4.706212649254175*xi*pow(eta,2.)
1666 - 2.6145625829189862*pow(xi,3.)
1667 + 7.843687748756959*pow(xi,3.)*pow(eta,2.);
1668 shape(123,2) = 9.41242529850835*xi*eta*nu
1669 - 15.687375497513917*pow(xi,3.)*eta*nu;
1670 shape(124,2) = 1.5687375497513918*xi
1671 - 4.706212649254175*xi*pow(eta,2.)
1672 - 2.6145625829189862*pow(xi,3.)
1673 + 7.843687748756959*pow(xi,3.)*pow(eta,2.);
1674 shape(125,0) = 0.6889189901577688*nu
1675 - 6.889189901577688*pow(xi,2.)*nu
1676 + 8.037388218507303*pow(xi,4.)*nu;
1677 shape(125,2) = -2.2963966338592297*xi
1678 + 6.889189901577688*xi*pow(nu,2.)
1679 + 5.358258812338202*pow(xi,3.)
1680 - 16.074776437014606*pow(xi,3.)*pow(nu,2.);
1681 shape(126,1) = 0.6889189901577688*nu
1682 - 6.889189901577688*pow(xi,2.)*nu
1683 + 8.037388218507303*pow(xi,4.)*nu;
1684 shape(127,0) = 0.6889189901577688*eta
1685 - 6.889189901577688*pow(xi,2.)*eta
1686 + 8.037388218507303*pow(xi,4.)*eta;
1687 shape(127,2) = 13.778379803155376*xi*eta*nu
1688 - 32.14955287402921*pow(xi,3.)*eta*nu;
1689 shape(128,1) = 0.6889189901577688*eta
1690 - 6.889189901577688*pow(xi,2.)*eta
1691 + 8.037388218507303*pow(xi,4.)*eta;
1692 shape(128,2) = -0.6889189901577688*nu
1693 + 6.889189901577688*pow(xi,2.)*nu
1694 - 8.037388218507303*pow(xi,4.)*nu;
1695 shape(129,2) = 0.6889189901577688*eta
1696 - 6.889189901577688*pow(xi,2.)*eta
1697 + 8.037388218507303*pow(xi,4.)*eta;
1698 shape(130,0) = 2.1986323874172324*xi - 10.260284474613751*pow(xi,3.)
1699 + 9.234256027152377*pow(xi,5.);
1700 shape(130,2) = -2.1986323874172324*nu
1701 + 30.780853423841258*pow(xi,2.)*nu
1702 - 46.17128013576188*pow(xi,4.)*nu;
1703 shape(131,1) = 2.1986323874172324*xi - 10.260284474613751*pow(xi,3.)
1704 + 9.234256027152377*pow(xi,5.);
1705 shape(132,2) = 2.1986323874172324*xi - 10.260284474613751*pow(xi,3.)
1706 + 9.234256027152377*pow(xi,5.);
1707 case 4:
1708 shape( 50,0) = 0.397747564417433 - 3.97747564417433*pow(nu,2.)
1709 + 4.640388251536718*pow(nu,4.);
1710 shape( 51,1) = 0.397747564417433 - 3.97747564417433*pow(nu,2.)
1711 + 4.640388251536718*pow(nu,4.);
1712 shape( 52,0) = -2.4302777619029476*eta*nu
1713 + 4.050462936504912*eta*pow(nu,3.);
1714 shape( 53,1) = -2.4302777619029476*eta*nu
1715 + 4.050462936504912*eta*pow(nu,3.);
1716 shape( 53,2) = -0.20252314682524564 + 1.2151388809514738*pow(nu,2.)
1717 - 1.0126157341262283*pow(nu,4.);
1718 shape( 54,0) = 0.4419417382415922 - 1.3258252147247767*pow(nu,2.)
1719 - 1.3258252147247767*pow(eta,2.)
1720 + 3.97747564417433*pow(eta,2.)*pow(nu,2.);
1721 shape( 55,1) = 0.4419417382415922 - 1.3258252147247767*pow(nu,2.)
1722 - 1.3258252147247767*pow(eta,2.)
1723 + 3.97747564417433*pow(eta,2.)*pow(nu,2.);
1724 shape( 55,2) = 2.6516504294495533*eta*nu
1725 - 2.6516504294495533*eta*pow(nu,3.);
1726 shape( 56,0) = -2.4302777619029476*eta*nu
1727 + 4.050462936504912*pow(eta,3.)*nu;
1728 shape( 57,1) = -2.4302777619029476*eta*nu
1729 + 4.050462936504912*pow(eta,3.)*nu;
1730 shape( 57,2) = -0.4050462936504913 + 1.2151388809514738*pow(nu,2.)
1731 + 2.025231468252456*pow(eta,2.)
1732 - 6.075694404757369*pow(eta,2.)*pow(nu,2.);
1733 shape( 58,0) = 0.397747564417433 - 3.97747564417433*pow(eta,2.)
1734 + 4.640388251536718*pow(eta,4.);
1735 shape( 59,1) = 0.397747564417433 - 3.97747564417433*pow(eta,2.)
1736 + 4.640388251536718*pow(eta,4.);
1737 shape( 59,2) = 7.95495128834866*eta*nu
1738 - 18.561553006146873*pow(eta,3.)*nu;
1739 shape( 60,2) = 0.397747564417433 - 3.97747564417433*pow(eta,2.)
1740 + 4.640388251536718*pow(eta,4.);
1741 shape( 61,0) = -2.4302777619029476*xi*nu
1742 + 4.050462936504912*xi*pow(nu,3.);
1743 shape( 61,2) = -0.20252314682524564 + 1.2151388809514738*pow(nu,2.)
1744 - 1.0126157341262283*pow(nu,4.);
1745 shape( 62,1) = -2.4302777619029476*xi*nu
1746 + 4.050462936504912*xi*pow(nu,3.);
1747 shape( 63,0) = -1.1858541225631423*xi*eta
1748 + 3.557562367689427*xi*eta*pow(nu,2.);
1749 shape( 63,2) = 1.1858541225631423*eta*nu
1750 - 1.1858541225631423*eta*pow(nu,3.);
1751 shape( 64,1) = -1.1858541225631423*xi*eta
1752 + 3.557562367689427*xi*eta*pow(nu,2.);
1753 shape( 64,2) = 1.1858541225631423*xi*nu
1754 - 1.1858541225631423*xi*pow(nu,3.);
1755 shape( 65,0) = -1.1858541225631423*xi*nu
1756 + 3.557562367689427*xi*pow(eta,2.)*nu;
1757 shape( 65,2) = -0.19764235376052372
1758 + 0.5929270612815711*pow(nu,2.)
1759 + 0.5929270612815711*pow(eta,2.)
1760 - 1.7787811838447134*pow(eta,2.)*pow(nu,2.);
1761 shape( 66,1) = -1.1858541225631423*xi*nu
1762 + 3.557562367689427*xi*pow(eta,2.)*nu;
1763 shape( 66,2) = 1.1858541225631423*xi*eta
1764 - 3.557562367689427*xi*eta*pow(nu,2.);
1765 shape( 67,0) = -2.4302777619029476*xi*eta
1766 + 4.050462936504912*xi*pow(eta,3.);
1767 shape( 67,2) = 2.4302777619029476*eta*nu
1768 - 4.050462936504912*pow(eta,3.)*nu;
1769 shape( 68,1) = -2.4302777619029476*xi*eta
1770 + 4.050462936504912*xi*pow(eta,3.);
1771 shape( 68,2) = 2.4302777619029476*xi*nu
1772 - 12.151388809514739*xi*pow(eta,2.)*nu;
1773 shape( 69,2) = -2.4302777619029476*xi*eta
1774 + 4.050462936504912*xi*pow(eta,3.);
1775 shape( 70,0) = 0.4419417382415922 - 1.3258252147247767*pow(nu,2.)
1776 - 1.3258252147247767*pow(xi,2.)
1777 + 3.97747564417433*pow(xi,2.)*pow(nu,2.);
1778 shape( 70,2) = 2.6516504294495533*xi*nu
1779 - 2.6516504294495533*xi*pow(nu,3.);
1780 shape( 71,1) = 0.4419417382415922 - 1.3258252147247767*pow(nu,2.)
1781 - 1.3258252147247767*pow(xi,2.)
1782 + 3.97747564417433*pow(xi,2.)*pow(nu,2.);
1783 shape( 72,0) = -1.1858541225631423*eta*nu
1784 + 3.557562367689427*pow(xi,2.)*eta*nu;
1785 shape( 72,2) = 1.1858541225631423*xi*eta
1786 - 3.557562367689427*xi*eta*pow(nu,2.);
1787 shape( 73,1) = -1.1858541225631423*eta*nu
1788 + 3.557562367689427*pow(xi,2.)*eta*nu;
1789 shape( 73,2) = -0.19764235376052372 + 0.5929270612815711*pow(nu,2.)
1790 + 0.5929270612815711*pow(xi,2.)
1791 - 1.7787811838447134*pow(xi,2.)*pow(nu,2.);
1792 shape( 74,0) = 0.4419417382415922 - 1.3258252147247767*pow(eta,2.)
1793 - 1.3258252147247767*pow(xi,2.)
1794 + 3.97747564417433*pow(xi,2.)*pow(eta,2.);
1795 shape( 74,2) = 2.6516504294495533*xi*nu
1796 - 7.95495128834866*xi*pow(eta,2.)*nu;
1797 shape( 75,1) = 0.4419417382415922 - 1.3258252147247767*pow(eta,2.)
1798 - 1.3258252147247767*pow(xi,2.)
1799 + 3.97747564417433*pow(xi,2.)*pow(eta,2.);
1800 shape( 75,2) = 2.6516504294495533*eta*nu
1801 - 7.95495128834866*pow(xi,2.)*eta*nu;
1802 shape( 76,2) = 0.4419417382415922 - 1.3258252147247767*pow(eta,2.)
1803 - 1.3258252147247767*pow(xi,2.)
1804 + 3.97747564417433*pow(xi,2.)*pow(eta,2.);
1805 shape( 77,0) = -2.4302777619029476*xi*nu
1806 + 4.050462936504912*pow(xi,3.)*nu;
1807 shape( 77,2) = -0.4050462936504913 + 1.2151388809514738*pow(nu,2.)
1808 + 2.025231468252456*pow(xi,2.)
1809 - 6.075694404757369*pow(xi,2.)*pow(nu,2.);
1810 shape( 78,1) = -2.4302777619029476*xi*nu
1811 + 4.050462936504912*pow(xi,3.)*nu;
1812 shape( 79,0) = -2.4302777619029476*xi*eta
1813 + 4.050462936504912*pow(xi,3.)*eta;
1814 shape( 79,2) = 2.4302777619029476*eta*nu
1815 - 12.151388809514739*pow(xi,2.)*eta*nu;
1816 shape( 80,1) = 2.4302777619029476*xi*eta
1817 + 4.050462936504912*pow(xi,3.)*eta;
1818 shape( 80,2) = 2.4302777619029476*xi*nu
1819 - 4.050462936504912*pow(xi,3.)*nu;
1820 shape( 81,2) = -2.4302777619029476*xi*eta
1821 + 4.050462936504912*pow(xi,3.)*eta;
1822 shape( 82,0) = 0.397747564417433 - 3.97747564417433*pow(xi,2.)
1823 + 4.640388251536718*pow(xi,4.);
1824 shape( 82,2) = 7.95495128834866*xi*nu
1825 - 18.561553006146873*pow(xi,3.)*nu;
1826 shape( 83,1) = 0.397747564417433 - 3.97747564417433*pow(xi,2.)
1827 + 4.640388251536718*pow(xi,4.);
1828 shape( 84,2) = 0.397747564417433 - 3.97747564417433*pow(xi,2.)
1829 + 4.640388251536718*pow(xi,4.);
1830 case 3:
1831 shape( 26,0) = -1.403121520040228*nu + 2.3385358667337135*pow(nu,3.);
1832 shape( 27,1) = -1.403121520040228*nu + 2.3385358667337135*pow(nu,3.);
1833 shape( 28,0) = -0.6846531968814576*eta
1834 + 2.053959590644373*eta*pow(nu,2.);
1835 shape( 29,1) = -0.6846531968814576*eta
1836 + 2.053959590644373*eta*pow(nu,2.);
1837 shape( 29,2) = 0.6846531968814576*nu - 0.6846531968814576*pow(nu,3.);
1838 shape( 30,0) = -0.6846531968814576*nu
1839 + 2.053959590644373*pow(eta,2.)*nu;
1840 shape( 31,1) = -0.6846531968814576*nu
1841 + 2.053959590644373*pow(eta,2.)*nu;
1842 shape( 31,2) = 0.6846531968814576*eta
1843 - 2.053959590644373*eta*pow(nu,2.);
1844 shape( 32,0) = -1.403121520040228*eta
1845 + 2.3385358667337135*pow(eta,3.);
1846 shape( 33,1) = -1.403121520040228*eta
1847 + 2.3385358667337135*pow(eta,3.);
1848 shape( 33,2) = 1.403121520040228*nu
1849 - 7.0156076002011405*pow(eta,2.)*nu;
1850 shape( 34,2) = -1.403121520040228*eta
1851 + 2.3385358667337135*pow(eta,3.);
1852 shape( 35,0) = -0.6846531968814576*xi
1853 + 2.053959590644373*xi*pow(nu,2.);
1854 shape( 35,2) = 0.6846531968814576*nu - 0.6846531968814576*pow(nu,3.);
1855 shape( 36,1) = -0.6846531968814576*xi
1856 + 2.053959590644373*xi*pow(nu,2.);
1857 shape( 37,0) = 1.8371173070873836*xi*eta*nu;
1858 shape( 37,2) = 0.30618621784789724*eta
1859 - 0.9185586535436918*eta*pow(nu,2.);
1860 shape( 38,1) = 1.8371173070873836*xi*eta*nu;
1861 shape( 38,2) = 0.30618621784789724*xi
1862 - 0.9185586535436918*xi*pow(nu,2.);
1863 shape( 39,0) = -0.6846531968814576*xi
1864 + 2.053959590644373*xi*pow(eta,2.);
1865 shape( 39,2) = 0.6846531968814576*nu
1866 - 2.053959590644373*pow(eta,2.)*nu;
1867 shape( 40,1) = -0.6846531968814576*xi
1868 + 2.053959590644373*xi*pow(eta,2.);
1869 shape( 40,2) = -4.107919181288746*xi*eta*nu;
1870 shape( 41,2) = -0.6846531968814576*xi
1871 + 2.053959590644373*xi*pow(eta,2.);
1872 shape( 42,0) = -0.6846531968814576*nu
1873 + 2.053959590644373*pow(xi,2.)*nu;
1874 shape( 42,2) = 0.6846531968814576*xi
1875 - 2.053959590644373*xi*pow(nu,2.);
1876 shape( 43,1) = -0.6846531968814576*nu
1877 + 2.053959590644373*pow(xi,2.)*nu;
1878 shape( 44,0) = -0.6846531968814576*eta
1879 + 2.053959590644373*pow(xi,2.)*eta;
1880 shape( 44,2) = -4.107919181288746*xi*eta*nu;
1881 shape( 45,1) = -0.6846531968814576*eta
1882 + 2.053959590644373*pow(xi,2.)*eta;
1883 shape( 45,2) = 0.6846531968814576*nu
1884 - 2.053959590644373*pow(xi,2.)*nu;
1885 shape( 46,2) = -0.6846531968814576*eta
1886 + 2.053959590644373*pow(xi,2.)*eta;
1887 shape( 47,0) = -1.403121520040228*xi + 2.3385358667337135*pow(xi,3.);
1888 shape( 47,2) = 1.403121520040228*nu
1889 - 7.0156076002011405*pow(xi,2.)*nu;
1890 shape( 48,1) = -1.403121520040228*xi + 2.3385358667337135*pow(xi,3.);
1891 shape( 49,2) = -1.403121520040228*xi + 2.3385358667337135*pow(xi,3.);
1892 case 2:
1893 shape( 11,0) = -0.39528470752104744 + 1.1858541225631423*pow(nu,2.);
1894 shape( 12,1) = -0.39528470752104744 + 1.1858541225631423*pow(nu,2.);
1895 shape( 13,0) = 1.0606601717798212*eta*nu;
1896 shape( 14,1) = 1.0606601717798212*eta*nu;
1897 shape( 14,2) = 0.1767766952966369 - 0.5303300858899106*pow(nu,2.);
1898 shape( 15,0) = -0.39528470752104744 + 1.1858541225631423*pow(eta,2.);
1899 shape( 16,1) = -0.39528470752104744 + 1.1858541225631423*pow(eta,2.);
1900 shape( 16,2) = -2.3717082451262845*eta*nu;
1901 shape( 17,2) = -0.39528470752104744 + 1.1858541225631423*pow(eta,2.);
1902 shape( 18,0) = 1.0606601717798212*xi*nu;
1903 shape( 18,2) = 0.1767766952966369 - 0.5303300858899106*pow(nu,2.);
1904 shape( 19,1) = 1.0606601717798212*xi*nu;
1905 shape( 20,0) = 1.0606601717798212*xi*eta;
1906 shape( 20,2) = -1.0606601717798212*eta*nu;
1907 shape( 21,1) = 1.0606601717798212*xi*eta;
1908 shape( 21,2) = -1.0606601717798212*xi*nu;
1909 shape( 22,2) = 1.0606601717798212*xi*eta;
1910 shape( 23,0) = -0.39528470752104744 + 1.1858541225631423*pow(xi,2.);
1911 shape( 23,2) = -2.3717082451262845*xi*nu;
1912 shape( 24,1) = -0.39528470752104744 + 1.1858541225631423*pow(xi,2.);
1913 shape( 25,2) = -0.39528470752104744 + 1.1858541225631423*pow(xi,2.);
1914 case 1:
1915 shape( 3,0) = 0.6123724356957945*nu;
1916 shape( 4,1) = 0.6123724356957945*nu;
1917 shape( 5,0) = 0.6123724356957945*eta;
1918 shape( 6,1) = 0.6123724356957945*eta;
1919 shape( 6,2) = -0.6123724356957945*nu;
1920 shape( 7,2) = 0.6123724356957945*eta;
1921 shape( 8,0) = 0.6123724356957945*xi;
1922 shape( 8,2) = -0.6123724356957945*nu;
1923 shape( 9,1) = 0.6123724356957945*xi;
1924 shape( 10,2) = 0.6123724356957945*xi;
1925 case 0:
1926 shape( 0,0) = 0.3535533905932738;
1927 shape( 1,1) = 0.3535533905932738;
1928 shape( 2,2) = 0.3535533905932738;
1929 }
1930}
1931
1932}
1933
1934}
1935
1936#endif
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Abstract class for construction of IntegrationRules in cut elements.
virtual void SetLevelSetProjectionOrder(int order)
Coefficient * LvlSet
The zero level set of this Coefficient defines the cut surface.
virtual ~CutIntegrationRules()
Destructor of CutIntegrationRules.
virtual void SetLevelSetCoefficient(Coefficient &ls)
Change the Coefficient whose zero level set specifies the cut.
int lsOrder
Space order for the LS projection.
virtual void GetSurfaceIntegrationRule(ElementTransformation &Tr, IntegrationRule &result)=0
Construct a cut-surface IntegrationRule.
int Order
Order of the IntegrationRule.
virtual void GetVolumeIntegrationRule(ElementTransformation &Tr, IntegrationRule &result, const IntegrationRule *sir=NULL)=0
Construct a cut-volume IntegrationRule.
virtual void SetOrder(int order)
Change the order of the constructed IntegrationRule.
virtual void GetSurfaceWeights(ElementTransformation &Tr, const IntegrationRule &sir, Vector &weights)=0
Compute transformation quadrature weights for surface integration.
CutIntegrationRules(int order, Coefficient &lvlset, int lsO=2)
Constructor to set up the generated cut IntegrationRules.
Class for Singular Value Decomposition of a DenseMatrix.
Definition densemat.hpp:956
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
Rank 3 tensor (array of matrices)
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Class for subdomain IntegrationRules by means of moment-fitting.
void BasisAD2D(const IntegrationPoint &ip, DenseMatrix &shape)
Antiderivatives of the monomial basis on the element [-1,1]^2.
void mGSStep(DenseMatrix &shape, DenseTensor &shapeMFN, int step)
A step of the modified Gram-Schmidt algorithm.
void OrthoBasis3D(const IntegrationPoint &ip, DenseMatrix &shape)
A orthogonalized divergence free basis on the element [-1,1]^3.
Vector FaceWeightsComp
Indicates the already computed face IntegrationRules.
Array< IntegrationPoint > FaceIP
Array of face integration points.
int nBasis
Number of divergence-free basis functions for surface integration.
void ComputeSurfaceWeights2D(ElementTransformation &Tr)
Compute 2D quadrature weights.
void Init(int order, Coefficient &levelset, int lsO)
Initialize the MomentFittingIntRules.
~MomentFittingIntRules() override
Destructor of MomentFittingIntRules.
void Basis2D(const IntegrationPoint &ip, Vector &shape)
Monomial basis on the element [-1,1]^2.
void InitSurface(int order, Coefficient &levelset, int lsO, ElementTransformation &Tr)
Initialization for surface IntegrationRule.
void InitVolume(int order, Coefficient &levelset, int lsO, ElementTransformation &Tr)
Initialization for volume IntegrationRule.
DenseMatrix FaceWeights
Column-wise Matrix for the face quadrature weights.
void ComputeSurfaceWeights1D(ElementTransformation &Tr)
Compute 1D quadrature weights.
void SetOrder(int order) override
Change the order of the constructed IntegrationRule.
int dim
Space Dimension of the element.
void OrthoBasis2D(const IntegrationPoint &ip, DenseMatrix &shape)
A orthogonalized divergence free basis on the element [-1,1]^2.
void ComputeVolumeWeights1D(ElementTransformation &Tr, const IntegrationRule *sir)
Compute the 1D quadrature weights.
void ComputeVolumeWeights2D(ElementTransformation &Tr, const IntegrationRule *sir)
Compute the 2D quadrature weights.
DenseMatrixSVD * VolumeSVD
SVD of the matrix for volumetric IntegrationRules.
void Basis3D(const IntegrationPoint &ip, Vector &shape)
Monomial basis on the element [-1,1]^3.
void ComputeFaceWeights(ElementTransformation &Tr)
Compute the IntegrationRules on the faces.
void GetSurfaceIntegrationRule(ElementTransformation &Tr, IntegrationRule &result) override
Construct a cut-surface IntegrationRule.
void GetSurfaceWeights(ElementTransformation &Tr, const IntegrationRule &sir, Vector &weights) override
Compute transformation quadrature weights for surface integration.
void ComputeSurfaceWeights3D(ElementTransformation &Tr)
Compute 2D quadrature weights.
void DivFreeBasis2D(const IntegrationPoint &ip, DenseMatrix &shape)
A divergence free basis on the element [-1,1]^2.
IntegrationRule ir
IntegrationRule representing the reused IntegrationPoints.
int nBasisVolume
Number of basis functions for volume integration.
void Clear()
Clear stored data of the MomentFittingIntRules.
void BasisAD3D(const IntegrationPoint &ip, DenseMatrix &shape)
Antiderivatives of the monomial basis on the element [-1,1]^3.
MomentFittingIntRules(int order, Coefficient &lvlset, int lsO)
Constructor to set up the generated cut IntegrationRules.
void GetVolumeIntegrationRule(ElementTransformation &Tr, IntegrationRule &result, const IntegrationRule *sir=nullptr) override
Construct a cut-volume IntegrationRule.
void ComputeVolumeWeights3D(ElementTransformation &Tr, const IntegrationRule *sir)
Compute the 3D quadrature weights.
Vector data type.
Definition vector.hpp:80
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
real_t lvlset(const Vector &X)
Level-set function defining the implicit interface.
Definition ex38.cpp:49
void GetDivFree3DBasis(const Vector &X, DenseMatrix &shape, int Order)
3 dimensional divergence free basis functions on [-1,1]^3
float real_t
Definition config.hpp:43