OOFEM 3.0
Loading...
Searching...
No Matches
quad10_2d_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 "quad10_2d_supg.h"
36#include "fei2dquadlin.h"
37#include "fei2dquadconst.h"
38#include "dofmanager.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 "materialinterface.h"
52#include "contextioerr.h"
53#include "crosssection.h"
54#include "classfactory.h"
55
56#ifdef __OOFEG
57 #include "oofeggraphiccontext.h"
58#endif
59
60namespace oofem {
62
63FEI2dQuadLin Quad10_2D_SUPG :: velocityInterpolation(1, 2);
64FEI2dQuadConst Quad10_2D_SUPG :: pressureInterpolation(1, 2);
65
66
67Quad10_2D_SUPG :: Quad10_2D_SUPG(int n, Domain *aDomain) :
68 SUPGElement2(n, aDomain), ZZNodalRecoveryModelInterface(this), pressureNode(1, aDomain, this)
69{
71}
72
74Quad10_2D_SUPG :: giveInterpolation() const
75{
76 return & this->velocityInterpolation;
77}
78
80Quad10_2D_SUPG :: giveInterpolation(DofIDItem id) const
81{
82 if ( id == P_f ) {
83 return & this->pressureInterpolation;
84 } else {
85 return & this->velocityInterpolation;
86 }
87}
88
90Quad10_2D_SUPG :: giveInternalDofManager(int i) const
91{
92 return const_cast< ElementDofManager * >(& pressureNode);
93}
94
95
96int
97Quad10_2D_SUPG :: computeNumberOfDofs()
98{
99 return 9;
100}
101
102
103void
104Quad10_2D_SUPG :: giveDofManDofIDMask(int inode, IntArray &answer) const
105{
106 answer = {V_u, V_v};
107}
108
109
110void
111Quad10_2D_SUPG :: giveInternalDofManDofIDMask(int i, IntArray &answer) const
112{
113 answer = {P_f};
114}
115
116
117void
118Quad10_2D_SUPG :: initializeFrom(InputRecord &ir, int priority)
119{
120 this->pressureNode.initializeFrom(ir, priority);
121
122 SUPGElement2 :: initializeFrom(ir, priority);
123}
124
125
126void
127Quad10_2D_SUPG :: giveInputRecord(DynamicInputRecord &input)
128{
129 SUPGElement2 :: giveInputRecord(input);
130 this->pressureNode.giveInputRecord(input);
131}
132
133
134void
135Quad10_2D_SUPG :: computeGaussPoints()
136// Sets up the array containing the four Gauss points of the receiver.
137{
138 if ( integrationRulesArray.size() == 0 ) {
139 integrationRulesArray.resize(3);
140
141
142 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 3);
143 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], 4, this);
144
145 //seven point Gauss integration
146 integrationRulesArray [ 1 ] = std::make_unique<GaussIntegrationRule>(2, this, 1, 3);
147 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 1 ], 4, this);
148
149 integrationRulesArray [ 2 ] = std::make_unique<GaussIntegrationRule>(3, this, 1, 3);
150 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 3 ], 4, this);
151 }
152}
153
154
155void
156Quad10_2D_SUPG :: computeDeviatoricStress(FloatArray &answer, const FloatArray &eps, GaussPoint *gp, TimeStep *tStep)
157{
158 answer = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->computeDeviatoricStress2D(eps, gp, tStep);
159}
160
161
162void
163Quad10_2D_SUPG :: computeTangent(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
164{
165 answer = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->computeTangent2D(mode, gp, tStep);
166}
167
168
169void
170Quad10_2D_SUPG :: computeNuMatrix(FloatMatrix &answer, GaussPoint *gp)
171{
172 FloatArray n;
174 answer.beNMatrixOf(n, 2);
175}
176
177
178void
179Quad10_2D_SUPG :: computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
180{
181 FloatMatrix n, dn;
182 FloatArray u, un;
184 this->computeNuMatrix(n, gp);
185 this->computeVectorOfVelocities(VM_Total, tStep, un);
186
187 u.beProductOf(n, un);
188
189 answer.resize(2, 8);
190 answer.zero();
191 for ( int i = 1; i <= 4; i++ ) {
192 answer.at(1, 2 * i - 1) = dn.at(i, 1) * u.at(1) + dn.at(i, 2) * u.at(2);
193 answer.at(2, 2 * i) = dn.at(i, 1) * u.at(1) + dn.at(i, 2) * u.at(2);
194 }
195}
196
197void
198Quad10_2D_SUPG :: computeBMatrix(FloatMatrix &answer, GaussPoint *gp)
199{
200 FloatMatrix dn;
202
203 answer.resize(3, 8);
204 answer.zero();
205
206 for ( int i = 1; i <= 4; i++ ) {
207 answer.at(1, 2 * i - 1) = dn.at(i, 1);
208 answer.at(2, 2 * i) = dn.at(i, 2);
209 answer.at(3, 2 * i - 1) = dn.at(i, 2);
210 answer.at(3, 2 * i) = dn.at(i, 1);
211 }
212}
213
214void
215Quad10_2D_SUPG :: computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp)
216{
217 FloatMatrix dn;
219
220 answer.resize(1, 8);
221 answer.zero();
222
223 for ( int i = 1; i <= 4; i++ ) {
224 answer.at(1, 2 * i - 1) = dn.at(i, 1);
225 answer.at(1, 2 * i) = dn.at(i, 2);
226 }
227}
228
229void
230Quad10_2D_SUPG :: computeNpMatrix(FloatMatrix &answer, GaussPoint *gp)
231{
232 FloatArray n;
234
235 answer.resize(1, 1);
236 answer.zero();
237
238 answer.at(1, 1) = n.at(1);
239}
240
241
242void
243Quad10_2D_SUPG :: computeGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
244{
245 FloatArray dnx(4), dny(4), u, u1(4), u2(4);
246 FloatMatrix dn;
247
248 answer.resize(2, 2);
249 answer.zero();
250
251 this->computeVectorOfVelocities(VM_Total, tStep, u);
252
254 for ( int i = 1; i <= 4; i++ ) {
255 dnx.at(i) = dn.at(i, 1);
256 dny.at(i) = dn.at(i, 2);
257
258 u1.at(i) = u.at(2 * i - 1);
259 u2.at(i) = u.at(2 * i);
260 }
261
262
263 answer.at(1, 1) = u1.dotProduct(dnx);
264 answer.at(1, 2) = u1.dotProduct(dny);
265 answer.at(2, 1) = u2.dotProduct(dnx);
266 answer.at(2, 2) = u2.dotProduct(dny);
267}
268
269void
270Quad10_2D_SUPG :: computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp)
271{
272 FloatMatrix dn;
274
275 answer.beTranspositionOf(dn);
276}
277
278
279void
280Quad10_2D_SUPG :: computeDivTauMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
281{
282 answer.resize(2, 8);
283 answer.zero();
284}
285
286
287void
288Quad10_2D_SUPG :: updateStabilizationCoeffs(TimeStep *tStep)
289{
290 double Re, norm_un, mu, mu_min, nu, norm_N, norm_N_d, norm_M_d, norm_LSIC, norm_G_c, norm_M_c, norm_N_c, t_p1, t_p2, t_p3, t_s1, t_s2, t_s3, rho;
291 FloatMatrix N, N_d, M_d, LSIC, G_c, M_c, N_c;
292 FloatArray u;
293
294 mu_min = 1;
295 rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity( integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
296 for ( GaussPoint *gp: *integrationRulesArray [ 1 ] ) {
297 mu = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial()->giveEffectiveViscosity(gp, tStep);
298 if ( mu_min > mu ) {
299 mu_min = mu;
300 }
301 }
302
303 nu = mu_min / rho;
304
305 //this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
306 this->computeVectorOfVelocities(VM_Total, tStep, u);
307
308 norm_un = u.computeNorm();
309
310 this->computeAdvectionTerm(N, tStep);
311 this->computeAdvectionDeltaTerm(N_d, tStep);
312 this->computeMassDeltaTerm(M_d, tStep);
313 this->computeLSICTerm(LSIC, tStep);
314 this->computeLinearAdvectionTerm_MC(G_c, tStep);
315 this->computeMassEpsilonTerm(M_c, tStep);
316 this->computeAdvectionEpsilonTerm(N_c, tStep);
317
318
319 norm_N = N.computeFrobeniusNorm();
320 norm_N_d = N_d.computeFrobeniusNorm();
321 norm_M_d = M_d.computeFrobeniusNorm();
322 norm_LSIC = LSIC.computeFrobeniusNorm();
323 norm_G_c = G_c.computeFrobeniusNorm();
324 norm_M_c = M_c.computeFrobeniusNorm();
325 norm_N_c = N_c.computeFrobeniusNorm();
326
327 if ( ( norm_N == 0 ) || ( norm_N_d == 0 ) || ( norm_M_d == 0 ) ) {
328 t_supg = 0;
329 } else {
330 Re = ( norm_un / nu ) * ( norm_N / norm_N_d );
331
332 t_s1 = norm_N / norm_N_d;
333
334 t_s2 = tStep->giveTimeIncrement() * ( norm_N / norm_M_d ) * 0.5;
335
336 t_s3 = t_s1 * Re;
337
338 t_supg = 1. / sqrt( 1. / ( t_s1 * t_s1 ) + 1. / ( t_s2 * t_s2 ) + 1. / ( t_s3 * t_s3 ) );
339 // t_supg = 0;
340 }
341
342 if ( norm_LSIC == 0 ) {
343 t_lsic = 0;
344 } else {
345 t_lsic = norm_N / norm_LSIC;
346
347 // t_lsic = 0;
348 }
349
350 if ( ( norm_G_c == 0 ) || ( norm_N_c == 0 ) || ( norm_M_c == 0 ) ) {
351 t_pspg = 0;
352 } else {
353 Re = ( norm_un / nu ) * ( norm_N / norm_N_d );
354
355 t_p1 = norm_G_c / norm_N_c;
356
357 t_p2 = tStep->giveTimeIncrement() * ( norm_G_c / norm_M_c ) * 0.5;
358
359 t_p3 = t_p1 * Re;
360
361 t_pspg = 1. / sqrt( 1. / ( t_p1 * t_p1 ) + 1. / ( t_p2 * t_p2 ) + 1. / ( t_p3 * t_p3 ) );
362 // t_pspg = 0;
363 }
364
365 //t_pspg = 0;
366}
367
368
369void
370Quad10_2D_SUPG :: computeAdvectionTerm(FloatMatrix &answer, TimeStep *tStep)
371{
372 FloatMatrix n, b;
373
374 answer.clear();
375
376 /* consistent part + supg stabilization term */
377 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
378 this->computeNuMatrix(n, gp);
379 this->computeUDotGradUMatrix(b, gp, tStep);
380 double dV = this->computeVolumeAround(gp);
381 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
382 answer.plusProductUnsym(n, b, rho * dV);
383 }
384}
385
386
387void
388Quad10_2D_SUPG :: computeAdvectionDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
389{
390 FloatMatrix n, b;
391
392 answer.clear();
393
394 /* consistent part + supg stabilization term */
395 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
396 this->computeNuMatrix(n, gp);
397 this->computeUDotGradUMatrix(b, gp, tStep);
398 double dV = this->computeVolumeAround(gp);
399 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
400
401 answer.plusProductUnsym(b, b, rho * dV);
402 }
403}
404
405
406void
407Quad10_2D_SUPG :: computeMassDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
408{
409 FloatMatrix n, b;
410
411 answer.clear();
412
413 /* mtrx for computing t_supg, norm of this mtrx is computed */
414 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
415 this->computeNuMatrix(n, gp);
416 this->computeUDotGradUMatrix(b, gp, tStep);
417 double dV = this->computeVolumeAround(gp);
418 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
419
420 answer.plusProductUnsym(b, n, rho * dV);
421 }
422}
423
424void
425Quad10_2D_SUPG :: computeLSICTerm(FloatMatrix &answer, TimeStep *tStep)
426{
427 FloatMatrix b;
428
429 answer.clear();
430
431 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
432 double dV = this->computeVolumeAround(gp);
433 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
434 this->computeDivUMatrix(b, gp);
435
436 answer.plusProductSymmUpper(b, b, dV * rho);
437 }
438
439 answer.symmetrized();
440}
441
442
443void
444Quad10_2D_SUPG :: computeAdvectionEpsilonTerm(FloatMatrix &answer, TimeStep *tStep)
445{
446 // to compute t_pspg
447 FloatMatrix g, b;
448
449 answer.clear();
450
451 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
452 this->computeGradPMatrix(g, gp);
453 this->computeUDotGradUMatrix(b, gp, tStep);
454 double dV = this->computeVolumeAround(gp);
455
456 answer.plusProductUnsym(g, b, dV);
457 }
458}
459
460void
461Quad10_2D_SUPG :: computeMassEpsilonTerm(FloatMatrix &answer, TimeStep *tStep)
462{
463 // to compute t_pspg
464 FloatMatrix g, n;
465
466 answer.clear();
467
468 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
469 this->computeGradPMatrix(g, gp);
470 this->computeNuMatrix(n, gp);
471 double dV = this->computeVolumeAround(gp);
472
473 answer.plusProductUnsym(g, n, dV);
474 }
475}
476
477int
478Quad10_2D_SUPG :: giveNumberOfSpatialDimensions()
479{
480 return 2;
481}
482
483
484double
485Quad10_2D_SUPG :: LS_PCS_computeF(LevelSetPCS *ls, TimeStep *tStep)
486{
487#if 0
488 FloatArray fi(3), un;
489
490 this->computeVectorOfVelocities(VM_Total, tStep, un);
491 for ( int i = 1; i <= 3; i++ ) {
492 fi.at(i) = ls->giveLevelSetDofManValue( dofManArray.at(i) );
493 }
494
495 double fix = b [ 0 ] * fi.at(1) + b [ 1 ] * fi.at(2) + b [ 2 ] * fi.at(3);
496 double fiy = c [ 0 ] * fi.at(1) + c [ 1 ] * fi.at(2) + c [ 2 ] * fi.at(3);
497 double norm = sqrt(fix * fix + fiy * fiy);
498
499 return ( 1. / 3. ) * ( fix * ( un.at(1) + un.at(3) + un.at(5) ) + fiy * ( un.at(2) + un.at(4) + un.at(6) ) ) / norm;
500
501#endif
502 return 0.0;
503}
504
505
506double
507Quad10_2D_SUPG :: LS_PCS_computeS(LevelSetPCS *ls, TimeStep *tStep)
508{
509 return 0.0;
510}
511
512
513void
514Quad10_2D_SUPG :: LS_PCS_computedN(FloatMatrix &answer)
515{ }
516
517
518void
519Quad10_2D_SUPG :: LS_PCS_computeVOFFractions(FloatArray &answer, FloatArray &fi)
520{ }
521
522
523double
524Quad10_2D_SUPG :: computeCriticalTimeStep(TimeStep *tStep)
525{
526 return 1.e6;
527}
528
529
530void
531Quad10_2D_SUPG :: NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node,
532 InternalStateType type, TimeStep *tStep)
533{
534 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
535 this->giveIPValue(answer, gp, type, tStep);
536}
537
538
539int
540Quad10_2D_SUPG :: checkConsistency()
541{
542 return SUPGElement :: checkConsistency();
543}
544
545
546void
547Quad10_2D_SUPG :: updateYourself(TimeStep *tStep)
548{
549 SUPGElement :: updateYourself(tStep);
550}
551
552int
553Quad10_2D_SUPG :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
554{
555 if ( type == IST_VOFFraction ) {
556 MaterialInterface *mi = domain->giveEngngModel()->giveMaterialInterface( domain->giveNumber() );
557 if ( mi ) {
558 FloatArray val;
560 answer = FloatArray{val.at(1)};
561 return 1;
562 } else {
563 answer = FloatArray{1.0};
564 return 1;
565 }
566 } else {
567 return SUPGElement :: giveIPValue(answer, gp, type, tStep);
568 }
569}
570
571
572void
573Quad10_2D_SUPG :: saveContext(DataStream &stream, ContextMode mode)
574{
575 SUPGElement :: saveContext(stream, mode);
576}
577
578
579
580void
581Quad10_2D_SUPG :: restoreContext(DataStream &stream, ContextMode mode)
582{
583 SUPGElement :: restoreContext(stream, mode);
584}
585
586
587double
588Quad10_2D_SUPG :: computeVolumeAround(GaussPoint *gp)
589// Returns the portion of the receiver which is attached to gp.
590{
591 double determinant, weight, volume;
592
593 determinant = fabs( this->velocityInterpolation.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
594
595
596 weight = gp->giveWeight();
597 volume = determinant * weight;
598
599 return volume;
600}
601
602#if 0
603double
604Quad10_2D_SUPG :: computeVolumeAroundPressure(FEInterpolation2d &interpol, GaussPoint *gp)
605// Returns the portion of the receiver which is attached to gp.
606{
607 double determinant, weight, volume;
608
609 determinant = fabs( interpol.giveTransformationJacobian(domain, pressureDofManArray, gp->giveNaturalCoordinates(), 0.0) );
610
611 weight = gp->giveWeight();
612 volume = determinant * weight;
613
614 return volume;
615}
616#endif
617
618Interface *
619Quad10_2D_SUPG :: giveInterface(InterfaceType interface)
620{
621 if ( interface == LevelSetPCSElementInterfaceType ) {
622 return static_cast< LevelSetPCSElementInterface * >(this);
623 } else if ( interface == ZZNodalRecoveryModelInterfaceType ) {
624 return static_cast< ZZNodalRecoveryModelInterface * >(this);
625 } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
626 return static_cast< NodalAveragingRecoveryModelInterface * >(this);
627 }
628
629 return NULL;
630}
631
632
633void
634Quad10_2D_SUPG :: giveLocalVelocityDofMap(IntArray &map)
635{
636 map.enumerate(8);
637}
638
639
640void
641Quad10_2D_SUPG :: giveLocalPressureDofMap(IntArray &map)
642{
643 map = {9};
644}
645
646
647#ifdef __OOFEG
648int
649Quad10_2D_SUPG :: giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode,
650 int node, TimeStep *tStep)
651{
652 return SUPGElement :: giveInternalStateAtNode(answer, type, mode, node, tStep);
653}
654
655void
656Quad10_2D_SUPG :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
657{
658 WCRec p [ 3 ];
659 GraphicObj *go;
660
661 if ( !gc.testElementGraphicActivity(this) ) {
662 return;
663 }
664
665 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
666 EASValsSetColor( gc.getElementColor() );
667 EASValsSetEdgeColor( gc.getElementEdgeColor() );
668 EASValsSetEdgeFlag(true);
669 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
670 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
671 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
672 p [ 0 ].z = 0.;
673 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
674 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
675 p [ 1 ].z = 0.;
676 p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
677 p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
678 p [ 2 ].z = 0.;
679
680 go = CreateTriangle3D(p);
681 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
682 EGAttachObject(go, ( EObjectP ) this);
683 EMAddGraphicsToModel(ESIModel(), go);
684}
685
686void
687Quad10_2D_SUPG :: drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
688{
689 int i, indx, result = 0;
690 WCRec p [ 3 ];
691 GraphicObj *tr;
692 FloatArray v1, v2, v3;
693 double s [ 3 ];
694
695 if ( !gc.testElementGraphicActivity(this) ) {
696 return;
697 }
698
699 EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
700
701 // if ((gc.giveIntVarMode() == ISM_local) && (gc.giveIntVarType() == IST_VOFFraction)) {
702 if ( ( gc.giveIntVarType() == IST_VOFFraction ) && ( gc.giveIntVarMode() == ISM_local ) ) {
703 Polygon matvolpoly;
704 //this->formMaterialVolumePoly(matvolpoly, NULL, temp_normal, temp_p, false);
705 EASValsSetColor( gc.getStandardSparseProfileColor() );
706 //GraphicObj *go = matvolpoly.draw(gc,true,OOFEG_VARPLOT_PATTERN_LAYER);
707 matvolpoly.draw(gc, true, OOFEG_VARPLOT_PATTERN_LAYER);
708 return;
709 }
710
711 if ( gc.giveIntVarMode() == ISM_recovered ) {
712 result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
713 result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
714 result += this->giveInternalStateAtNode(v3, gc.giveIntVarType(), gc.giveIntVarMode(), 3, tStep);
715 } else if ( gc.giveIntVarMode() == ISM_local ) {
716 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
717 result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
718 v2 = v1;
719 v3 = v1;
720 result *= 3;
721 }
722
723 if ( result != 3 ) {
724 return;
725 }
726
727 indx = gc.giveIntVarIndx();
728
729 s [ 0 ] = v1.at(indx);
730 s [ 1 ] = v2.at(indx);
731 s [ 2 ] = v3.at(indx);
732
733 EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
734
735 if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
736 for ( i = 0; i < 3; i++ ) {
737 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
738 p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
739 p [ i ].z = 0.;
740 }
741
742 //EASValsSetColor(gc.getYieldPlotColor(ratio));
743 gc.updateFringeTableMinMax(s, 3);
744 tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
745 EGWithMaskChangeAttributes(LAYER_MASK, tr);
746 EMAddGraphicsToModel(ESIModel(), tr);
747 } else if ( ( gc.getScalarAlgo() == SA_ZPROFILE ) || ( gc.getScalarAlgo() == SA_COLORZPROFILE ) ) {
748 double landScale = gc.getLandScale();
749
750 for ( i = 0; i < 3; i++ ) {
751 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
752 p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
753 p [ i ].z = s [ i ] * landScale;
754 }
755
756 if ( gc.getScalarAlgo() == SA_ZPROFILE ) {
757 EASValsSetColor( gc.getDeformedElementColor() );
758 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
759 EASValsSetFillStyle(FILL_SOLID);
760 tr = CreateTriangle3D(p);
761 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | FILL_MASK | LAYER_MASK, tr);
762 } else {
763 gc.updateFringeTableMinMax(s, 3);
764 EASValsSetFillStyle(FILL_SOLID);
765 tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
766 EGWithMaskChangeAttributes(FILL_MASK | LAYER_MASK, tr);
767 }
768
769 EMAddGraphicsToModel(ESIModel(), tr);
770 }
771}
772
773#endif
774} // end namespace oofem
#define N(a, b)
#define REGISTER_Element(class)
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
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
Definition feinterpol.C:81
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
int giveNumber() const
Definition femcmpnn.h:104
void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition fmelement.C:44
double computeNorm() const
Definition floatarray.C:861
double & at(Index i)
Definition floatarray.h:202
double dotProduct(const FloatArray &x) const
Definition floatarray.C:524
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
void beNMatrixOf(const FloatArray &n, int nsd)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
double computeFrobeniusNorm() const
double at(std::size_t i, std::size_t j) const
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
Element * giveElement()
Returns corresponding element to receiver.
Definition gausspoint.h:187
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
void enumerate(int maxVal)
Definition intarray.C:85
double giveLevelSetDofManValue(int i)
Returns level set value in specific node.
virtual void giveElementMaterialMixture(FloatArray &answer, int ielem)=0
GraphicObj * draw(oofegGraphicContext &, bool filled, int layer=OOFEG_DEBUG_LAYER)
Definition geotoolbox.C:152
int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) override
void computeMassDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
double computeVolumeAround(GaussPoint *gp) override
void computeAdvectionDeltaTerm(FloatMatrix &answer, TimeStep *tStep)
void computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp) override
void computeAdvectionEpsilonTerm(FloatMatrix &answer, TimeStep *tStep)
void computeAdvectionTerm(FloatMatrix &answer, TimeStep *tStep)
ElementDofManager pressureNode
static FEI2dQuadLin velocityInterpolation
void computeNuMatrix(FloatMatrix &answer, GaussPoint *gp) override
void computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) override
void computeLSICTerm(FloatMatrix &answer, TimeStep *tStep)
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
static FEI2dQuadConst pressureInterpolation
void computeMassEpsilonTerm(FloatMatrix &answer, TimeStep *tStep)
void computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp) override
SUPGElement2(int n, Domain *aDomain)
void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep) override
double giveTimeIncrement()
Returns solution step associated time increment.
Definition timestep.h:168
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
long ContextMode
Definition contextmode.h:43
double norm(const FloatArray &x)
InternalStateMode
Determines the mode of internal variable.
@ ZZNodalRecoveryModelInterfaceType
@ LevelSetPCSElementInterfaceType
@ NodalAveragingRecoveryModelInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_VARPLOT_PATTERN_LAYER
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
#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