OOFEM 3.0
Loading...
Searching...
No Matches
supgelement2.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 "supgelement2.h"
36#include "domain.h"
37#include "timestep.h"
38#include "bodyload.h"
39#include "intarray.h"
40#include "floatarray.h"
41#include "floatmatrix.h"
42#include "fluidmodel.h"
44#include "fluidcrosssection.h"
45#include "integrationrule.h"
46#include "node.h"
47#include "dof.h"
48
49#ifdef __OOFEG
50 #include "oofeggraphiccontext.h"
51 #include "connectivitytable.h"
52#endif
53
54namespace oofem {
55SUPGElement2 :: SUPGElement2(int n, Domain *aDomain) :
56 SUPGElement(n, aDomain)
57{ }
58
59
60void
61SUPGElement2 :: initializeFrom(InputRecord &ir, int priority)
62{
63 SUPGElement :: initializeFrom(ir, priority);
64}
65
66
67void
68SUPGElement2 :: giveInputRecord(DynamicInputRecord &input)
69{
70 SUPGElement :: giveInputRecord(input);
71}
72
73
74void
75SUPGElement2 :: giveCharacteristicMatrix(FloatMatrix &answer,
76 CharType mtrx, TimeStep *tStep)
77//
78// returns characteristics matrix of receiver according to mtrx
79//
80{
81#if 0
82 if ( mtrx == TangentStiffnessMatrix ) {
83 // support for stokes solver
84 IntArray vloc, ploc;
86 int size = this->computeNumberOfDofs();
87 this->giveLocalVelocityDofMap(vloc);
88 this->giveLocalPressureDofMap(ploc);
89 answer.resize(size, size);
90 answer.zero();
92 answer.assemble(h, vloc);
93 this->computeDiffusionDerivativeTerm_MB(h, TangentStiffness, tStep);
94 answer.assemble(h, vloc);
95 this->computePressureTerm_MB(h, tStep);
96 answer.assemble(h, vloc, ploc);
97 this->computeLinearAdvectionTerm_MC(h, tStep);
98 answer.assemble(h, ploc, vloc);
100 answer.assemble(h, ploc, vloc);
101 this->computeDiffusionDerivativeTerm_MC(h, tStep);
102 answer.assemble(h, ploc, vloc);
103 this->computePressureTerm_MC(h, tStep);
104 answer.assemble(h, ploc);
105 this->computeLSICStabilizationTerm_MB(h, tStep);
106 answer.assemble(h, vloc);
107 } else
108#endif
109 {
110 OOFEM_ERROR("Unknown Type of characteristic mtrx.");
111 }
112}
113
114
115void
116SUPGElement2 :: giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode,
117 TimeStep *tStep)
118//
119// returns characteristics vector of receiver according to requested type
120//
121{
122 if ( mtrx == ExternalForcesVector ) {
123 // stokes flow
124 IntArray vloc, ploc;
125 FloatArray h;
126 int size = this->computeNumberOfDofs();
127 this->giveLocalVelocityDofMap(vloc);
128 this->giveLocalPressureDofMap(ploc);
129 answer.resize(size);
130 answer.zero();
131 this->computeBCRhsTerm_MB(h, tStep);
132 answer.assemble(h, vloc);
133 this->computeBCRhsTerm_MC(h, tStep);
134 answer.assemble(h, ploc);
135 } else
136#if 0
137 if ( mtrx == InternalForcesVector ) {
138 // stokes flow
139 IntArray vloc, ploc;
140 FloatArray h;
141 int size = this->computeNumberOfDofs();
142 this->giveLocalVelocityDofMap(vloc);
143 this->giveLocalPressureDofMap(ploc);
144 answer.resize(size);
145 answer.zero();
146 this->computeAdvectionTerm_MB(h, tStep);
147 answer.assemble(h, vloc);
148 this->computeAdvectionTerm_MC(h, tStep);
149 answer.assemble(h, ploc);
150 this->computeDiffusionTerm_MB(h, tStep);
151 answer.assemble(h, vloc);
152 this->computeDiffusionTerm_MC(h, tStep);
153 answer.assemble(h, ploc);
154
155 FloatMatrix m1;
156 FloatArray v, p;
157 // add lsic stabilization term
158 this->giveCharacteristicMatrix(m1, LSICStabilizationTerm_MB, tStep);
159 //m1.times( lscale / ( dscale * uscale * uscale ) );
160 this->computeVectorOfVelocities(VM_Total, tStep, v);
161 h.beProductOf(m1, v);
162 answer.assemble(h, vloc);
163 this->giveCharacteristicMatrix(m1, LinearAdvectionTerm_MC, tStep);
164 //m1.times( 1. / ( dscale * uscale ) );
165 h.beProductOf(m1, v);
166 answer.assemble(h, ploc);
167
168
169 // add pressure term
170 this->giveCharacteristicMatrix(m1, PressureTerm_MB, tStep);
171 this->computeVectorOfPressures(VM_Total, tStep, v);
172 h.beProductOf(m1, v);
173 answer.assemble(h, vloc);
174
175 // pressure term
176 this->giveCharacteristicMatrix(m1, PressureTerm_MC, tStep);
177 this->computeVectorOfPressures(VM_Total, tStep, v);
178 h.beProductOf(m1, v);
179 answer.assemble(h, ploc);
180 } else
181#endif
182 {
183 OOFEM_ERROR("Unknown Type of characteristic mtrx.");
184 }
185}
186
187
188int
189SUPGElement2 :: checkConsistency()
190//
191// check internal consistency
192// mainly tests, whether material and crossSection data
193// are safe for conversion to "Structural" versions
194//
195{
196 int result = 1;
197 return result;
198}
199
200
201void
202SUPGElement2 :: updateInternalState(TimeStep *tStep)
203{
204 FloatArray stress, eps;
205
206 // force updating strains & stresses
207 for ( auto &iRule: integrationRulesArray ) {
208 for ( auto &gp: *iRule ) {
209 computeDeviatoricStrain(eps, gp, tStep);
210 computeDeviatoricStress(stress, eps, gp, tStep);
211 }
212 }
213}
214
215
216#ifdef __OOFEG
217int
218SUPGElement2 :: giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode,
219 int node, TimeStep *tStep)
220{
221 int indx = 1;
222 Node *n = this->giveNode(node);
223
224 if ( type == IST_Velocity ) {
225 answer.resize( this->giveSpatialDimension() );
226 std::vector< Dof* >::const_iterator dofindx;
227 if ( ( dofindx = n->findDofWithDofId(V_u) ) != n->end() ) {
228 answer.at(indx++) = (*dofindx)->giveUnknown(VM_Total, tStep);
229 } else if ( ( dofindx = n->findDofWithDofId(V_v) ) != n->end() ) {
230 answer.at(indx++) = (*dofindx)->giveUnknown(VM_Total, tStep);
231 } else if ( ( dofindx = n->findDofWithDofId(V_w) ) != n->end() ) {
232 answer.at(indx++) = (*dofindx)->giveUnknown(VM_Total, tStep);
233 }
234
235 return 1;
236 } else if ( type == IST_Pressure ) {
237 auto dofindx = n->findDofWithDofId(P_f);
238 if ( dofindx != n->end() ) {
239 answer.resize(1);
240 answer.at(1) = (*dofindx)->giveUnknown(VM_Total, tStep);
241 return 1;
242 } else {
243 return 0;
244 }
245 } else {
246 return Element :: giveInternalStateAtNode(answer, type, mode, node, tStep);
247 }
248}
249
250#endif
251
252void
253SUPGElement2 :: computeAccelerationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
254{
255 FloatMatrix n, b;
256
257 answer.clear();
258
259 int rule = 2;
260 /* consistent part + supg stabilization term */
261 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
262 this->computeNuMatrix(n, gp);
263 this->computeUDotGradUMatrix( b, gp, tStep->givePreviousStep() );
264 double dV = this->computeVolumeAround(gp);
265 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
266 /* consistent part */
267 answer.plusProductUnsym(n, n, rho * dV);
268 /* supg stabilization */
269 answer.plusProductUnsym(b, n, rho * t_supg * dV);
270 }
271}
272
273void
274SUPGElement2 :: computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep)
275{
276 FloatMatrix n, b, bn;
277 FloatArray u, v;
278
279 answer.clear();
280
281 this->computeVectorOfVelocities(VM_Total, tStep, u);
282
283 int rule = 2;
284 /* consistent part + supg stabilization term */
285 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
286 this->computeNuMatrix(n, gp);
287 this->computeUDotGradUMatrix( bn, gp, tStep->givePreviousStep() );
288 this->computeUDotGradUMatrix(b, gp, tStep);
289 v.beProductOf(b, u);
290 double dV = this->computeVolumeAround(gp);
291 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
292 /* consistent part */
293 answer.plusProduct(n, v, rho * dV);
294
295 /* supg stabilization */
296 answer.plusProduct(bn, v, t_supg * rho * dV);
297 }
298}
299
300void
301SUPGElement2 :: computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep)
302{
303 FloatMatrix n, b, bn, grad_u, grad_uN;
304
305 answer.clear();
306
307 int rule = 2;
308 /* consistent part + supg stabilization term */
309 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
310 this->computeNuMatrix(n, gp);
311 this->computeUDotGradUMatrix( bn, gp, tStep->givePreviousStep() );
312 this->computeUDotGradUMatrix(b, gp, tStep);
313 double dV = this->computeVolumeAround(gp);
314 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
315
316 this->computeGradUMatrix(grad_u, gp, tStep);
317
318 /* consistent part */
319 answer.plusProductUnsym(n, b, rho * dV);
320
321 grad_uN.beProductOf(grad_u, n);
322 answer.plusProductUnsym(n, grad_uN, rho * dV);
323 /* supg stabilization */
324 answer.plusProductUnsym(bn, b, t_supg * rho * dV);
325 answer.plusProductUnsym(bn, grad_uN, t_supg * rho * dV);
326 }
327}
328
329void
330SUPGElement2 :: computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)
331{
332 FloatArray u, eps, stress, dDB_u;
333 FloatMatrix b, un_gu, dDB;
334 double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
335
336 answer.clear();
337
338 this->computeVectorOfVelocities(VM_Total, tStep, u);
339
340 int rule = 1;
341 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
342 double dV = this->computeVolumeAround(gp);
343 this->computeBMatrix(b, gp);
344 this->computeDivTauMatrix(dDB, gp, tStep);
345 this->computeUDotGradUMatrix( un_gu, gp, tStep->givePreviousStep() );
346 eps.beProductOf(b, u);
347 this->computeDeviatoricStress(stress, eps, gp, tStep);
348 dDB_u.beProductOf(dDB, u);
349 /* consistent part */
350 answer.plusProduct(b, stress, dV / Re);
351
352 /* SUPG term */
353 answer.plusProduct(un_gu, dDB_u, ( -1.0 ) * t_supg * dV);
354 }
355}
356
357void
358SUPGElement2 :: computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
359{
360 FloatMatrix _db, _d, _b, dDB, un_gu;
361 double Re = static_cast< FluidModel * >( domain->giveEngngModel() )->giveReynoldsNumber();
362
363 answer.clear();
364
365 int rule = 1;
366 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
367 double dV = this->computeVolumeAround(gp);
368 this->computeBMatrix(_b, gp);
369 this->computeTangent(_d, mode, gp, tStep);
370
371 this->computeDivTauMatrix(dDB, gp, tStep);
372 this->computeUDotGradUMatrix( un_gu, gp, tStep->givePreviousStep() );
373 /* standard term */
374 _db.beProductOf(_d, _b);
375 answer.plusProductUnsym(_b, _db, dV / Re); //answer.plusProduct (_b,_db,area);
376
377 /* SUPG term */
378 answer.plusProductUnsym( un_gu, dDB, t_supg * dV * ( -1.0 ) * ( 1. / Re ) );
379 }
380}
381
382
383void
384SUPGElement2 :: computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
385{
386 FloatMatrix gu, np, b;
387
388 answer.clear();
389
390 int rule = 1;
391 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
392 double dV = this->computeVolumeAround(gp);
393 this->computeDivUMatrix(gu, gp);
394 this->computeNpMatrix(np, gp);
395
396 /*alternative computing*/
397 //this->computeNuMatrix(gu, gp);
398 //this->computeGradPMatrix(np, gp);
399
400 /* standard term */
401 answer.plusProductUnsym(gu, np, ( -1.0 ) * dV);
402 }
403
404 for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
405 double dV = this->computeVolumeAround(gp);
406 this->computeUDotGradUMatrix( b, gp, tStep->givePreviousStep() );
407 this->computeGradPMatrix(np, gp);
408
409 // supg term
410 answer.plusProductUnsym(b, np, t_supg * dV);
411 }
412}
413void
414SUPGElement2 :: computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep)
415{
416 FloatMatrix b;
417
418 answer.clear();
419
420 int rule = 0;
421 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
422 double dV = this->computeVolumeAround(gp);
423 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
424 this->computeDivUMatrix(b, gp);
425
426 answer.plusProductSymmUpper(b, b, dV * rho * t_lsic);
427 }
428
429 answer.symmetrized();
430}
431
432void
433SUPGElement2 :: computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)
434{
435 FloatMatrix gu, np;
436
437 answer.clear();
438
439 int rule = 0;
440 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
441 double dV = this->computeVolumeAround(gp);
442 this->computeDivUMatrix(gu, gp);
443 this->computeNpMatrix(np, gp);
444
445 /* standard term */
446 answer.plusProductUnsym(np, gu, dV);
447 }
448}
449
450void
451SUPGElement2 :: computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)
452{
453 // N_epsilon (due to PSPG stabilization)
454 FloatMatrix g, b;
455 FloatArray u, v;
456
457 answer.clear();
458
459 int rule = 1;
460 this->computeVectorOfVelocities(VM_Total, tStep, u);
461
462 /* pspg stabilization term */
463 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
464 this->computeGradPMatrix(g, gp);
465 this->computeUDotGradUMatrix(b, gp, tStep);
466 v.beProductOf(b, u);
467 double dV = this->computeVolumeAround(gp);
468 answer.plusProduct(g, v, t_pspg * dV);
469 }
470}
471
472
473void
474SUPGElement2 :: computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
475{
476 FloatMatrix g, b;
477
478 answer.clear();
479
480 int rule = 1;
481 /* pspg stabilization term */
482 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
483 this->computeGradPMatrix(g, gp);
484 this->computeUDotGradUMatrix(b, gp, tStep);
485 double dV = this->computeVolumeAround(gp);
486 answer.plusProductUnsym(g, b, dV * t_pspg);
487 }
488}
489
490void
491SUPGElement2 :: computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)
492{
493 FloatMatrix dDB, g;
494
495 answer.clear();
496
497 int rule = 1;
498 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
499 double dV = this->computeVolumeAround(gp);
500 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
501
502 this->computeDivTauMatrix(dDB, gp, tStep);
503 this->computeGradPMatrix(g, gp);
504
505 answer.plusProductUnsym(g, dDB, ( -1.0 ) * dV * t_pspg / rho);
506 }
507
508 answer.zero();
509}
510
511void
512SUPGElement2 :: computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)
513{
514 answer.clear();
515
516 /*
517 * FloatMatrix dDB, _d, g;
518 * double sum, dV, coeff, rho, isd, nsd = this->giveNumberOfSpatialDimensions();
519 * GaussPoint *gp;
520 * int i,k;
521 * FloatArray u, dDB_u;
522 *
523 * this->computeVectorOfVelocities(VM_Total, tStep, u);
524 *
525 * for ( GaussPoint *gp: *this->integrationRulesArray [ 1 ] ) {
526 * dV = this->computeVolumeAround(gp);
527 * rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
528 *
529 * coeff = (-1.0) * dV * t_pspg / rho;
530 *
531 * this->computeGradPMatrix(g, gp);
532 * this->computeDivTauMatrix(dDB, gp, tStep);
533 *
534 * dDB_u.beProductOf(dDB, u);
535 *
536 * for ( i = 1; i <= pndofs; i++ ) {
537 * for ( sum = 0, isd = 1; isd <= nsd; isd++ ) {
538 * sum += g.at(isd, i) * dDB_u.at(isd);
539 * }
540 *
541 * answer.at(i) += coeff * sum;
542 * }
543 *
544 *
545 *
546 * }
547 */
548}
549
550
551void
552SUPGElement2 :: computeAccelerationTerm_MC(FloatMatrix &answer, TimeStep *tStep)
553{
554 FloatMatrix g, n;
555
556 answer.clear();
557 // pspg stabilization term: M_\epsilon term
558
559 int rule = 0;
560 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
561 this->computeGradPMatrix(g, gp);
562 this->computeNuMatrix(n, gp);
563 double dV = this->computeVolumeAround(gp);
564 answer.plusProductUnsym(g, n, dV * t_pspg);
565 }
566}
567
568void
569SUPGElement2 :: computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
570{
571 FloatMatrix g;
572
573 answer.clear();
574
575 int rule = 0;
576 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
577 this->computeGradPMatrix(g, gp);
578 double dV = this->computeVolumeAround(gp);
579 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
580 answer.plusProductSymmUpper(g, g, dV * t_pspg / rho);
581 }
582
583 answer.symmetrized();
584}
585
586
587void
588SUPGElement2 :: computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep)
589{
590 int nLoads;
591
592 answer.clear();
593
594 int rule = 0;
595 FloatArray gVector, helpLoadVector;
596 FloatMatrix b, nu;
597
598 // add body load (gravity) termms
599 nLoads = this->giveBodyLoadArray()->giveSize();
600 for ( int i = 1; i <= nLoads; i++ ) {
601 Load *load = domain->giveLoad( bodyLoadArray.at(i) );
602 bcGeomType ltype = load->giveBCGeoType();
603 if ( ( ltype == BodyLoadBGT ) && ( load->giveBCValType() == ForceLoadBVT ) ) {
604 load->computeComponentArrayAt(gVector, tStep, VM_Total);
605 if ( gVector.giveSize() ) {
606 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
607 this->computeUDotGradUMatrix( b, gp, tStep->givePreviousStep() );
608 this->computeNuMatrix(nu, gp);
609 double dV = this->computeVolumeAround(gp);
610 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
611 answer.plusProduct(b, gVector, t_supg * rho * dV);
612 answer.plusProduct(nu, gVector, rho * dV);
613 }
614 }
615 }
616 }
617
618 // integrate tractions
619 // if no traction bc applied but side marked as with traction load
620 // then zero traction is assumed !!!
621
622 // loop over boundary load array
623 nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
624 for ( int i = 1; i <= nLoads; i++ ) {
625 int n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
626 int id = boundaryLoadArray.at(i * 2);
627 Load *load = domain->giveLoad(n);
628 bcGeomType ltype = load->giveBCGeoType();
629 if ( ltype == EdgeLoadBGT ) {
630 this->computeEdgeLoadVector_MB(helpLoadVector, load, id, tStep);
631 if ( helpLoadVector.giveSize() ) {
632 answer.add(helpLoadVector);
633 }
634 } else if ( ltype == SurfaceLoadBGT ) {
635 this->computeSurfaceLoadVector_MB(helpLoadVector, load, id, tStep);
636 if ( helpLoadVector.giveSize() ) {
637 answer.add(helpLoadVector);
638 }
639 } else {
640 OOFEM_ERROR("unsupported load type class");
641 }
642 }
643}
644
645
646void
647SUPGElement2 :: computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep)
648{
649 int nLoads;
650 FloatArray gVector, helpLoadVector;
651 FloatMatrix g;
652
653 int rule = 1;
654
655 answer.clear();
656
657 nLoads = this->giveBodyLoadArray()->giveSize();
658 for ( int i = 1; i <= nLoads; i++ ) {
659 Load *load = domain->giveLoad( bodyLoadArray.at(i) );
660 bcGeomType ltype = load->giveBCGeoType();
661 if ( ( ltype == BodyLoadBGT ) && ( load->giveBCValType() == ForceLoadBVT ) ) {
662 load->computeComponentArrayAt(gVector, tStep, VM_Total);
663 if ( gVector.giveSize() ) {
664 for ( GaussPoint *gp: *this->integrationRulesArray [ rule ] ) {
665 this->computeGradPMatrix(g, gp);
666 double dV = this->computeVolumeAround(gp);
667 answer.plusProduct(g, gVector, t_pspg * dV);
668 }
669 }
670 }
671 }
672
673 // integrate tractions
674 // if no traction bc applied but side marked as with traction load
675 // then zero traction is assumed !!!
676
677 // loop over boundary load array
678 nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
679 for ( int i = 1; i <= nLoads; i++ ) {
680 int n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
681 int id = boundaryLoadArray.at(i * 2);
682 Load *load = domain->giveLoad(n);
683 bcGeomType ltype = load->giveBCGeoType();
684 if ( ltype == EdgeLoadBGT ) {
685 this->computeEdgeLoadVector_MC(helpLoadVector, load, id, tStep);
686 if ( helpLoadVector.giveSize() ) {
687 answer.add(helpLoadVector);
688 }
689 } else if ( ltype == SurfaceLoadBGT ) {
690 this->computeSurfaceLoadVector_MC(helpLoadVector, load, id, tStep);
691 if ( helpLoadVector.giveSize() ) {
692 answer.add(helpLoadVector);
693 }
694 } else {
695 OOFEM_ERROR("unsupported load type class");
696 }
697 }
698}
699
700
701void
702SUPGElement2 :: computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
703{
704 if ( type != ExternalForcesVector ) {
705 answer.clear();
706 return;
707 }
708
709 int rule = 1;
710 FloatArray un, nV;
711 FloatArray mb, mc;
712
713 this->computeVectorOfVelocities(VM_Total, tStep->givePreviousStep(), un);
714
715 if ( load->giveBCValType() == ForceLoadBVT ) {
716 FloatMatrix b, g, nu;
717 FloatArray gVector;
718 load->computeComponentArrayAt(gVector, tStep, VM_Total);
719
720 for ( auto &gp : *integrationRulesArray [ rule ] ) {
721 double dV = this->computeVolumeAround(gp);
722 double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
723
724 this->computeNuMatrix(nu, gp);
725 mb.plusProduct(nu, gVector, rho * dV);
726 if ( t_supg != 0 ) {
727 this->computeUDotGradUMatrix( b, gp, tStep->givePreviousStep() );
728 mb.plusProduct(b, gVector, t_supg * rho * dV);
729 }
730
731 if ( t_pspg != 0. ) {
732 this->computeGradPMatrix(g, gp);
733 mc.plusProduct(g, gVector, t_pspg * dV);
734 }
735 }
736 }
737
738 IntArray vloc, ploc;
739 this->giveLocalVelocityDofMap(vloc);
740 this->giveLocalPressureDofMap(ploc);
741 answer.resize(this->computeNumberOfDofs());
742 answer.zero();
743 answer.assemble(mb, vloc);
744 answer.assemble(mc, ploc);
745}
746
747void
748SUPGElement2 :: computeEdgeLoadVector_MB(FloatArray &answer, Load *load, int id, TimeStep *tStep)
749{
750 OOFEM_ERROR("not implemented");
751}
752
753void
754SUPGElement2 :: computeSurfaceLoadVector_MB(FloatArray &answer, Load *load, int id, TimeStep *tStep)
755{
756 OOFEM_ERROR("not implemented");
757}
758
759void
760SUPGElement2 :: computeEdgeLoadVector_MC(FloatArray &answer, Load *load, int id, TimeStep *tStep)
761{
762 OOFEM_ERROR("not implemented");
763}
764
765void
766SUPGElement2 :: computeSurfaceLoadVector_MC(FloatArray &answer, Load *load, int id, TimeStep *tStep)
767{
768 OOFEM_ERROR("not implemented");
769}
770
771
772void
773SUPGElement2 :: computeDeviatoricStrain(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
774{
775 FloatArray u;
776 FloatMatrix b;
777 this->computeVectorOfVelocities(VM_Total, tStep, u);
778
779 this->computeBMatrix(b, gp);
780 answer.beProductOf(b, u);
781}
782
783} // end namespace oofem
std::vector< Dof * >::const_iterator findDofWithDofId(DofIDItem dofID) const
Definition dofmanager.C:274
std::vector< Dof * >::iterator end()
Definition dofmanager.h:162
Node * giveNode(int i) const
Definition element.h:629
IntArray boundaryLoadArray
Definition element.h:147
virtual int giveSpatialDimension()
Definition element.C:1391
virtual int computeNumberOfDofs()
Definition element.h:436
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
Definition element.C:411
IntArray * giveBoundaryLoadArray()
Returns array containing load numbers of boundary loads acting on element.
Definition element.C:420
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
virtual double computeVolumeAround(GaussPoint *gp)
Definition element.h:530
IntArray bodyLoadArray
Definition element.h:147
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
Definition fmelement.C:51
void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition fmelement.C:44
void assemble(const FloatArray &fe, const IntArray &loc)
Definition floatarray.C:616
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Definition floatarray.C:288
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void add(const FloatArray &src)
Definition floatarray.C:218
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 beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void zero()
Zeroes all coefficient of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
virtual double giveReynoldsNumber()=0
virtual bcGeomType giveBCGeoType() const
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition load.C:84
void computeAdvectionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep) override
void computeLSICStabilizationTerm_MB(FloatMatrix &answer, TimeStep *tStep) override
virtual void computeEdgeLoadVector_MB(FloatArray &answer, Load *load, int id, TimeStep *tStep)
void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep) override
virtual void computeDivTauMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)=0
void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep) override
virtual void computeUDotGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)=0
void computeAdvectionDerivativeTerm_MB(FloatMatrix &answer, TimeStep *tStep) override
virtual void computeNpMatrix(FloatMatrix &answer, GaussPoint *gp)=0
void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep) override
virtual void computeGradPMatrix(FloatMatrix &answer, GaussPoint *gp)=0
void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep) override
void computeDeviatoricStrain(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override
virtual void computeBMatrix(FloatMatrix &anwer, GaussPoint *gp)=0
void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep) override
virtual void computeGradUMatrix(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeSurfaceLoadVector_MC(FloatArray &answer, Load *load, int id, TimeStep *tStep)
void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep) override
void computeAdvectionTerm_MB(FloatArray &answer, TimeStep *tStep) override
virtual void computeEdgeLoadVector_MC(FloatArray &answer, Load *load, int id, TimeStep *tStep)
virtual void computeDivUMatrix(FloatMatrix &answer, GaussPoint *gp)=0
void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *tStep) override
void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep) override
void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep) override
virtual void computeSurfaceLoadVector_MB(FloatArray &answer, Load *load, int id, TimeStep *tStep)
virtual void computeNuMatrix(FloatMatrix &answer, GaussPoint *gp)=0
void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep) override
void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep) override
virtual void giveLocalPressureDofMap(IntArray &map)
SUPGElement(int n, Domain *aDomain)
Definition supgelement.C:61
virtual void computeDeviatoricStress(FloatArray &answer, const FloatArray &eps, GaussPoint *gp, TimeStep *tStep)=0
virtual void giveLocalVelocityDofMap(IntArray &map)
virtual void computeTangent(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
TimeStep * givePreviousStep()
Returns pointer to previous solution step.
Definition timestep.C:132
#define OOFEM_ERROR(...)
Definition error.h:79
bcGeomType
Type representing the geometric character of loading.
Definition bcgeomtype.h:40
@ SurfaceLoadBGT
Distributed surface load.
Definition bcgeomtype.h:45
@ BodyLoadBGT
Distributed body load.
Definition bcgeomtype.h:43
@ EdgeLoadBGT
Distributed edge load.
Definition bcgeomtype.h:44
InternalStateMode
Determines the mode of internal variable.
@ ForceLoadBVT
Definition bcvaltype.h:43

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