OOFEM 3.0
Loading...
Searching...
No Matches
tet1_3d_supg.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35#include "tet1_3d_supg.h"
36#include "fei3dtetlin.h"
37#include "fluidmodel.h"
38#include "node.h"
39#include "material.h"
40#include "gausspoint.h"
43#include "fluidcrosssection.h"
44#include "floatmatrix.h"
45#include "floatarray.h"
46#include "intarray.h"
47#include "domain.h"
48#include "mathfem.h"
49#include "engngm.h"
50#include "timestep.h"
51#include "crosssection.h"
52#include "classfactory.h"
53
54#ifdef __OOFEG
55 #include "oofeggraphiccontext.h"
56#endif
57
58namespace oofem {
60
61FEI3dTetLin Tet1_3D_SUPG :: interpolation;
62
63Tet1_3D_SUPG :: Tet1_3D_SUPG(int n, Domain *aDomain) :
64 SUPGElement2(n, aDomain)
65{
67}
68
69int
70Tet1_3D_SUPG :: computeNumberOfDofs()
71{
72 return 16;
73}
74
75void
76Tet1_3D_SUPG :: giveDofManDofIDMask(int inode, IntArray &answer) const
77{
78 answer = {V_u, V_v, V_w, P_f};
79}
80
81
82void
83Tet1_3D_SUPG :: computeGaussPoints()
84// Sets up the array containing the integration points of the receiver.
85{
86 if ( integrationRulesArray.size() == 0 ) {
87 integrationRulesArray.resize(3);
88 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 3);
89 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], 1, this);
90
91 integrationRulesArray [ 1 ] = std::make_unique<GaussIntegrationRule>(2, this, 1, 3);
92 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 1 ], 4, this);
93
94 integrationRulesArray [ 2 ] = std::make_unique<GaussIntegrationRule>(3, this, 1, 3);
95 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 2 ], 4, this);
96 }
97}
98
99
100Interface *
101Tet1_3D_SUPG :: giveInterface(InterfaceType interface)
102{
103 if ( interface == LevelSetPCSElementInterfaceType ) {
104 return static_cast< LevelSetPCSElementInterface * >(this);
105 }
106
107 return NULL;
108}
109
110
111void
112Tet1_3D_SUPG :: computeNuMatrix(FloatMatrix &answer, GaussPoint *gp)
113{
114 FloatArray n;
116 answer.beNMatrixOf(n, 3);
117}
118
119void
120Tet1_3D_SUPG :: computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
121{
122 FloatMatrix n, dn;
123 FloatArray u, un;
125 this->computeNuMatrix(n, gp);
126 this->computeVectorOfVelocities(VM_Total, tStep, un);
127
128 u.beProductOf(n, un);
129
130 answer.resize(3, 12);
131 answer.zero();
132 for ( int i = 1; i <= 4; i++ ) {
133 answer.at(1, 3 * i - 2) = u.at(1) * dn.at(i, 1) + u.at(2) * dn.at(i, 2) + u.at(3) * dn.at(i, 3);
134 answer.at(2, 3 * i - 1) = u.at(1) * dn.at(i, 1) + u.at(2) * dn.at(i, 2) + u.at(3) * dn.at(i, 3);
135 answer.at(3, 3 * i - 0) = u.at(1) * dn.at(i, 1) + u.at(2) * dn.at(i, 2) + u.at(3) * dn.at(i, 3);
136 }
137}
138
139void
140Tet1_3D_SUPG :: computeDivTauMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
141{
142 answer.resize(3, 12);
143 answer.zero();
144}
145
146
147void
148Tet1_3D_SUPG :: computeGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
149{
150 FloatArray u;
151 FloatMatrix dn, um(3, 4);
152
153 this->computeVectorOfVelocities(VM_Total, tStep, u);
154
156 for ( int i = 1; i <= 4; i++ ) {
157 um.at(1, i) = u.at(3 * i - 2);
158 um.at(2, i) = u.at(3 * i - 1);
159 um.at(3, i) = u.at(3 * i - 0);
160 }
161
162 answer.beProductOf(um, dn);
163}
164
165
166void
167Tet1_3D_SUPG :: computeBMatrix(FloatMatrix &answer, GaussPoint *gp)
168{
169 FloatMatrix dn;
171
172 answer.resize(6, 12);
173 answer.zero();
174
175 for ( int i = 1; i <= 4; i++ ) {
176 answer.at(1, 3 * i - 2) = dn.at(i, 1);
177 answer.at(2, 3 * i - 1) = dn.at(i, 2);
178 answer.at(3, 3 * i - 0) = dn.at(i, 3);
179
180 answer.at(4, 3 * i - 1) = dn.at(i, 3);
181 answer.at(4, 3 * i - 0) = dn.at(i, 2);
182
183 answer.at(5, 3 * i - 2) = dn.at(i, 3);
184 answer.at(5, 3 * i - 0) = dn.at(i, 1);
185
186 answer.at(6, 3 * i - 2) = dn.at(i, 2);
187 answer.at(6, 3 * i - 1) = dn.at(i, 1);
188 }
189}
190
191void
192Tet1_3D_SUPG :: computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp)
193{
194 FloatMatrix dn;
196
197 answer.resize(1, 12);
198 answer.zero();
199
200 for ( int i = 1; i <= 4; i++ ) {
201 answer.at(1, 3 * i - 2) = dn.at(i, 1);
202 answer.at(1, 3 * i - 1) = dn.at(i, 2);
203 answer.at(1, 3 * i - 0) = dn.at(i, 3);
204 }
205}
206
207void
208Tet1_3D_SUPG :: computeNpMatrix(FloatMatrix &answer, GaussPoint *gp)
209{
210 FloatArray n;
212
213 answer.resize(1, 4);
214 answer.zero();
215
216 answer.at(1, 1) = n.at(1);
217 answer.at(1, 2) = n.at(2);
218 answer.at(1, 3) = n.at(3);
219 answer.at(1, 4) = n.at(4);
220}
221
222
223void
224Tet1_3D_SUPG :: computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp)
225{
226 FloatMatrix dn;
228
229 answer.beTranspositionOf(dn);
230}
231
232
233
234void
235Tet1_3D_SUPG :: updateStabilizationCoeffs(TimeStep *tStep)
236{
237 //TR1_2D_SUPG :: updateStabilizationCoeffs (tStep);
238 /* UGN-Based Stabilization */
239 double h_ugn, sum = 0.0, vnorm, t_sugn1, t_sugn2, t_sugn3, u_1, u_2, u_3, z, Re_ugn;
240 double dscale, uscale, lscale, tscale, dt;
241 //bool zeroFlag = false;
242 int im1;
243 FloatArray u, divu;
244 FloatMatrix du;
245
246 uscale = domain->giveEngngModel()->giveVariableScale(VST_Velocity);
247 lscale = domain->giveEngngModel()->giveVariableScale(VST_Length);
248 tscale = domain->giveEngngModel()->giveVariableScale(VST_Time);
249 dscale = domain->giveEngngModel()->giveVariableScale(VST_Density);
250
251 this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), u);
252
253 u.times(uscale);
254 double nu;
255
256 // compute averaged viscosity based on rule of mixture
257
258 dt = tStep->giveTimeIncrement() * tscale;
259
260 std :: unique_ptr< IntegrationRule > &iRule = this->integrationRulesArray [ 1 ];
261 GaussPoint *gp_first = iRule->getIntegrationPoint(0);
262 nu = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->giveEffectiveViscosity(gp_first, tStep->givePreviousStep() );
263 nu *= domain->giveEngngModel()->giveVariableScale(VST_Viscosity);
264
265 for ( auto *gp: *iRule ) {
266 this->computeDivUMatrix(du, gp);
267 divu.beProductOf(du, u);
268 sum += divu.at(1);
269 }
270
271 sum *= ( 1. / lscale / iRule->giveNumberOfIntegrationPoints() );
272
273 /*
274 * for (i=1; i<=3;i++) {
275 * im1=i-1;
276 * sum+= fabs(u.at((im1)*2+1)*b[im1]/lscale + u.at(im1*2+2)*c[im1]/lscale);
277 * }
278 */
279 vnorm = 0.;
280 int nsd = this->giveNumberOfSpatialDimensions();
281 for ( int i = 1; i <= numberOfDofMans; i++ ) {
282 im1 = i - 1;
283 u_1 = u.at( ( im1 ) * nsd + 1 );
284 u_2 = u.at( ( im1 ) * nsd + 2 );
285 if ( nsd > 2 ) {
286 u_3 = u.at( ( im1 ) * nsd + 3 );
287 } else {
288 u_3 = 0.;
289 }
290
291 vnorm = max( vnorm, sqrt(u_1 * u_1 + u_2 * u_2 + u_3 * u_3) );
292 }
293
294 if ( ( vnorm == 0.0 ) || ( sum < vnorm * 1e-10 ) ) {
295 //t_sugn1 = inf;
296 t_sugn2 = dt / 2.0;
297 //t_sugn3 = inf;
298 this->t_supg = 1. / sqrt( 1. / ( t_sugn2 * t_sugn2 ) );
299 this->t_pspg = this->t_supg;
300 this->t_lsic = 0.0;
301 } else {
302 h_ugn = 2.0 * vnorm / sum;
303 t_sugn1 = 1. / sum;
304 t_sugn2 = dt / 2.0;
305 t_sugn3 = h_ugn * h_ugn / 4.0 / nu;
306
307 this->t_supg = 1. / sqrt( 1. / ( t_sugn1 * t_sugn1 ) + 1. / ( t_sugn2 * t_sugn2 ) + 1. / ( t_sugn3 * t_sugn3 ) );
308 this->t_pspg = this->t_supg;
309
310 Re_ugn = vnorm * h_ugn / ( 2. * nu );
311 z = ( Re_ugn <= 3. ) ? Re_ugn / 3. : 1.0;
312 this->t_lsic = h_ugn * vnorm * z / 2.0;
313 }
314
315 // if (this->number == 1) {
316 // printf ("t_supg %e t_pspg %e t_lsic %e\n", t_supg, t_pspg, t_lsic);
317 // }
318
319
320 this->t_supg *= uscale / lscale;
321 this->t_pspg *= 1. / ( lscale * dscale );
322 this->t_lsic *= ( dscale * uscale ) / ( lscale * lscale );
323
324 this->t_lsic = 0.0;
325
326 //this->t_lsic=0.0;
327 //this->t_pspg=0.0;
328}
329
330int
331Tet1_3D_SUPG :: giveNumberOfSpatialDimensions()
332{
333 return 3;
334}
335
336double
337Tet1_3D_SUPG :: computeCriticalTimeStep(TimeStep *tStep)
338{
339 FloatArray u;
340 double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
341
342 this->computeVectorOfVelocities(VM_Total, tStep, u);
343
344 double vn1 = sqrt( u.at(1) * u.at(1) + u.at(2) * u.at(2) + u.at(3) * u.at(3) );
345 double vn2 = sqrt( u.at(4) * u.at(4) + u.at(5) * u.at(5) + u.at(6) * u.at(6) );
346 double vn3 = sqrt( u.at(7) * u.at(7) + u.at(8) * u.at(8) + u.at(9) * u.at(9) );
347 double vn4 = sqrt( u.at(10) * u.at(10) + u.at(11) * u.at(11) + u.at(12) * u.at(12) );
348 double veln = max( vn1, max( vn2, max(vn3, vn4) ) );
349
350 double ln = 1.e6;
351 Node *inode, *jnode, *knode, *lnode;
352 FloatArray t1, t2, n3, n;
353 for ( int l = 1; l <= 4; l++ ) {
354 int i = ( l > 3 ) ? 1 : l + 1;
355 int j = ( i > 3 ) ? 1 : i + 1;
356 int k = ( j > 3 ) ? 1 : j + 1;
357
358 inode = this->giveNode(i);
359 jnode = this->giveNode(j);
360 knode = this->giveNode(k);
361 lnode = this->giveNode(l);
362 t1.beDifferenceOf(inode->giveCoordinates(), jnode->giveCoordinates());
363
364 t2.beDifferenceOf(knode->giveCoordinates(), jnode->giveCoordinates());
365
366 n.beVectorProductOf(t1, t2);
367 n.normalize();
368
369 n3.beDifferenceOf(lnode->giveCoordinates(), jnode->giveCoordinates());
370
371 ln = min( ln, sqrt( fabs( n.dotProduct(n3) ) ) );
372 }
373
374 if ( veln != 0.0 ) {
375 return ln / veln;
376 } else {
377 return 0.5 * ln * ln * Re;
378 }
379}
380
381
382double
383Tet1_3D_SUPG :: computeVolumeAround(GaussPoint *gp)
384// Returns the portion of the receiver which is attached to gp.
385{
386 double determinant, weight, volume;
387 determinant = fabs( this->interpolation.giveTransformationJacobian( gp->giveNaturalCoordinates(),
389
390 weight = gp->giveWeight();
391 volume = determinant * weight;
392
393 return volume;
394}
395
396
397void
398Tet1_3D_SUPG :: computeDeviatoricStress(FloatArray &answer, const FloatArray &eps, GaussPoint *gp, TimeStep *tStep)
399{
400 answer = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->computeDeviatoricStress3D(eps, gp, tStep);
401}
402
403
404void
405Tet1_3D_SUPG :: computeTangent(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
406{
407 answer = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->computeTangent3D(mode, gp, tStep);
408}
409
410
411double
412Tet1_3D_SUPG :: LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep)
413{
414 double answer = 0.0, norm, dV, vol = 0.0;
415 FloatMatrix n, dn;
416 FloatArray fi(4), u, un, gfi;
417
418 this->computeVectorOfVelocities(VM_Total, tStep, un);
419
420 for ( int i = 1; i <= 4; i++ ) {
421 fi.at(i) = ls->giveLevelSetDofManValue( dofManArray.at(i) );
422 }
423
424 for ( GaussPoint *gp: *this->integrationRulesArray [ 0 ] ) {
425 dV = this->computeVolumeAround(gp);
426 interpolation.evaldNdx( dn, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
427 this->computeNuMatrix(n, gp);
428 u.beProductOf(n, un);
429 gfi.beTProductOf(dn, fi);
430 norm = gfi.computeNorm();
431
432 vol += dV;
433 answer += dV * u.dotProduct(gfi) / norm;
434 }
435
436 return answer / vol;
437}
438
439void
440Tet1_3D_SUPG :: LS_PCS_computedN(FloatMatrix &answer)
441{
442 GaussPoint *gp = this->integrationRulesArray [ 0 ]->getIntegrationPoint(0);
443 interpolation.evaldNdx( answer, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
444}
445
446
447double
448Tet1_3D_SUPG :: LS_PCS_computeVolume()
449{
450 double answer = 0.0;
451
452 for ( GaussPoint *gp: *this->integrationRulesArray [ 0 ] ) {
453 answer += this->computeVolumeAround(gp);
454 }
455
456 return answer;
457}
458
459double
460Tet1_3D_SUPG :: LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep)
461{
462 FloatArray voff(2), fi(4);
463 for ( int i = 1; i <= 4; i++ ) {
464 fi.at(i) = ls->giveLevelSetDofManValue( dofManArray.at(i) );
465 }
466
467 this->LS_PCS_computeVOFFractions(voff, fi);
468 return ( voff.at(1) - voff.at(2) );
469}
470
471
472
473
474void
475Tet1_3D_SUPG :: LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi)
476{
477 int neg = 0, pos = 0, zero = 0, si = 0;
478 double x1, y1, z1;
479 answer.resize(2);
480
481 for ( double ifi: fi ) {
482 if ( ifi >= 0. ) {
483 pos++;
484 } else if ( ifi < 0. ) {
485 neg++;
486 } else {
487 zero++;
488 }
489 }
490
491 if ( neg == 0 ) { // all level set values positive
492 answer.at(1) = 1.0;
493 answer.at(2) = 0.0;
494 return; //return area;
495 } else if ( pos == 0 ) { // all level set values negative
496 answer.at(1) = 0.0;
497 answer.at(2) = 1.0;
498 return; //return -area;
499 } else if ( zero == 4 ) {
500 // ???????
501 answer.at(1) = 1.0;
502 answer.at(2) = 0.0;
503 return;
504 } else {
505 // zero level set inside
506 // distinguish two poosible cases
507 if ( max(pos, neg) == 3 ) {
508 // find the vertex vith level set sign different from other three
509 for ( int i = 1; i <= 4; i++ ) {
510 if ( ( pos > neg ) && ( fi.at(i) < 0.0 ) ) {
511 si = i;
512 break;
513 }
514
515 if ( ( pos < neg ) && ( fi.at(i) >= 0.0 ) ) {
516 si = i;
517 break;
518 }
519 }
520
521 if ( si ) {
522 x1 = this->giveNode(si)->giveCoordinate(1);
523 y1 = this->giveNode(si)->giveCoordinate(2);
524 z1 = this->giveNode(si)->giveCoordinate(3);
525
526
527 int ii;
528 double t, xi [ 3 ], yi [ 3 ], zi [ 3 ];
529 // compute intersections with element sides originating from this vertex
530 for ( int i = 0; i < 3; i++ ) {
531 ii = ( si + i ) % 4 + 1;
532 t = fi.at(si) / ( fi.at(si) - fi.at(ii) );
533 xi [ i ] = x1 + t * ( this->giveNode(ii)->giveCoordinate(1) - x1 );
534 yi [ i ] = y1 + t * ( this->giveNode(ii)->giveCoordinate(2) - y1 );
535 zi [ i ] = z1 + t * ( this->giveNode(ii)->giveCoordinate(3) - z1 );
536 }
537
538 // compute volume of this pyramid (si, xi[0], xi[1], xi[2])
539 double __vol = fabs( ( 1. / 6. ) * ( ( x1 - xi [ 0 ] ) * ( ( yi [ 1 ] - yi [ 0 ] ) * ( zi [ 2 ] - zi [ 0 ] ) - ( zi [ 1 ] - zi [ 0 ] ) * ( yi [ 2 ] - yi [ 0 ] ) ) +
540 ( y1 - yi [ 0 ] ) * ( ( zi [ 1 ] - zi [ 0 ] ) * ( xi [ 2 ] - xi [ 0 ] ) - ( xi [ 1 ] - xi [ 0 ] ) * ( zi [ 2 ] - zi [ 0 ] ) ) +
541 ( z1 - zi [ 0 ] ) * ( ( xi [ 1 ] - xi [ 0 ] ) * ( yi [ 2 ] - yi [ 0 ] ) - ( yi [ 1 ] - yi [ 0 ] ) * ( xi [ 2 ] - xi [ 0 ] ) ) ) );
542
543
544 double vol = LS_PCS_computeVolume();
545 if ( ( fabs(__vol) - vol ) < 0.0000001 ) {
546 __vol = sgn(__vol) * vol;
547 }
548
549 if ( ( __vol < 0 ) || ( fabs(__vol) / vol > 1.0000001 ) ) {
550 OOFEM_ERROR("internal consistency error");
551 }
552
553 if ( pos > neg ) {
554 // negative vol computed
555 answer.at(2) = min(fabs(__vol) / vol, 1.0);
556 answer.at(1) = 1.0 - answer.at(2);
557 } else {
558 // postive vol computed
559 answer.at(1) = min(fabs(__vol) / vol, 1.0);
560 answer.at(2) = 1.0 - answer.at(1);
561 }
562 } else {
563 OOFEM_ERROR("internal consistency error");
564 }
565 } else if ( max(pos, neg) == 2 ) {
566 // two vertices positive; two negative; compute positive volume
567 int p1 = 0, p2 = 0;
568 // find the vertex vith level set sign different from other three
569 for ( int i = 1; i <= 4; i++ ) {
570 if ( fi.at(i) >= 0.0 ) {
571 if ( p1 ) {
572 p2 = i;
573 break;
574 } else {
575 p1 = i;
576 }
577 }
578 }
579
580 int _ind, ii;
581 double t;
582 double p1i_x [ 3 ], p1i_y [ 3 ], p1i_z [ 3 ];
583 double p2i_x [ 3 ], p2i_y [ 3 ], p2i_z [ 3 ];
584
585 if ( p1 && p2 ) {
586 // find the two intersections sharing edge with p1 and p2
587 p1i_x [ 0 ] = this->giveNode(p1)->giveCoordinate(1);
588 p1i_y [ 0 ] = this->giveNode(p1)->giveCoordinate(2);
589 p1i_z [ 0 ] = this->giveNode(p1)->giveCoordinate(3);
590
591 p2i_x [ 0 ] = this->giveNode(p2)->giveCoordinate(1);
592 p2i_y [ 0 ] = this->giveNode(p2)->giveCoordinate(2);
593 p2i_z [ 0 ] = this->giveNode(p2)->giveCoordinate(3);
594
595 _ind = 1;
596 for ( int i = 0; i < 3; i++ ) {
597 ii = ( p1 + i ) % 4 + 1;
598 if ( ( ii == p2 ) || ( ii == p1 ) ) {
599 continue;
600 }
601
602 t = fi.at(p1) / ( fi.at(p1) - fi.at(ii) );
603 p1i_x [ _ind ] = p1i_x [ 0 ] + t * ( this->giveNode(ii)->giveCoordinate(1) - p1i_x [ 0 ] );
604 p1i_y [ _ind ] = p1i_y [ 0 ] + t * ( this->giveNode(ii)->giveCoordinate(2) - p1i_y [ 0 ] );
605 p1i_z [ _ind ] = p1i_z [ 0 ] + t * ( this->giveNode(ii)->giveCoordinate(3) - p1i_z [ 0 ] );
606
607 t = fi.at(p2) / ( fi.at(p2) - fi.at(ii) );
608 p2i_x [ _ind ] = p2i_x [ 0 ] + t * ( this->giveNode(ii)->giveCoordinate(1) - p2i_x [ 0 ] );
609 p2i_y [ _ind ] = p2i_y [ 0 ] + t * ( this->giveNode(ii)->giveCoordinate(2) - p2i_y [ 0 ] );
610 p2i_z [ _ind++ ] = p2i_z [ 0 ] + t * ( this->giveNode(ii)->giveCoordinate(3) - p2i_z [ 0 ] );
611 }
612
613 // compute volume of this wedge as a sum of volumes of three
614 // pyramids
615 double __v1 = ( ( p2i_x [ 0 ] - p1i_x [ 0 ] ) * ( p1i_y [ 1 ] - p1i_y [ 0 ] ) * ( p1i_z [ 2 ] - p1i_z [ 0 ] ) -
616 ( p2i_x [ 0 ] - p1i_x [ 0 ] ) * ( p1i_y [ 2 ] - p1i_y [ 0 ] ) * ( p1i_z [ 1 ] - p1i_z [ 0 ] ) +
617 ( p1i_x [ 2 ] - p1i_x [ 0 ] ) * ( p2i_y [ 0 ] - p1i_y [ 0 ] ) * ( p1i_z [ 1 ] - p1i_z [ 0 ] ) -
618 ( p1i_x [ 1 ] - p1i_x [ 0 ] ) * ( p2i_y [ 0 ] - p1i_y [ 0 ] ) * ( p1i_z [ 2 ] - p1i_z [ 0 ] ) +
619 ( p1i_x [ 1 ] - p1i_x [ 0 ] ) * ( p1i_y [ 2 ] - p1i_y [ 0 ] ) * ( p2i_z [ 0 ] - p1i_z [ 0 ] ) -
620 ( p1i_x [ 2 ] - p1i_x [ 0 ] ) * ( p1i_y [ 1 ] - p1i_y [ 0 ] ) * ( p2i_z [ 0 ] - p1i_z [ 0 ] ) ) / 6.0;
621
622 double __v2 = ( ( p2i_x [ 0 ] - p1i_x [ 1 ] ) * ( p1i_y [ 2 ] - p1i_y [ 1 ] ) * ( p2i_z [ 1 ] - p1i_z [ 1 ] ) -
623 ( p2i_x [ 0 ] - p1i_x [ 1 ] ) * ( p2i_y [ 1 ] - p1i_y [ 1 ] ) * ( p1i_z [ 2 ] - p1i_z [ 1 ] ) +
624 ( p2i_x [ 1 ] - p1i_x [ 1 ] ) * ( p2i_y [ 0 ] - p1i_y [ 1 ] ) * ( p1i_z [ 2 ] - p1i_z [ 1 ] ) -
625 ( p1i_x [ 2 ] - p1i_x [ 1 ] ) * ( p2i_y [ 0 ] - p1i_y [ 1 ] ) * ( p2i_z [ 1 ] - p1i_z [ 1 ] ) +
626 ( p1i_x [ 2 ] - p1i_x [ 1 ] ) * ( p2i_y [ 1 ] - p1i_y [ 1 ] ) * ( p2i_z [ 0 ] - p1i_z [ 1 ] ) -
627 ( p2i_x [ 1 ] - p1i_x [ 1 ] ) * ( p1i_y [ 2 ] - p1i_y [ 1 ] ) * ( p2i_z [ 0 ] - p1i_z [ 1 ] ) ) / 6.0;
628
629 double __v3 = ( ( p1i_x [ 2 ] - p2i_x [ 0 ] ) * ( p2i_y [ 1 ] - p2i_y [ 0 ] ) * ( p2i_z [ 2 ] - p2i_z [ 0 ] ) -
630 ( p1i_x [ 2 ] - p2i_x [ 0 ] ) * ( p2i_y [ 2 ] - p2i_y [ 0 ] ) * ( p2i_z [ 1 ] - p2i_z [ 0 ] ) +
631 ( p2i_x [ 2 ] - p2i_x [ 0 ] ) * ( p1i_y [ 2 ] - p2i_y [ 0 ] ) * ( p2i_z [ 1 ] - p2i_z [ 0 ] ) -
632 ( p2i_x [ 1 ] - p2i_x [ 0 ] ) * ( p1i_y [ 2 ] - p2i_y [ 0 ] ) * ( p2i_z [ 2 ] - p2i_z [ 0 ] ) +
633 ( p2i_x [ 1 ] - p2i_x [ 0 ] ) * ( p2i_y [ 2 ] - p2i_y [ 0 ] ) * ( p1i_z [ 2 ] - p2i_z [ 0 ] ) -
634 ( p2i_x [ 2 ] - p2i_x [ 0 ] ) * ( p2i_y [ 1 ] - p2i_y [ 0 ] ) * ( p1i_z [ 2 ] - p2i_z [ 0 ] ) ) / 6.0;
635
636 double __vol = fabs(__v1) + fabs(__v2) + fabs(__v3);
637 double vol = LS_PCS_computeVolume();
638
639 if ( ( __vol < 0 ) || ( fabs(__vol) / vol > 1.0000001 ) ) {
640 OOFEM_ERROR("internal consistency error");
641 }
642
643 answer.at(1) = min(fabs(__vol) / vol, 1.0);
644 answer.at(2) = 1.0 - answer.at(1);
645 } else {
646 OOFEM_ERROR("internal consistency error");
647 }
648 } else {
649 OOFEM_ERROR("internal consistency error");
650 }
651 }
652}
653
654
655#ifdef __OOFEG
656
657void Tet1_3D_SUPG :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
658{
659 WCRec p [ 4 ];
660 GraphicObj *go;
661
662 if ( !gc.testElementGraphicActivity(this) ) {
663 return;
664 }
665
666 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
667 EASValsSetColor( gc.getElementColor() );
668 EASValsSetEdgeColor( gc.getElementEdgeColor() );
669 EASValsSetEdgeFlag(true);
670 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
671 EASValsSetFillStyle(FILL_SOLID);
672 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
673 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
674 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
675 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
676 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
677 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
678 p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
679 p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
680 p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveCoordinate(3);
681 p [ 3 ].x = ( FPNum ) this->giveNode(4)->giveCoordinate(1);
682 p [ 3 ].y = ( FPNum ) this->giveNode(4)->giveCoordinate(2);
683 p [ 3 ].z = ( FPNum ) this->giveNode(4)->giveCoordinate(3);
684
685 go = CreateTetra(p);
686 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
687 EGAttachObject(go, ( EObjectP ) this);
688 EMAddGraphicsToModel(ESIModel(), go);
689}
690
691#endif
692} // end namespace oofem
#define REGISTER_Element(class)
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
Node * giveNode(int i) const
Definition element.h:629
IntArray dofManArray
Array containing dofmanager numbers.
Definition element.h:138
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition fmelement.C:44
double computeNorm() const
Definition floatarray.C:861
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
double dotProduct(const FloatArray &x) const
Definition floatarray.C:524
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Definition floatarray.C:403
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
void times(double s)
Definition floatarray.C:834
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void beNMatrixOf(const FloatArray &n, int nsd)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
virtual double giveReynoldsNumber()=0
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
double giveLevelSetDofManValue(int i)
Returns level set value in specific node.
SUPGElement2(int n, Domain *aDomain)
static FEI3dTetLin interpolation
void computeNuMatrix(FloatMatrix &answer, GaussPoint *gp) override
void computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp) override
void LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi) override
double computeVolumeAround(GaussPoint *gp) override
int giveNumberOfSpatialDimensions() override
double LS_PCS_computeVolume() override
Returns receiver's volume.
double giveTimeIncrement()
Returns solution step associated time increment.
Definition timestep.h:168
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
Definition timestep.C:132
#define OOFEM_ERROR(...)
Definition error.h:79
double norm(const FloatArray &x)
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ VST_Length
@ VST_Viscosity
@ VST_Velocity
@ VST_Density
FloatArrayF< N > max(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
double sum(const FloatArray &x)
FloatMatrixF< N, M > zero()
Constructs a zero matrix (this is the default behavior when constructing a matrix,...
double sgn(double i)
Returns the signum of given value (if value is < 0 returns -1, otherwise returns 1).
Definition mathfem.h:104
@ LevelSetPCSElementInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_RAW_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_LAYER

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011