OOFEM 3.0
Loading...
Searching...
No Matches
truss3dnl.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#include "parametermanager.h"
49#include "paramkey.h"
50
51
52namespace oofem {
55
56Truss3dnl :: Truss3dnl(int n, Domain *aDomain) : Truss3d(n, aDomain)
57{
59}
60
61
62void
63Truss3dnl :: initializeFrom(InputRecord &ir, int priority)
64{
65 Truss3d :: initializeFrom(ir, priority);
66 ParameterManager &ppm = this->giveDomain()->dofmanPPM;
68}
69
70
71void
72Truss3dnl :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
73{
74 FloatMatrix B, Be;
75 FloatArray vStress, vStrain, u;
76
77 // This function can be quite costly to do inside the loops when one has many slave dofs.
78 this->computeVectorOf(VM_Total, tStep, u);
79 // subtract initial displacements, if defined
82 }
83
84 // zero answer will resize accordingly when adding first contribution
85 answer.clear();
86
87 for ( auto &gp: *this->giveDefaultIntegrationRulePtr() ) {
88 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
89 this->_computeBmatrixAt(gp, B, tStep, true);
90 this->_computeBmatrixAt(gp, Be, tStep);
91 if ( useUpdatedGpRecord == 1 ) {
92 vStress = matStat->giveStressVector();
93 } else {
95 if ( !this->isActivated(tStep) ) {
96 vStrain.resize( StructuralMaterial :: giveSizeOfVoigtSymVector( gp->giveMaterialMode() ) );
97 vStrain.zero();
98 }
99 vStrain.beProductOf(Be, u);
100 // add influence of initial stress/stretch
101 double l2 = initialStretch*initialStretch;
102 vStrain.times(l2);
103 FloatArray E0(1);
104 E0.at(1) = (l2-1.)/2.;
105 vStrain.add(E0);
106 //
107 this->computeStressVector(vStress, vStrain, gp, tStep);
108 }
109
110 if ( vStress.giveSize() == 0 ) {
111 break;
112 }
113
114 // Compute nodal internal forces at nodes as f = B^T*Stress dV
115 double dV = this->computeVolumeAround(gp);
116
117 if ( vStress.giveSize() == 6 ) {
118 // It may happen that e.g. plane strain is computed
119 // using the default 3D implementation. If so,
120 // the stress needs to be reduced.
121 // (Note that no reduction will take place if
122 // the simulation is actually 3D.)
123 FloatArray stressTemp;
124 StructuralMaterial :: giveReducedSymVectorForm( stressTemp, vStress, gp->giveMaterialMode() );
125 answer.plusProduct(B, stressTemp, dV);
126 } else {
127 answer.plusProduct(B, vStress, dV);
128 }
129
130
131 // If inactive: update fields but do not give any contribution to the internal forces
132 if ( !this->isActivated(tStep) ) {
133 answer.zero();
134 return;
135 }
136 }
137}
138
139
140
141
142void
143Truss3dnl :: computeStiffnessMatrix(FloatMatrix &answer,
144 MatResponseMode rMode, TimeStep *tStep)
145{
147 bool matStiffSymmFlag = cs->isCharacteristicMtrxSymmetric(rMode);
148
149 answer.clear();
150
151 if ( !this->isActivated(tStep) ) {
152 return;
153 }
154
155 // Compute matrix from material stiffness (total stiffness for small def.) - B^T * dS/dE * B
156 if ( integrationRulesArray.size() == 1 ) {
157 FloatMatrix B, D, DB, Ksigma;
158 for ( auto &gp : *this->giveDefaultIntegrationRulePtr() ) {
159 this->_computeBmatrixAt(gp, B, tStep, true);
160 this->computeConstitutiveMatrixAt(D, rMode, gp, tStep);
161 double dV = this->computeVolumeAround(gp);
162 DB.beProductOf(D, B);
163 if ( matStiffSymmFlag ) {
164 answer.plusProductSymmUpper(B, DB, dV);
165 } else {
166 answer.plusProductUnsym(B, DB, dV);
167 }
168 this->computeInitialStressStiffness(Ksigma, rMode, gp, tStep);
169 Ksigma.times(dV);
170 answer.add(Ksigma);
171
172 }
173
174 if ( matStiffSymmFlag ) {
175 answer.symmetrized();
176 }
177 }
178}
179
180
181void
182Truss3dnl :: _computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, TimeStep *tStep, bool lin)
183{
184 FloatMatrix Bl, Bnl;
185 this->computeBlMatrixAt(gp, Bl);
186 this->computeBnlMatrixAt(gp, Bnl, tStep, lin);
187 answer = Bl;
188 answer.add(Bnl);
189}
190
191
192
193void
194Truss3dnl :: computeBlMatrixAt(GaussPoint *gp, FloatMatrix &answer)
195//
196// Returns linear part of geometrical equations of the receiver at gp.
197// Returns the linear part of the B matrix
198//
199{
200 Truss3d::computeBmatrixAt(gp, answer);
201}
202
203
204
205void
206Truss3dnl :: computeBnlMatrixAt(GaussPoint *gp, FloatMatrix &answer, TimeStep *tStep, bool lin)
207//
208// Returns linear part of geometrical equations of the receiver at gp.
209// Returns the linear part of the B matrix
210//
211{
212 FloatArray d;
213 this->computeVectorOf(VM_Total, tStep, d);
214
215 FloatMatrix Bnl, A(6,6);
216 A.at(1,1) = A.at(2,2) = A.at(3,3) = A.at(4,4) = A.at(5,5) = A.at(6,6) = 1.0;
217 A.at(1,4) = A.at(2,5) = A.at(3,6) = A.at(4,1) = A.at(5,2) = A.at(6,3) = -1.0;
218 double l0 = this->computeLength();
219 double factor = 1/l0/l0;
220 if(!lin) {
221 factor /= 2;
222 }
224 Bnl.times(factor);
225 answer.beTranspositionOf(Bnl);
226
227}
228
229
230void
231Truss3dnl :: computeInitialStressStiffness(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
232{
233 answer.resize(6,6);
234 answer.at(1,1) = answer.at(2,2) = answer.at(3,3) = answer.at(4,4) = answer.at(5,5) = answer.at(6,6) = 1.0;
235 answer.at(1,4) = answer.at(2,5) = answer.at(3,6) = answer.at(4,1) = answer.at(5,2) = answer.at(6,3) = -1.0;
236
237 FloatArray d, strain;
238 FloatMatrix B;
239 this->computeVectorOf(VM_Total, tStep, d);
240 this->_computeBmatrixAt(gp, B, tStep);
241 strain.beProductOf(B, d);
242 // add influence of initial stress/stretch
243 double l2 = initialStretch*initialStretch;
244 strain.times(l2);
245 FloatArray E0(1);
246 E0.at(1) = (l2-1.)/2;
247 strain.add(E0);
249 auto stress = this->giveStructuralCrossSection()->giveRealStress_1d(strain, gp, tStep);
250 double l0 = this->computeLength();
251 double factor = 1/l0/l0;
252 // prevent zero initial stress stiffness
253 if ( stress.at(1) == 0 ) {
254 FloatMatrix D;
255 this->computeConstitutiveMatrixAt(D, rMode, gp, tStep);
256 stress.at(1) = D.at(1,1) * 1.e-8;
257 }
258 answer.times(stress.at(1));
259 answer.times(factor);
260}
261
262
263void
264Truss3dnl :: computeGaussPoints()
265// Sets up the array of Gauss Points of the receiver.
266{
267 if ( integrationRulesArray.size() == 0 ) {
268 integrationRulesArray.resize( 1 );
269 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 2);
270 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
271 }
272}
273
274
275
276
277
278} // end namespace oofem
279
#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
Domain * giveDomain() const
Definition femcmpnn.h:97
int number
Component number.
Definition femcmpnn.h:77
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 subtract(const FloatArray &src)
Definition floatarray.C:320
void times(double s)
Definition floatarray.C:834
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)
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 beTranspositionOf(const FloatMatrix &src)
double at(std::size_t i, std::size_t j) const
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 & giveStressVector() const
Returns the const pointer to receiver's stress vector.
double computeVolumeAround(GaussPoint *gp) override
Definition truss3d.C:192
double computeLength() override
Definition truss3d.C:129
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition truss3d.C:90
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition truss3d.C:243
Truss3d(int n, Domain *d)
Definition truss3d.C:57
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss3d.C:267
static ParamKey IPK_Truss3dnl_initialStretch
[optional] Initial stretch of the truss element.
Definition truss3dnl.h:52
void computeInitialStressStiffness(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Definition truss3dnl.C:231
double initialStretch
Definition truss3dnl.h:51
void computeBnlMatrixAt(GaussPoint *gp, FloatMatrix &answer, TimeStep *tStep, bool lin=false)
Definition truss3dnl.C:206
void _computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, TimeStep *tStep, bool lin=false)
Definition truss3dnl.C:182
void computeBlMatrixAt(GaussPoint *gp, FloatMatrix &answer)
Definition truss3dnl.C:194
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)

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