OOFEM 3.0
Loading...
Searching...
No Matches
truss3dnl2.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
38#include "fei3dlinelin.h"
39#include "node.h"
40#include "material.h"
41#include "gausspoint.h"
43#include "floatmatrix.h"
44#include "floatarray.h"
45#include "intarray.h"
46#include "mathfem.h"
47#include "classfactory.h"
48
49
50namespace oofem {
52
53Truss3dnl2 :: Truss3dnl2(int n, Domain *aDomain) : Truss3d(n, aDomain)
54{
56}
57
58void
59Truss3dnl2 :: postInitialize()
60{
61 Truss3d :: postInitialize();
62 FloatArray vc1 = this-> giveCellGeometryWrapper()->giveVertexCoordinates( 1 );
63 FloatArray vc2 = this-> giveCellGeometryWrapper()->giveVertexCoordinates( 2 );
65 L = this->computeLength();
66}
67
68void
69Truss3dnl2 :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
70{
71 FloatMatrix B, Be;
72 FloatArray vStress, vStrain, u;
73
74 // This function can be quite costly to do inside the loops when one has many slave dofs.
75 this->computeVectorOf(VM_Total, tStep, u);
76 // subtract initial displacements, if defined
79 }
80
81 // zero answer will resize accordingly when adding first contribution
82 answer.clear();
83 for ( auto &gp: *this->giveDefaultIntegrationRulePtr() ) {
84 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
85 this->_computeBmatrixAt(gp, B, tStep, u);
86 if ( useUpdatedGpRecord == 1 ) {
87 vStress = matStat->givePVector();
88 } else {
90 if ( !this->isActivated(tStep) ) {
91 vStrain.resize( StructuralMaterial :: giveSizeOfVoigtSymVector( gp->giveMaterialMode() ) );
92 vStrain.zero();
93 }
94 // compute strain tensor, i.e., Biot strain
95 auto vStrain = this->_computeStrainVector(gp, u);
96 // compute stress tensor, i.e., firt Piola-Kirchhoff
97 vStress = this->giveStructuralCrossSection()->giveFirstPKStresses(vStrain, gp, tStep);
98 }
99
100 if ( vStress.giveSize() == 0 ) {
101 break;
102 }
103
104 // Compute nodal internal forces at nodes as f = B^T*Stress dV
105 double dV = this->computeVolumeAround(gp);
106
107 if ( vStress.giveSize() == 6 ) {
108 // It may happen that e.g. plane strain is computed
109 // using the default 3D implementation. If so,
110 // the stress needs to be reduced.
111 // (Note that no reduction will take place if
112 // the simulation is actually 3D.)
113 FloatArray stressTemp;
114 StructuralMaterial :: giveReducedSymVectorForm( stressTemp, vStress, gp->giveMaterialMode() );
115 answer.plusProduct(B, stressTemp, dV);
116 } else {
117 answer.plusProduct(B, vStress, dV);
118 }
119
120
121 // If inactive: update fields but do not give any contribution to the internal forces
122 if ( !this->isActivated(tStep) ) {
123 answer.zero();
124 return;
125 }
126 }
127}
128
129
130double
131Truss3dnl2 :: computeLength()
132{
133 FloatArray X12;
134 X12.beProductOf(this->givePmatrix(), X);
135 return X12.computeNorm();
136}
137
138double
139Truss3dnl2 :: computeDeformedLength(const FloatArray &d)
140{
141 FloatArray x12, x(X);
142 x.add(d);
143 x12.beProductOf(this->givePmatrix(), x);
144 return x12.computeNorm();
145}
146
147
148
149
150
152Truss3dnl2 :: _computeStrainVector(GaussPoint *gp, const FloatArray &d)
153{
154 FloatArray answer(1);
155 auto l = this->computeDeformedLength(d);
156 answer.at(1) = l/this->computeLength();
157 return answer;
158}
159
160
161void
162Truss3dnl2 :: computeStiffnessMatrix(FloatMatrix &answer,
163 MatResponseMode rMode, TimeStep *tStep)
164{
166 bool matStiffSymmFlag = cs->isCharacteristicMtrxSymmetric(rMode);
167
168 answer.clear();
169
170 if ( !this->isActivated(tStep) ) {
171 return;
172 }
173 FloatArray u;
174 this->computeVectorOf(VM_Total, tStep, u);
175
176 // Compute matrix from material stiffness (total stiffness for small def.) - B^T * dS/dE * B
177 if ( integrationRulesArray.size() == 1 ) {
178 FloatMatrix B, D, DB, Ksigma;
179 for ( auto &gp : *this->giveDefaultIntegrationRulePtr() ) {
180 this->_computeBmatrixAt(gp, B, tStep, u);
181 this->computeConstitutiveMatrixAt(D, rMode, gp, tStep);
182 double dV = this->computeVolumeAround(gp);
183 DB.beProductOf(D, B);
184 if ( matStiffSymmFlag ) {
185 answer.plusProductSymmUpper(B, DB, dV);
186 } else {
187 answer.plusProductUnsym(B, DB, dV);
188 }
189
190 this->computeInitialStressStiffness(Ksigma, rMode, gp, tStep,B,u);
191 Ksigma.times(dV);
192 answer.add(Ksigma);
193 }
194
195 if ( matStiffSymmFlag ) {
196 answer.symmetrized();
197 }
198 }
199}
200
201
202void
203Truss3dnl2 :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
204{
205 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_dPdF_1d(rMode, gp, tStep);
206}
207
208
209
210void
211Truss3dnl2 :: _computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, TimeStep *tStep, const FloatArray &d)
212{
213 double L = computeLength();
214 double l = computeDeformedLength(d);
215 FloatMatrixF<6,6> A = this->giveAmatrix();
216 FloatArray x(X);
217 x.add(d);
219 answer.times(1./l/L);
220}
221
222
224Truss3dnl2 :: giveAmatrix()
225{
226 FloatMatrix A;
227 A = FloatMatrix({{1, 0, 0, -1, 0, 0}, {0, 1, 0, 0, -1, 0}, {0, 0, 1, 0, 0, -1}, {-1, 0, 0, 1, 0, 0}, {0, -1, 0, 0, 1, 0}, {0, 0, -1, 0, 0, 1}});
228 return A;
229}
230
232Truss3dnl2 :: givePmatrix()
233{
234 FloatMatrix P;
235 P = FloatMatrix({{1,0,0},{0,1,0},{0,0,1},{-1,0,0},{0,-1,0},{0,0,-1}});
236 return P;
237}
238
239
240
241
242
243void
244Truss3dnl2 :: computeInitialStressStiffness(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep, const FloatMatrix &B, const FloatArray &d)
245//
246{
247
248 FloatMatrix BB, A = this->giveAmatrix();
249 double L = computeLength();
250 double l = computeDeformedLength(d);
251 FloatArray x(X);
252 FloatMatrix xx, Axx, AxxA;
253 x.add(d);
254
256 Axx.beProductOf(A,xx);
257 AxxA.beProductOf(Axx,A);
258 AxxA.times(1./l/l);
259
260 answer = A;
261 answer.subtract(AxxA);
262 answer.times(1./l/L);
263 auto stress = this->giveStructuralCrossSection()->giveFirstPKStresses(this->_computeStrainVector(gp, d), gp, tStep);
264
265 // prevent zero initial stress stiffness
266 if ( stress.at(1) == 0 ) {
267 FloatMatrix D;
268 this->computeConstitutiveMatrixAt(D, rMode, gp, tStep);
269 stress.at(1) = D.at(1,1) * initStressFactor;
270
271 }
272 answer.times(stress.at(1));
273}
274
275
276
277void
278Truss3dnl2 :: computeGaussPoints()
279// Sets up the array of Gauss Points of the receiver.
280{
281 if ( integrationRulesArray.size() == 0 ) {
282 integrationRulesArray.resize( 1 );
283 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 2);
284 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
285 }
286}
287
288
290Truss3dnl2 :: giveCellGeometryWrapper()
291{
292 if ( !cellGeometryWrapper ) {
294 }
295
296 return cellGeometryWrapper;
297}
298
299} // end namespace oofem
300
#define REGISTER_Element(class)
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
double computeNorm() const
Definition floatarray.C:861
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
static FloatArray fromConcatenated(std::initializer_list< FloatArray > ini)
Definition floatarray.C:146
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void add(const FloatArray &src)
Definition floatarray.C:218
void subtract(const FloatArray &src)
Definition floatarray.C:320
void times(double f)
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
static FloatMatrix fromArray(const FloatArray &vector, bool transpose=false)
void add(const FloatMatrix &a)
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void beProductTOf(const FloatMatrix &a, const FloatMatrix &b)
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
double at(std::size_t i, std::size_t j) const
void subtract(const FloatMatrix &a)
void beTProductOf(const FloatMatrix &a, const FloatMatrix &b)
bool isCharacteristicMtrxSymmetric(MatResponseMode mode) const override=0
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
const FloatArray & givePVector() const
Returns the const pointer to receiver's first Piola-Kirchhoff stress vector.
double computeVolumeAround(GaussPoint *gp) override
Definition truss3d.C:192
Truss3d(int n, Domain *d)
Definition truss3d.C:57
FloatMatrixF< 3, 6 > givePmatrix()
Definition truss3dnl2.C:232
FloatMatrixF< 6, 6 > giveAmatrix()
Definition truss3dnl2.C:224
double computeDeformedLength(const FloatArray &d)
Definition truss3dnl2.C:139
virtual FEICellGeometry * giveCellGeometryWrapper()
Definition truss3dnl2.C:290
FEICellGeometry * cellGeometryWrapper
Definition truss3dnl2.h:52
double initStressFactor
Definition truss3dnl2.h:58
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss3dnl2.C:203
double computeLength() override
Definition truss3dnl2.C:131
void _computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, TimeStep *tStep, const FloatArray &u)
Definition truss3dnl2.C:211
FloatArray _computeStrainVector(GaussPoint *gp, const FloatArray &u)
Definition truss3dnl2.C:152
void computeInitialStressStiffness(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep, const FloatMatrix &B, const FloatArray &d)
Definition truss3dnl2.C:244

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