OOFEM 3.0
Loading...
Searching...
No Matches
pfemelement.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 "pfemelement.h"
36#include "domain.h"
37#include "timestep.h"
38#include "node.h"
39#include "dof.h"
40#include "load.h"
41#include "boundaryload.h"
42#include "gausspoint.h"
44#include "intarray.h"
45#include "floatarray.h"
46#include "floatmatrix.h"
47#include "verbose.h"
48#include "material.h"
49#include "elementside.h"
50#include "mathfem.h"
51
52#ifdef __OOFEG
53 #include "oofeggraphiccontext.h"
54 #include "connectivitytable.h"
55#endif
56
57namespace oofem {
58PFEMElement :: PFEMElement(int n, Domain *aDomain) :
59 FMElement(n, aDomain)
60{ }
61
62
63PFEMElement :: ~PFEMElement()
64{ }
65
66void
67PFEMElement :: initializeFrom(InputRecord &ir, int priority)
68{
69 FMElement :: initializeFrom(ir, priority);
70
71 this->computeGaussPoints();
72}
73
74
75void
76PFEMElement :: giveCharacteristicMatrix(FloatMatrix &answer,
77 CharType mtrx, TimeStep *tStep)
78//
79// returns characteristics matrix of receiver according to mtrx
80//
81{
82 FloatMatrix reducedAnswer;
83 int size = this->computeNumberOfDofs();
84 answer.resize(size, size);
85 answer.zero();
86
87 if ( mtrx == LumpedMassMatrix ) {
88 this->computeDiagonalMassMtrx(reducedAnswer, tStep);
89 answer.assemble(reducedAnswer, this->giveVelocityDofMask(), this->giveVelocityDofMask());
90 } else if ( mtrx == StiffnessMatrix ) { //Tangent stiffness matrix set as default
91 this->computeStiffnessMatrix(reducedAnswer, TangentStiffness, tStep);
92 answer.assemble(reducedAnswer, this->giveVelocityDofMask(), this->giveVelocityDofMask());
93 } else if ( mtrx == TangentStiffnessMatrix ) {
94 this->computeStiffnessMatrix(reducedAnswer, TangentStiffness, tStep);
95 answer.assemble(reducedAnswer, this->giveVelocityDofMask(), this->giveVelocityDofMask());
96 } else if ( mtrx == PressureGradientMatrix ) {
97 this->computeGradientMatrix(reducedAnswer, tStep);
98 answer.assemble(reducedAnswer, this->giveVelocityDofMask(), this->givePressureDofMask());
99 } else if ( mtrx == DivergenceMatrix ) {
100 this->computeDivergenceMatrix(reducedAnswer, tStep);
101 answer.assemble(reducedAnswer, this->givePressureDofMask(), this->giveVelocityDofMask());
102 } else if ( mtrx == PressureLaplacianMatrix ) {
103 this->computePressureLaplacianMatrix(reducedAnswer, tStep);
104 answer.assemble(reducedAnswer, this->givePressureDofMask(), this->givePressureDofMask());
105 } else {
106 OOFEM_ERROR("giveCharacteristicMatrix: Unknown Type of characteristic mtrx.");
107 }
108}
109
110
111void
112PFEMElement :: giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode,
113 TimeStep *tStep)
114//
115// returns characteristics vector of receiver according to requested type
116//
117{
118 FloatArray reducedVector;
119 answer.resize(this->computeNumberOfDofs());
120 answer.zero();
121
122 if ( mtrx == ExternalForcesVector ) {
123 //this->computeForceVector(reducedVector, tStep);
124 //answer.assemble(reducedVector, this->giveVelocityDofMask());
125 } else if ( mtrx == PressureGradientVector ) {
126 FloatMatrix g;
127 this->giveCharacteristicMatrix(g, PressureGradientMatrix, tStep);
128 FloatArray p;
129 this->computeVectorOfPressures(VM_Total, tStep, reducedVector);
130 p.resize(this->computeNumberOfDofs());
131 p.assemble(reducedVector, this->givePressureDofMask());
132 answer.beProductOf(g, p);
133 } else if ( mtrx == MassVelocityVector ) {
134 FloatMatrix m;
135 FloatArray u;
136 this->giveCharacteristicMatrix(m, LumpedMassMatrix, tStep);
137 this->computeVectorOfVelocities(VM_Total, tStep, reducedVector);
138 u.resize(this->computeNumberOfDofs());
139 u.assemble(reducedVector, this->giveVelocityDofMask());
140 answer.beProductOf(m, u);
141 } else if ( mtrx == LaplaceVelocityVector ) {
142 FloatMatrix l;
143 FloatArray u;
144 this->giveCharacteristicMatrix(l, StiffnessMatrix, tStep);
145 this->computeVectorOfVelocities(VM_Total, tStep, reducedVector);
146 u.resize(this->computeNumberOfDofs());
147 u.assemble(reducedVector, this->giveVelocityDofMask());
148 answer.beProductOf(l, u);
149 } else if ( mtrx == MassAuxVelocityVector ) {
150 FloatMatrix m;
151 FloatArray u;
152 this->giveCharacteristicMatrix(m, LumpedMassMatrix, tStep);
153 this->computeVectorOf(VM_Intermediate, tStep, u);
154 answer.beProductOf(m, u);
155 } else if ( mtrx == DivergenceAuxVelocityVector ) {
156 FloatMatrix d;
157 this->giveCharacteristicMatrix(d, DivergenceMatrix, tStep);
158 FloatArray u_star;
159 this->computeVectorOf(VM_Intermediate, tStep, u_star);
160 answer.beProductOf(d, u_star);
161 } else if ( mtrx == DivergenceVelocityVector ) {
162 FloatMatrix d;
163 this->giveCharacteristicMatrix(d, DivergenceMatrix, tStep);
164 FloatArray u;
165 this->computeVectorOfVelocities(VM_Total, tStep, u);
166 answer.beProductOf(d, u);
167 } else if ( mtrx == LumpedMassMatrix ) {
168 this->computeDiagonalMassMtrx(reducedVector, tStep);
169 answer.assemble(reducedVector, this->giveVelocityDofMask());
170 } else {
171 OOFEM_ERROR("giveCharacteristicVector: Unknown Type of characteristic mtrx.");
172 }
173}
174
175
176void PFEMElement :: computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
177{
178 if ( type != ExternalForcesVector ) {
179 answer.clear();
180 return;
181 }
182 answer.resize(this->computeNumberOfDofs());
183 answer.zero();
184 FloatArray reducedVector;
185 this->computeBodyLoadVectorAt(reducedVector, load, tStep, mode);
186 answer.assemble(reducedVector, this->giveVelocityDofMask());
187}
188
189
190
191int
192PFEMElement :: checkConsistency()
193// no check at the moment
194{
195 int result = 1;
196
197 return result;
198}
199
200void
201PFEMElement :: updateInternalState(TimeStep *stepN)
202{
203 FloatArray stress;
204 // force updating strains & stresses
205 for ( auto &iRule: integrationRulesArray ) {
206 for ( GaussPoint *gp: *iRule ) {
207 computeDeviatoricStress(stress, gp, stepN);
208 }
209 }
210}
211
212void
213PFEMElement :: printOutputAt(FILE *file, TimeStep *tStep)
214// Performs end-of-step operations.
215{
216#ifdef __MPI_PARALLEL_MODE
217 fprintf( file, "element %d [%8d] :\n", this->giveNumber(), this->giveGlobalNumber() );
218#else
219 fprintf(file, "element %d :\n", number);
220#endif
221 for ( auto &iRule: integrationRulesArray ) {
222 iRule->printOutputAt(file, tStep);
223 }
224}
225
226
227#ifdef __OOFEG
228int
229PFEMElement :: giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode,
230 int node, TimeStep *atTime)
231{
232 int indx = 1;
233 Node *n = this->giveNode(node);
234
235 if ( type == IST_Velocity ) {
236 answer.resize( this->giveSpatialDimension() );
237 int dofindx;
238 if ( ( dofindx = n->findDofWithDofId(V_u) ) ) {
239 answer.at(indx++) = n->giveDof(dofindx)->giveUnknown(VM_Total, atTime);
240 }
241
242 if ( ( dofindx = n->findDofWithDofId(V_v) ) ) {
243 answer.at(indx++) = n->giveDof(dofindx)->giveUnknown(VM_Total, atTime);
244 }
245
246 if ( ( dofindx = n->findDofWithDofId(V_w) ) ) {
247 answer.at(indx++) = n->giveDof(dofindx)->giveUnknown(VM_Total, atTime);
248 }
249
250 return 1;
251 } else if ( type == IST_Pressure ) {
252 int dofindx;
253 if ( ( dofindx = n->findDofWithDofId(P_f) ) ) {
254 answer.resize(1);
255 answer.at(1) = n->giveDof(dofindx)->giveUnknown(VM_Total, atTime);
256 return 1;
257 } else {
258 return 0;
259 }
260 } else {
261 return Element :: giveInternalStateAtNode(answer, type, mode, node, atTime);
262 }
263}
264
265#endif
266
267// might be needed for OOFEG-plot
268// int
269// PFEMElement :: giveIntVarCompFullIndx(IntArray &answer, InternalStateType type)
270// {
271// if ( ( type == IST_Velocity ) ) {
272// IntArray mask;
273// int indx = 1;
274// answer.resize(3);
275// this->giveElementDofIDMask(EID_MomentBalance, mask);
276// if ( mask.findFirstIndexOf(V_u) ) {
277// answer.at(1) = indx++;
278// }
279//
280// if ( mask.findFirstIndexOf(V_v) ) {
281// answer.at(2) = indx++;
282// }
283//
284// if ( mask.findFirstIndexOf(V_w) ) {
285// answer.at(3) = indx++;
286// }
287//
288// return 1;
289// } else if ( type == IST_Pressure ) {
290// answer.resize(1);
291// answer.at(1) = 1;
292// return 1;
293// } else {
294// return Element :: giveIntVarCompFullIndx(answer, type);
295// }
296// }
297} // end namespace oofem
std::vector< Dof * >::const_iterator findDofWithDofId(DofIDItem dofID) const
Definition dofmanager.C:274
int giveGlobalNumber() const
Definition element.h:1129
Node * giveNode(int i) const
Definition element.h:629
virtual void computeGaussPoints()
Definition element.h:1255
virtual int giveSpatialDimension()
Definition element.C:1391
virtual int computeNumberOfDofs()
Definition element.h:436
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
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 resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
void giveCharacteristicMatrix(FloatMatrix &answer, CharType, TimeStep *) override
Definition pfemelement.C:76
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *atTime)=0
Calculates the stiffness matrix.
virtual const IntArray & givePressureDofMask() const =0
Returns mask of pressure Dofs.
virtual void computeDiagonalMassMtrx(FloatArray &answer, TimeStep *)=0
virtual void computeDivergenceMatrix(FloatMatrix &answer, TimeStep *atTime)=0
Calculates the velocity divergence matrix.
virtual void computePressureLaplacianMatrix(FloatMatrix &answer, TimeStep *atTime)=0
Calculates the pressure laplacian matrix.
virtual const IntArray & giveVelocityDofMask() const =0
Returns mask of velocity Dofs.
virtual void computeDeviatoricStress(FloatArray &answer, GaussPoint *gp, TimeStep *)=0
Computes deviatoric stress vector in given integration point and solution step from given total strai...
virtual void computeGradientMatrix(FloatMatrix &answer, TimeStep *atTime)=0
Calculates the pressure gradient matrix.
virtual void computeBodyLoadVectorAt(FloatArray &answer, BodyLoad *load, TimeStep *tStep, ValueModeType mode)=0
#define OOFEM_ERROR(...)
Definition error.h:79
InternalStateMode
Determines the mode of internal variable.

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