OOFEM 3.0
Loading...
Searching...
No Matches
supgelement.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 "supgelement.h"
36#include "domain.h"
37#include "load.h"
38#include "gausspoint.h"
39#include "intarray.h"
40#include "floatarray.h"
41#include "floatmatrix.h"
43#include "fluidcrosssection.h"
44#include "dynamicinputrecord.h"
45#include "engngm.h"
46#include "node.h"
47#include "dof.h"
48#include "paramkey.h"
49
50#ifdef __OOFEG
51 #include "oofeggraphiccontext.h"
52 #include "connectivitytable.h"
53#endif
54
55namespace oofem {
56
59
60
61SUPGElement :: SUPGElement(int n, Domain *aDomain) :
62 FMElement(n, aDomain)
63{ }
64
65
66void
67SUPGElement :: initializeFrom(InputRecord &ir, int priority)
68{
69 FMElement :: initializeFrom(ir, priority);
70
71 ParameterManager &ppm = this->giveDomain()->elementPPM;
72 bool flag;
74 if (flag && !boundarySides.isEmpty()) {
76 if (!flag) {
77 OOFEM_ERROR("Boundary codes not provided for element %d", this->giveNumber());
78 }
79 }
80}
81
82
83void
84SUPGElement :: giveInputRecord(DynamicInputRecord &input)
85{
86 FMElement :: giveInputRecord(input);
87 if ( !boundarySides.isEmpty() ) {
88 input.setField(this->boundarySides, IPK_SUPGElement_bsides.getNameCStr());
89 input.setField(this->boundaryCodes, IPK_SUPGElement_bcodes.getNameCStr());
90 }
91}
92
93
94void
95SUPGElement :: giveCharacteristicMatrix(FloatMatrix &answer,
96 CharType type, TimeStep *tStep)
97//
98// returns characteristics matrix of receiver according to mtrx
99//
100{
101
102 if ( type == TangentStiffnessMatrix ) {
103 // stokes flow only
104 double dscale = this->giveDomain()->giveEngngModel()->giveVariableScale(VST_Density);
105 double uscale = this->giveDomain()->giveEngngModel()->giveVariableScale(VST_Velocity);
106
107 IntArray vloc, ploc;
108 FloatMatrix h;
109 int size = this->computeNumberOfDofs();
110 this->giveLocalVelocityDofMap(vloc);
111 this->giveLocalPressureDofMap(ploc);
112 answer.resize(size, size);
113 answer.zero();
114
115 //this->computeAccelerationTerm_MB(h, tStep);
116 //answer.assemble(h, vloc);
117 this->computeDiffusionDerivativeTerm_MB(h, TangentStiffness, tStep);
118 answer.assemble(h, vloc);
119 this->computePressureTerm_MB(h, tStep);
120 answer.assemble(h, vloc, ploc);
121 //this->computeLSICStabilizationTerm_MB(h, tStep);
122 //h.times( alpha * tStep->giveTimeIncrement() * lscale / ( dscale * uscale * uscale ) );
123 //answer.assemble(h, vloc);
124 this->computeBCLhsTerm_MB(h, tStep);
125 if ( h.isNotEmpty() ) {
126 answer.assemble(h, vloc);
127 }
128
129 this->computeBCLhsPressureTerm_MB(h, tStep);
130 if ( h.isNotEmpty() ) {
131 answer.assemble(h, vloc, ploc);
132 }
133
134 // conservation eq part
135 this->computeLinearAdvectionTerm_MC(h, tStep);
136 h.times( 1.0 / ( dscale * uscale ) );
137 answer.assemble(h, ploc, vloc);
138 this->computeBCLhsPressureTerm_MC(h, tStep);
139 if ( h.isNotEmpty() ) {
140 answer.assemble(h, ploc, vloc);
141 }
142
143 this->computeDiffusionDerivativeTerm_MC(h, tStep);
144 answer.assemble(h, ploc, vloc);
145 this->computePressureTerm_MC(h, tStep);
146 answer.assemble(h, ploc);
147 } else {
148 OOFEM_ERROR("giveCharacteristicMatrix: Unknown Type of characteristic mtrx.");
149 }
150
151}
152
153
154void
155SUPGElement :: giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode,
156 TimeStep *tStep)
157//
158// returns characteristics vector of receiver according to requested type
159//
160{
161 if ( mtrx == ExternalForcesVector ) {
162 answer.clear();
163#if 0
164 // assembled from assembler from loads directly
165 // stokes flow
166 IntArray vloc, ploc;
167 FloatArray h;
168 int size = this->computeNumberOfDofs();
169 this->giveLocalVelocityDofMap(vloc);
170 this->giveLocalPressureDofMap(ploc);
171 answer.resize(size);
172 answer.zero();
173 this->computeBCRhsTerm_MB(h, tStep);
174 answer.assemble(h, vloc);
175 this->computeBCRhsTerm_MC(h, tStep);
176 answer.assemble(h, ploc);
177#endif
178 }
179
180#if 1
181 else if ( mtrx == InternalForcesVector ) {
182 // stokes flow
183 IntArray vloc, ploc;
184 FloatArray h;
185 int size = this->computeNumberOfDofs();
186 this->giveLocalVelocityDofMap(vloc);
187 this->giveLocalPressureDofMap(ploc);
188 answer.resize(size);
189 answer.zero();
190 //this->computeAdvectionTerm_MB(h, tStep); answer.assemble(h, vloc);
191 this->computeAdvectionTerm_MC(h, tStep);
192 answer.assemble(h, ploc);
193 this->computeDiffusionTerm_MB(h, tStep);
194 answer.assemble(h, vloc);
195 this->computeDiffusionTerm_MC(h, tStep);
196 answer.assemble(h, ploc);
197
198 FloatMatrix m1;
199 FloatArray v;
200 // add lsic stabilization term
201 //this->giveCharacteristicMatrix(m1, LSICStabilizationTerm_MB, tStep);
202 //m1.times( lscale / ( dscale * uscale * uscale ) );
203 this->computeVectorOfVelocities(VM_Total, tStep, v);
204 //h.beProductOf(m1, v);
205 //answer.assemble(h, vloc);
206 this->computeLinearAdvectionTerm_MC(m1, tStep);
207 //m1.times( 1. / ( dscale * uscale ) );
208 h.beProductOf(m1, v);
209 answer.assemble(h, ploc);
210
211 // add pressure term
212 this->computePressureTerm_MB(m1, tStep);
213 this->computeVectorOfPressures(VM_Total, tStep, v);
214 h.beProductOf(m1, v);
215 answer.assemble(h, vloc);
216
217 // pressure term
218 this->computePressureTerm_MC(m1, tStep);
219 this->computeVectorOfPressures(VM_Total, tStep, v);
220 h.beProductOf(m1, v);
221 answer.assemble(h, ploc);
222 }
223#endif
224 else {
225 OOFEM_ERROR("Unknown Type of characteristic mtrx.");
226 }
227}
228
229
230void
231SUPGElement :: computeBCLhsTerm_MB(FloatMatrix &answer, TimeStep *tStep)
232{
233 bcType boundarytype;
234 int nLoads = 0;
235 //bcType loadtype;
236 FloatMatrix helpMatrix;
237 // loop over boundary load array
238
239 answer.clear();
240
241 nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
242 if ( nLoads ) {
243 for ( int i = 1; i <= nLoads; i++ ) {
244 int n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
245 int side = boundaryLoadArray.at(i * 2);
246 Load *load = domain->giveLoad(n);
247 boundarytype = load->giveType();
248 if ( boundarytype == SlipWithFriction ) {
249 this->computeSlipWithFrictionBCTerm_MB(helpMatrix, load, side, tStep);
250 answer.add(helpMatrix);
251 } else if ( boundarytype == PenetrationWithResistance ) {
252 this->computePenetrationWithResistanceBCTerm_MB(helpMatrix, load, side, tStep);
253 answer.add(helpMatrix);
254 } else {
255 // OOFEM_ERROR("unsupported load type class");
256 }
257 }
258 }
259
260 nLoads = this->giveBodyLoadArray()->giveSize();
261
262 if ( nLoads ) {
263 bcGeomType ltype;
264 for ( int i = 1; i <= nLoads; i++ ) {
265 Load *load = domain->giveLoad( bodyLoadArray.at(i) );
266 ltype = load->giveBCGeoType();
267 if ( ( ltype == BodyLoadBGT ) && ( load->giveBCValType() == ReinforceBVT ) ) {
268 this->computeHomogenizedReinforceTerm_MB(helpMatrix, load, tStep);
269 answer.add(helpMatrix);
270 }
271 }
272 }
273}
274
275void
276SUPGElement :: computeBCLhsPressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
277{
278 bcType boundarytype;
279 int nLoads = 0;
280 //bcType loadtype;
281 FloatMatrix helpMatrix;
282 // loop over boundary load array
283 answer.clear();
284
285 nLoads = this->giveBoundaryLoadArray()->giveSize() / 2;
286
287 if ( nLoads ) {
288 for ( int i = 1; i <= nLoads; i++ ) {
289 int n = boundaryLoadArray.at(1 + ( i - 1 ) * 2);
290 int side = boundaryLoadArray.at(i * 2);
291 Load *load = domain->giveLoad(n);
292 boundarytype = load->giveType();
293 if ( boundarytype == OutFlowBC ) {
294 this->computeOutFlowBCTerm_MB(helpMatrix, side, tStep);
295 answer.add(helpMatrix);
296 } else {
297 //_warning("computeForceLoadVector : unsupported load type class");
298 }
299 }
300 }
301}
302
303void
304SUPGElement :: computeBCLhsPressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
305{
306 int nLoads = 0;
307 //bcType loadtype;
308 FloatMatrix helpMatrix;
309
310 nLoads = this->giveBodyLoadArray()->giveSize();
311 answer.clear();
312 if ( nLoads ) {
313 bcGeomType ltype;
314 for ( int i = 1; i <= nLoads; i++ ) {
315 Load *load = domain->giveLoad( bodyLoadArray.at(i) );
316 ltype = load->giveBCGeoType();
317 if ( ( ltype == BodyLoadBGT ) && ( load->giveBCValType() == ReinforceBVT ) ) {
318 this->computeHomogenizedReinforceTerm_MC(helpMatrix, load, tStep);
319 answer.add(helpMatrix);
320 }
321 }
322 }
323}
324
325
326int
327SUPGElement :: checkConsistency()
328//
329// check internal consistency
330// mainly tests, whether material and crossSection data
331// are safe for conversion to "Structural" versions
332//
333{
334 int result = 1;
335 return result;
336}
337
338void
339SUPGElement :: updateInternalState(TimeStep *tStep)
340{
341 FloatArray stress, eps;
342
343 // force updating strains & stresses
344 for ( auto &iRule: integrationRulesArray ) {
345 for ( auto &gp: *iRule ) {
346 computeDeviatoricStrain(eps, gp, tStep);
347 computeDeviatoricStress(stress, eps, gp, tStep);
348 }
349 }
350}
351
352
353#ifdef __OOFEG
354int
355SUPGElement :: giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode,
356 int node, TimeStep *tStep)
357{
358 int indx = 1;
359 Node *n = this->giveNode(node);
360
361 if ( type == IST_Velocity ) {
362 answer.resize( this->giveSpatialDimension() );
363 std::vector< Dof* >::const_iterator dofindx;
364 if ( ( dofindx = n->findDofWithDofId(V_u) ) != n->end() ) {
365 answer.at(indx++) = (*dofindx)->giveUnknown(VM_Total, tStep);
366 } else if ( ( dofindx = n->findDofWithDofId(V_v) ) != n->end() ) {
367 answer.at(indx++) = (*dofindx)->giveUnknown(VM_Total, tStep);
368 } else if ( ( dofindx = n->findDofWithDofId(V_w) ) != n->end() ) {
369 answer.at(indx++) = (*dofindx)->giveUnknown(VM_Total, tStep);
370 }
371
372 return 1;
373 } else if ( type == IST_Pressure ) {
374 auto dofindx = n->findDofWithDofId(P_f);
375 if ( dofindx != n->end() ) {
376 answer.resize(1);
377 answer.at(1) = (*dofindx)->giveUnknown(VM_Total, tStep);
378 return 1;
379 } else {
380 return 0;
381 }
382 } else {
383 return Element :: giveInternalStateAtNode(answer, type, mode, node, tStep);
384 }
385}
386
387#endif
388
389} // 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
void setField(int item, InputFieldType id)
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
IntArray bodyLoadArray
Definition element.h:147
Domain * giveDomain() const
Definition femcmpnn.h:97
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
int number
Component number.
Definition femcmpnn.h:77
int giveNumber() const
Definition femcmpnn.h:104
void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
Definition fmelement.C:51
void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition fmelement.C:44
FMElement(int n, Domain *aDomain)
Definition fmelement.C:38
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 zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void times(double f)
void add(const FloatMatrix &a)
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 zero()
Zeroes all coefficient of receiver.
bool isNotEmpty() const
Tests for empty matrix.
void assemble(const FloatMatrix &src, const IntArray &loc)
virtual bcGeomType giveBCGeoType() const
static ParamKey IPK_SUPGElement_bsides
Definition supgelement.h:67
virtual void computeDiffusionDerivativeTerm_MB(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)=0
virtual void giveLocalPressureDofMap(IntArray &map)
virtual void computeBCLhsPressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)
virtual void computeDeviatoricStrain(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeAdvectionTerm_MC(FloatArray &answer, TimeStep *tStep)=0
virtual void computeHomogenizedReinforceTerm_MC(FloatMatrix &answer, Load *load, TimeStep *tStep)
virtual void computeDiffusionTerm_MB(FloatArray &answer, TimeStep *tStep)=0
virtual void computeDiffusionDerivativeTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeLinearAdvectionTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeBCLhsTerm_MB(FloatMatrix &answer, TimeStep *tStep)
virtual void computeBCLhsPressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)
IntArray boundaryCodes
Boundary sides codes.
Definition supgelement.h:59
virtual void computeDeviatoricStress(FloatArray &answer, const FloatArray &eps, GaussPoint *gp, TimeStep *tStep)=0
virtual void computeDiffusionTerm_MC(FloatArray &answer, TimeStep *tStep)=0
static ParamKey IPK_SUPGElement_bcodes
Definition supgelement.h:68
virtual void giveLocalVelocityDofMap(IntArray &map)
virtual void computeBCRhsTerm_MB(FloatArray &answer, TimeStep *tStep)=0
virtual void computePressureTerm_MC(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computePressureTerm_MB(FloatMatrix &answer, TimeStep *tStep)=0
virtual void computeBCRhsTerm_MC(FloatArray &answer, TimeStep *tStep)=0
IntArray boundarySides
Array of boundary sides.
Definition supgelement.h:57
virtual void computeSlipWithFrictionBCTerm_MB(FloatMatrix &answer, Load *load, int side, TimeStep *tStep)
virtual void computeHomogenizedReinforceTerm_MB(FloatMatrix &answer, Load *load, TimeStep *tStep)
virtual void computePenetrationWithResistanceBCTerm_MB(FloatMatrix &answer, Load *load, int side, TimeStep *tStep)
virtual void computeOutFlowBCTerm_MB(FloatMatrix &answer, int side, TimeStep *tStep)
#define OOFEM_ERROR(...)
Definition error.h:79
@ VST_Velocity
@ VST_Density
bcGeomType
Type representing the geometric character of loading.
Definition bcgeomtype.h:40
@ BodyLoadBGT
Distributed body load.
Definition bcgeomtype.h:43
InternalStateMode
Determines the mode of internal variable.
@ ReinforceBVT
Definition bcvaltype.h:49
bcType
Type representing the type of bc.
Definition bctype.h:40
@ SlipWithFriction
Definition bctype.h:45
@ PenetrationWithResistance
Definition bctype.h:46
@ OutFlowBC
Definition bctype.h:47
#define PM_UPDATE_PARAMETER_AND_REPORT(_val, _pm, _ir, _componentnum, _paramkey, _prio, _flag)

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