OOFEM 3.0
Loading...
Searching...
No Matches
truss3d.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
37#include "fei3dlinelin.h"
38#include "node.h"
39#include "material.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "floatarray.h"
44#include "intarray.h"
45#include "mathfem.h"
46#include "classfactory.h"
47
48#ifdef __OOFEG
49 #include "oofeggraphiccontext.h"
50#endif
51
52namespace oofem {
54
56
57Truss3d::Truss3d(int n, Domain *aDomain) :
59{
61}
62
63
65
66
69{
70 if ( interface == ZZNodalRecoveryModelInterfaceType ) {
71 return static_cast< ZZNodalRecoveryModelInterface * >( this );
72 } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
73 return static_cast< NodalAveragingRecoveryModelInterface * >( this );
74 }
75
76 //OOFEM_LOG_INFO("Interface on Truss3d element not supported");
77 return NULL;
78}
79
80
81void
83{
84 auto gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
85 this->giveCrossSection()->giveIPValue(answer, gp, type, tStep);
86}
87
88
89void
91//
92// Returns linear part of geometrical equations of the receiver at gp.
93// Returns the linear part of the B matrix
94//
95{
96 FloatMatrix dN;
97 this->interp.evaldNdx(dN, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
98
99 answer.resize(1, 6);
100 answer.at(1, 1) = dN.at(1, 1);
101 answer.at(1, 2) = dN.at(1, 2);
102 answer.at(1, 3) = dN.at(1, 3);
103 answer.at(1, 4) = dN.at(2, 1);
104 answer.at(1, 5) = dN.at(2, 2);
105 answer.at(1, 6) = dN.at(2, 3);
106}
107
108
109void
111{
112 this->computeBmatrixAt(gp, answer);
113}
114
115
116void
118// Sets up the array of Gauss Points of the receiver.
119{
120 if ( integrationRulesArray.size() == 0 ) {
121 integrationRulesArray.resize(1);
122 integrationRulesArray [ 0 ] = std::make_unique< GaussIntegrationRule >(1, this, 1, 2);
124 }
125}
126
127
128double
130{
131 return this->interp.giveLength(FEIElementGeometryWrapper(this) );
132}
133
134
135void
137{
138 // computes initial stress matrix of receiver (or geometric stiffness matrix)
139 double area, NForce, contrib;
140 double l0 = this->computeLength();
141 FloatArray stress, strain;
142 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
143
144 this->computeStrainVector(strain, gp, tStep);
145 this->computeStressVector(stress, strain, gp, tStep);
146 area = this->giveCrossSection()->give(CS_Area, gp);
147 NForce = stress.at(1)*area;
148 contrib = NForce / l0;
149
150 answer.resize(6, 6);
151 answer.zero();
152
153 answer.at(2,2)=answer.at(3,3)=answer.at(5,5)=answer.at(6,6) = contrib;
154 answer.at(2,5)=answer.at(5,2)=answer.at(3,6)=answer.at(6,3) = -contrib;
155}
156
157void
159// Returns the lumped mass matrix of the receiver. This expression is
160// valid in both local and global axes.
161{
162 answer.resize(6, 6);
163 answer.zero();
164 if ( !this->isActivated(tStep) ) {
165 return;
166 }
167
168 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
169 double density = this->giveStructuralCrossSection()->give('d', gp);
170 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() * 0.5;
171 answer.at(1, 1) = halfMass;
172 answer.at(2, 2) = halfMass;
173 answer.at(3, 3) = halfMass;
174 answer.at(4, 4) = halfMass;
175 answer.at(5, 5) = halfMass;
176 answer.at(6, 6) = halfMass;
177}
178
179
180void
182// Returns the displacement interpolation matrix {N} of the receiver, eva-
183// luated at gp.
184{
185 FloatArray n;
186 this->interp.evalN(n, iLocCoord, FEIElementGeometryWrapper(this) );
187 answer.beNMatrixOf(n, 3);
188}
189
190
191double
193// Returns the length of the receiver. This method is valid only if 1
194// Gauss point is used.
195{
196 double detJ = this->interp.giveTransformationJacobian(gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
197 double weight = gp->giveWeight();
198 return detJ * weight * this->giveCrossSection()->give(CS_Area, gp);
199}
200
201
202int
204//
205// returns a unit vectors of local coordinate system at element
206// stored rowwise (mainly used by some materials with ortho and anisotrophy)
207//
208{
209 FloatArray lx, ly(3), lz;
210
211 lx.beDifferenceOf(this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
212 lx.normalize();
213
214 ly(0) = lx(1);
215 ly(1) = -lx(2);
216 ly(2) = lx(0);
217
218 // Construct orthogonal vector
219 double npn = ly.dotProduct(lx);
220 ly.add(-npn, lx);
221 ly.normalize();
222 lz.beVectorProductOf(ly, lx);
223
224 answer.resize(3, 3);
225 for ( int i = 1; i <= 3; i++ ) {
226 answer.at(1, i) = lx.at(i);
227 answer.at(2, i) = ly.at(i);
228 answer.at(3, i) = lz.at(i);
229 }
230
231 return 1;
232}
233
234
235void
237{
239}
240
241
242void
244{
245 answer = this->giveStructuralCrossSection()->giveRealStress_1d(strain, gp, tStep);
246}
247
248int
250{
251 if ( type == IST_BeamForceMomentTensor ) {
252 FloatArray stress, strain;
253 double area;
254 this->computeStrainVector(strain, gp, tStep);
255 this->computeStressVector(stress, strain, gp, tStep);
256 area = this->giveCrossSection()->give(CS_Area, gp);
257 answer.resize(6);
258 answer.at(1) = stress.at(1)*area;
259 return 1;
260 } else {
261 return NLStructuralElement :: giveIPValue(answer, gp, type, tStep);
262 }
263}
264
265
266void
267Truss3d::computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
268{
269 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_1d(rMode, gp, tStep);
270}
271
272void
274{
275 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_dPdF_1d(rMode, gp, tStep);
276}
277
278void
279Truss3d::giveDofManDofIDMask(int inode, IntArray &answer) const
280{
281 answer = { D_u, D_v, D_w };
282}
283
284
285void
286Truss3d::giveEdgeDofMapping(IntArray &answer, int iEdge) const
287{
288 /*
289 * provides dof mapping of local edge dofs (only nonzero are taken into account)
290 * to global element dofs
291 */
292
293 if ( iEdge != 1 ) {
294 OOFEM_ERROR("wrong edge number");
295 }
296
297 answer.resize(6);
298 answer.at(1) = 1;
299 answer.at(2) = 2;
300 answer.at(3) = 3;
301 answer.at(4) = 4;
302 answer.at(5) = 5;
303 answer.at(6) = 6;
304}
305
306
307double
309{
310 if ( iEdge != 1 ) { // edge between nodes 1 2
311 OOFEM_ERROR("wrong edge number");
312 }
313
314 double weight = gp->giveWeight();
315 return this->interp.giveTransformationJacobian(gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) * weight;
316}
317
318
319int
321{
322 // returns transformation matrix from
323 // edge local coordinate system
324 // to element local c.s
325 // (same as global c.s in this case)
326 //
327 // i.e. f(element local) = T * f(edge local)
328 //
329 FloatMatrix lcs;
330 this->giveLocalCoordinateSystem(lcs);
331 answer.beTranspositionOf(lcs);
332
333 return 1;
334}
335
336
337#ifdef __OOFEG
339{
340 GraphicObj *go;
341 // if (!go) { // create new one
342 WCRec p[ 2 ]; /* point */
343 if ( !gc.testElementGraphicActivity(this) ) {
344 return;
345 }
346
347 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
348 EASValsSetColor(gc.getElementColor() );
349 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
350 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
351 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
352 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
353 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
354 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
355 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
356 go = CreateLine3D(p);
357 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
358 EGAttachObject(go, ( EObjectP ) this);
359 EMAddGraphicsToModel(ESIModel(), go);
360}
361
362
364{
365 GraphicObj *go;
366 double defScale = gc.getDefScale();
367 // if (!go) { // create new one
368 WCRec p[ 2 ]; /* point */
369 if ( !gc.testElementGraphicActivity(this) ) {
370 return;
371 }
372
373 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
374 EASValsSetColor(gc.getDeformedElementColor() );
375 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
376 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
377 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
378 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
379
380 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
381 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
382 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
383 go = CreateLine3D(p);
384 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
385 EMAddGraphicsToModel(ESIModel(), go);
386}
387#endif
388} // end namespace oofem
#define REGISTER_Element(class)
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
virtual int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep)
double giveCoordinate(int i) const
Definition dofmanager.h:383
Node * giveNode(int i) const
Definition element.h:629
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Definition floatarray.C:403
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beNMatrixOf(const FloatArray &n, int nsd)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
void initializeFrom(InputRecord &ir, int priority) override
NLStructuralElement(int n, Domain *d)
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Definition node.C:234
virtual FloatMatrixF< 1, 1 > giveStiffnessMatrix_dPdF_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatArrayF< 1 > giveRealStress_1d(const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatMatrixF< 1, 1 > giveStiffnessMatrix_1d(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition truss3d.C:181
double computeVolumeAround(GaussPoint *gp) override
Definition truss3d.C:192
void giveEdgeDofMapping(IntArray &answer, int iEdge) const override
Definition truss3d.C:286
void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition truss3d.C:136
int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int, GaussPoint *gp) override
Definition truss3d.C:320
void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override
Definition truss3d.C:363
void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) override
Definition truss3d.C:82
double computeEdgeVolumeAround(GaussPoint *gp, int) override
Definition truss3d.C:308
void computeConstitutiveMatrix_dPdF_At(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss3d.C:273
FEInterpolation * giveInterpolation() const override
Definition truss3d.C:64
void initializeFrom(InputRecord &ir, int priority) override
Definition truss3d.C:236
double computeLength() override
Definition truss3d.C:129
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition truss3d.C:90
int giveLocalCoordinateSystem(FloatMatrix &answer) override
Definition truss3d.C:203
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition truss3d.C:249
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition truss3d.C:243
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition truss3d.C:158
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition truss3d.C:279
void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer) override
Definition truss3d.C:110
static FEI3dLineLin interp
Definition truss3d.h:57
Truss3d(int n, Domain *d)
Definition truss3d.C:57
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss3d.C:267
Interface * giveInterface(InterfaceType it) override
Definition truss3d.C:68
void computeGaussPoints() override
Definition truss3d.C:117
void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep) override
Definition truss3d.C:338
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
#define OOFEM_ERROR(...)
Definition error.h:79
@ CS_Area
Area.
@ ZZNodalRecoveryModelInterfaceType
@ NodalAveragingRecoveryModelInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_DEFORMED_GEOMETRY_LAYER
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_LAYER

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