OOFEM 3.0
Loading...
Searching...
No Matches
truss2d.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 "node.h"
38#include "material.h"
39#include "gausspoint.h"
41#include "floatmatrix.h"
42#include "floatarray.h"
43#include "intarray.h"
44#include "mathfem.h"
45#include "classfactory.h"
46#include "fei2dlinelin.h"
47#include "parametermanager.h"
48#include "paramkey.h"
49
50#ifdef __OOFEG
51 #include "oofeggraphiccontext.h"
52#endif
53
54namespace oofem {
57
58
60 FEI2dLineLin(1, 2),
61 FEI2dLineLin(2, 3) };
62
63Truss2d::Truss2d(int n, Domain *aDomain) :
64 NLStructuralElement(n, aDomain)
65{
67 length = 0.;
68 pitch = 10.; // a dummy value
69 cs_mode = 0;
70}
71
72
73void
75//
76// Returns linear part of geometrical equations of the receiver at gp.
77// Returns the linear part of the B matrix
78//
79{
80 // eps = B*a = [-cos -sin cos sin]/l *[u1 v1 u2 v2]^t
81 // cos = (x2 - x1) / l
82 // sin = (z2 - z1) / l
83
84 double l, x1, x2, z1, z2;
85 //determine in which plane the truss is defined
86 int c1 = 0, c2 = 0;
87 resolveCoordIndices(c1, c2);
88
89 x1 = this->giveNode(1)->giveCoordinate(c1);
90 z1 = this->giveNode(1)->giveCoordinate(c2);
91 x2 = this->giveNode(2)->giveCoordinate(c1);
92 z2 = this->giveNode(2)->giveCoordinate(c2);
93
94 answer.resize(1, 4);
95
96 answer.at(1, 1) = x1 - x2;
97 answer.at(1, 2) = z1 - z2;
98 answer.at(1, 3) = x2 - x1;
99 answer.at(1, 4) = z2 - z1;
100
101 l = this->computeLength();
102 answer.times(1.0 / l / l);
103}
104
105
106void
108{
109 // Will be the same as the regular B-matrix
110 this->computeBmatrixAt(gp, answer);
111}
112
113
114
116// Sets up the array of Gauss Points of the receiver.
117{
118 if ( integrationRulesArray.size() == 0 ) {
119 integrationRulesArray.resize(1);
120 integrationRulesArray [ 0 ] = std::make_unique< GaussIntegrationRule >(1, this, 1, 2);
122 }
123}
124
125
126void
128// Returns the lumped mass matrix of the receiver. This expression is
129// valid in both local and global axes.
130{
131 answer.resize(4, 4);
132 answer.zero();
133 if ( !isActivated(tStep) ) {
134 return;
135 }
136
137 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
138 double density = this->giveStructuralCrossSection()->give('d', gp);
139 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() * 0.5;
140 answer.at(1, 1) = halfMass;
141 answer.at(2, 2) = halfMass;
142 answer.at(3, 3) = halfMass;
143 answer.at(4, 4) = halfMass;
144}
145
146
147void
149// Returns the displacement interpolation matrix {N} of the receiver, eva-
150// luated at gp.
151{
152 double ksi, n1, n2;
153
154 ksi = iLocCoord.at(1);
155 n1 = ( 1. - ksi ) * 0.5;
156 n2 = ( 1. + ksi ) * 0.5;
157
158 answer.resize(2, 4);
159 answer.zero();
160 answer.at(1, 1) = n1;
161 answer.at(1, 3) = n2;
162 answer.at(2, 2) = n1;
163 answer.at(2, 4) = n2;
164}
165
166
167int
169{
170 //determine in which plane the truss is defined
171 int c1 = 0, c2 = 0;
172 resolveCoordIndices(c1, c2);
173
174 double ksi, n1, n2;
175
176 ksi = lcoords.at(1);
177 n1 = ( 1. - ksi ) * 0.5;
178 n2 = ( 1. + ksi ) * 0.5;
179
180 answer.resize(3);
181 answer.zero();
182 answer.at(c1) = n1 * this->giveNode(1)->giveCoordinate(c1) + n2 * this->giveNode(2)->giveCoordinate(c1);
183 answer.at(c2) = n1 * this->giveNode(1)->giveCoordinate(c2) + n2 * this->giveNode(2)->giveCoordinate(c2);
184
185 return 1;
186}
187
188
190// Returns the length of the receiver. This method is valid only if 1
191// Gauss point is used.
192{
193 double weight = gp->giveWeight();
194 return 0.5 * this->computeLength() * weight * this->giveCrossSection()->give(CS_Area, gp);
195}
196
197
199// Returns the length of the receiver.
200{
201 //determine in which plane the truss is defined
202 int c1 = 0, c2 = 0;
203 resolveCoordIndices(c1, c2);
204
205 double dx, dz;
206 Node *nodeA, *nodeB;
207
208 if ( length == 0. ) {
209 nodeA = this->giveNode(1);
210 nodeB = this->giveNode(2);
211 dx = nodeB->giveCoordinate(c1) - nodeA->giveCoordinate(c1);
212 dz = nodeB->giveCoordinate(c2) - nodeA->giveCoordinate(c2);
213 length = sqrt(dx * dx + dz * dz);
214 }
215
216 return length;
217}
218
219
221// Returns the pitch of the receiver.
222{
223 double xA, xB, zA, zB;
224 Node *nodeA, *nodeB;
225 //determine in which plane the truss is defined
226 int c1 = 0, c2 = 0;
227 resolveCoordIndices(c1, c2);
228
229 if ( pitch == 10. ) { // 10. : dummy initialization value
230 nodeA = this->giveNode(1);
231 nodeB = this->giveNode(2);
232 xA = nodeA->giveCoordinate(c1);
233 xB = nodeB->giveCoordinate(c1);
234 zA = nodeA->giveCoordinate(c2);
235 zB = nodeB->giveCoordinate(c2);
236 pitch = atan2(zB - zA, xB - xA);
237 }
238
239 return pitch;
240}
241
242
243int
245//
246// returns a unit vectors of local coordinate system at element
247// stored rowwise (mainly used by some materials with ortho and anisotrophy)
248//
249{
250 double sine, cosine;
251
252 answer.resize(3, 3);
253 answer.zero();
254
255 sine = sin(this->givePitch() );
256 cosine = cos(pitch);
257
258 answer.at(1, 1) = cosine;
259 answer.at(1, 2) = sine;
260 answer.at(2, 1) = -sine;
261 answer.at(2, 2) = cosine;
262 answer.at(3, 3) = 1.0;
263
264 return 1;
265}
266
267
268void
270{
271 if ( cs_mode == 0 ) {
272 //xz-plane
273 c1 = 1;
274 c2 = 3;
275 } else if ( cs_mode == 1 ) {
276 //xy-plane
277 c1 = 1;
278 c2 = 2;
279 } else if ( cs_mode == 2 ) {
280 //yz-plane
281 c1 = 2;
282 c2 = 3;
283 } else {
284 OOFEM_ERROR("Unknow cs_mode");
285 }
286}
287
288void
290{
292 ParameterManager &ppm = this->giveDomain()->dofmanPPM;
293 PM_UPDATE_PARAMETER(cs_mode, ppm, ir, this->number, IPK_Truss2d_cs, priority) ;
294}
295
296void
298{
300 if ( cs_mode != 0 && cs_mode != 1 && cs_mode != 2 ) {
301 throw ComponentInputException(IPK_Truss2d_cs.getName(), ComponentInputException::ComponentType::ctElement, this->number, "Unsupported mode");
302 }
303}
304
305void
306Truss2d::giveDofManDofIDMask(int inode, IntArray &answer) const
307{
308 if ( cs_mode == 0 ) {
309 answer = { D_u, D_w };
310 } else if ( cs_mode == 1 ) {
311 answer = { D_u, D_v };
312 } else if ( cs_mode == 2 ) {
313 answer = { D_v, D_w };
314 }
315}
316
317
318void
320{
321 answer = this->giveStructuralCrossSection()->giveRealStress_1d(strain, gp, tStep);
322}
323
324
325void
326Truss2d::computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
327{
328 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_1d(rMode, gp, tStep);
329}
330
331
332void
334{
335 answer = this->giveStructuralCrossSection()->giveStiffnessMatrix_dPdF_1d(rMode, gp, tStep);
336}
337
338
339void
340Truss2d::giveEdgeDofMapping(IntArray &answer, int iEdge) const
341{
342 /*
343 * provides dof mapping of local edge dofs (only nonzero are taken into account)
344 * to global element dofs
345 */
346
347 if ( iEdge != 1 ) {
348 OOFEM_ERROR("wrong edge number");
349 }
350
351
352 answer.resize(4);
353 answer.at(1) = 1;
354 answer.at(2) = 2;
355 answer.at(3) = 3;
356 answer.at(4) = 4;
357}
358
359double
361{
362 if ( iEdge != 1 ) { // edge between nodes 1 2
363 OOFEM_ERROR("wrong egde number");
364 }
365
366 double weight = gp->giveWeight();
367 return 0.5 * this->computeLength() * weight;
368}
369
370int
372{
373 // returns transformation matrix from
374 // edge local coordinate system
375 // to element local c.s
376 // (same as global c.s in this case)
377 //
378 // i.e. f(element local) = T * f(edge local)
379 //
380 //double dx,dy, length ;
381 double sine, cosine;
382
383 answer.resize(2, 2);
384 answer.zero();
385
386 sine = sin(this->givePitch() );
387 cosine = cos(pitch);
388
389 answer.at(1, 1) = cosine;
390 answer.at(1, 2) = -sine;
391 answer.at(2, 1) = sine;
392 answer.at(2, 2) = cosine;
393
394 return 1;
395}
396
397
399{
400 // return interpolator
401 return & interp [ cs_mode ];
402}
403
404
405#ifdef __OOFEG
407{
408 int c1, c2;
409 resolveCoordIndices(c1, c2);
410
411 GraphicObj *go;
412 // if (!go) { // create new one
413 WCRec p[ 2 ]; /* point */
414 if ( !gc.testElementGraphicActivity(this) ) {
415 return;
416 }
417
418 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
419 EASValsSetColor(gc.getElementColor() );
420 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
421 if ( cs_mode == 0 ) {
422 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(c1);
423 p [ 0 ].y = 0.;
424 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(c2);
425 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(c1);
426 p [ 1 ].y = 0.;
427 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(c2);
428 } else if ( cs_mode == 1 ) {
429 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(c1);
430 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(c2);
431 p [ 0 ].z = 0.;
432 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(c1);
433 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(c2);
434 p [ 1 ].z = 0.;
435 } else if ( cs_mode == 2 ) {
436 p [ 0 ].x = 0.;
437 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(c1);
438 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(c2);
439 p [ 1 ].x = 0.;
440 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(c1);
441 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(c2);
442 }
443
444 go = CreateLine3D(p);
445 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
446 EGAttachObject(go, ( EObjectP ) this);
447 EMAddGraphicsToModel(ESIModel(), go);
448}
449
450
452{
453 int c1, c2;
454 resolveCoordIndices(c1, c2);
455
456 GraphicObj *go;
457 double defScale = gc.getDefScale();
458 // if (!go) { // create new one
459 WCRec p[ 2 ]; /* point */
460 if ( !gc.testElementGraphicActivity(this) ) {
461 return;
462 }
463
464 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
465 EASValsSetColor(gc.getDeformedElementColor() );
466 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
467
468 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
469 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
470 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
471
472 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
473 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
474 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
475
476 go = CreateLine3D(p);
477 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
478 EMAddGraphicsToModel(ESIModel(), go);
479}
480#endif
481} // end namespace oofem
#define REGISTER_Element(class)
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
double giveCoordinate(int i) const
Definition dofmanager.h:383
ParameterManager dofmanPPM
Definition domain.h:134
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
void postInitialize() override
Performs post initialization steps.
Definition element.C:797
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
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 zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void times(double f)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
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.
double computeVolumeAround(GaussPoint *gp) override
Definition truss2d.C:189
void computeGaussPoints() override
Definition truss2d.C:115
double pitch
Definition truss2d.h:64
void computeConstitutiveMatrix_dPdF_At(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss2d.C:333
FEInterpolation * giveInterpolation() const override
Definition truss2d.C:398
void giveEdgeDofMapping(IntArray &answer, int) const override
Definition truss2d.C:340
int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *gp) override
Definition truss2d.C:371
int giveLocalCoordinateSystem(FloatMatrix &answer) override
Definition truss2d.C:244
void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep) override
Definition truss2d.C:406
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition truss2d.C:319
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition truss2d.C:326
double computeLength() override
Definition truss2d.C:198
void postInitialize() override
Performs post initialization steps.
Definition truss2d.C:297
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition truss2d.C:127
double givePitch()
Definition truss2d.C:220
Truss2d(int n, Domain *d)
Definition truss2d.C:63
double computeEdgeVolumeAround(GaussPoint *gp, int) override
Definition truss2d.C:360
double length
Definition truss2d.h:63
static FEI2dLineLin interp[3]
Definition truss2d.h:67
void initializeFrom(InputRecord &ir, int priority) override
Definition truss2d.C:289
void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override
Definition truss2d.C:451
void resolveCoordIndices(int &c1, int &c2)
Definition truss2d.C:269
int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) override
Definition truss2d.C:168
void giveDofManDofIDMask(int inode, IntArray &answer) const override
Definition truss2d.C:306
void computeBHmatrixAt(GaussPoint *gp, FloatMatrix &) override
Definition truss2d.C:107
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &, int=1, int=ALL_STRAINS) override
Definition truss2d.C:74
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &) override
Definition truss2d.C:148
static ParamKey IPK_Truss2d_cs
Definition truss2d.h:68
#define OOFEM_ERROR(...)
Definition error.h:79
@ CS_Area
Area.
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
#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