OOFEM 3.0
Loading...
Searching...
No Matches
latticebeam3d.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 program is free software; you can redistribute it and/or modify
21 * it under the terms of the GNU General Public License as published by
22 * the Free Software Foundation; either version 2 of the License, or
23 * (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
28 * GNU General Public License for more details.
29 *
30 * You should have received a copy of the GNU General Public License
31 * along with this program; if not, write to the Free Software
32 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
33 */
34
35#include "domain.h"
39#include "node.h"
40#include "material.h"
41#include "gausspoint.h"
43#include "floatmatrix.h"
44#include "intarray.h"
45#include "floatarray.h"
46#include "mathfem.h"
48#include "contextioerr.h"
49#include "datastream.h"
50#include "classfactory.h"
52#include "parametermanager.h"
53#include "paramkey.h"
54
55#ifdef __OOFEG
56 #include "oofeggraphiccontext.h"
58#endif
59
60namespace oofem {
63
64LatticeBeam3d :: LatticeBeam3d(int n, Domain *aDomain) : LatticeStructuralElement(n, aDomain)
65{
66 geometryFlag = 0;
68 myPi = 3.14159;
69}
70
71LatticeBeam3d :: ~LatticeBeam3d()
72{}
73
74
75void
76LatticeBeam3d :: giveGPCoordinates(FloatArray &coords)
77{
78 if ( geometryFlag == 0 ) {
80 }
81 coords.resize(3);
82 coords = this->globalCentroid;
83
84 return;
85}
86
87double
88LatticeBeam3d :: giveLength()
89{
90 if ( geometryFlag == 0 ) {
92 }
93
94 return this->length;
95}
96
97void
98LatticeBeam3d :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
99{
100 answer = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->give3dStiffnessMatrix(rMode, gp, tStep);
101}
102
103
104void
105LatticeBeam3d :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
106 TimeStep *tStep)
107// Computes numerically the stiffness matrix of the receiver.
108{
109 FloatMatrix d, bi, bj, bjt, dbj, dij;
110
111 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
112
113 answer.resize(12, 12);
114 answer.zero();
115
116 //Stiffness matrix of Bernoulli beam from McGuire
117 double a = pow(this->diameter / 2., 2.) * myPi;
118 double iz = pow(this->diameter / 2., 4.) * myPi / 4.;
119 double iy = iz;
120 double j = iz + iy;
121 double l = this->giveLength();
122
123 double nu = ( static_cast< LatticeLinearElastic * >( this->giveMaterial() ) )->give('n', gp);
124 double e = ( static_cast< LatticeLinearElastic * >( this->giveMaterial() ) )->give('E', gp);
125
126 answer.at(1, 1) = a / l;
127 answer.at(1, 7) = -a / l;
128
129 answer.at(2, 2) = 12. * iz / pow(l, 3.);
130 answer.at(2, 6) = 6. * iz / pow(l, 2.);
131 answer.at(2, 8) = -12. * iz / pow(l, 3.);
132 answer.at(2, 12) = 6. * iz / pow(l, 2.);
133
134 answer.at(3, 3) = 12 * iy / pow(l, 3.);
135 answer.at(3, 5) = -6. * iy / pow(l, 2.);
136 answer.at(3, 9) = -12. * iy / pow(l, 3.);
137 answer.at(3, 11) = -6. * iy / pow(l, 2.);
138
139 answer.at(4, 4) = j / ( 2. * ( 1. + nu ) * l );
140 answer.at(4, 10) = -j / ( 2. * ( 1. + nu ) * l );
141
142 answer.at(5, 3) = answer.at(3, 5);
143 answer.at(5, 5) = 4 * iy / l;
144 answer.at(5, 9) = 6. * iy / pow(l, 2.);
145 answer.at(5, 11) = 2. * iy / l;
146
147 answer.at(6, 2) = answer.at(2, 6);
148 answer.at(6, 6) = 4. * iz / l;
149 answer.at(6, 8) = -6. * iz / pow(l, 2.);
150 answer.at(6, 12) = 2. * iz / l;
151
152 answer.at(7, 1) = answer.at(1, 7);
153 answer.at(7, 7) = a / l;
154
155 answer.at(8, 2) = answer.at(2, 8);
156 answer.at(8, 6) = answer.at(6, 8);
157 answer.at(8, 8) = 12. * iz / pow(l, 3.);
158 answer.at(8, 12) = -6. * iz / pow(l, 2.);
159
160 answer.at(9, 3) = answer.at(3, 9);
161 answer.at(9, 5) = answer.at(5, 9);
162 answer.at(9, 9) = 12 * iy / pow(l, 3.);
163 answer.at(9, 11) = 6 * iy / pow(l, 2.);
164
165 answer.at(10, 4) = answer.at(4, 10.);
166 answer.at(10, 10) = j / ( 2 * ( 1 + nu ) * l );
167
168 answer.at(11, 3) = answer.at(3, 11);
169 answer.at(11, 5) = answer.at(5, 11);
170 answer.at(11, 9) = answer.at(9, 11);
171 answer.at(11, 11) = 4 * iy / l;
172
173 answer.at(12, 2) = answer.at(2, 12);
174 answer.at(12, 6) = answer.at(6, 12);
175 answer.at(12, 8) = answer.at(8, 12);
176 answer.at(12, 12) = 4 * iz / l;
177
178 answer.times(e);
179
180 return;
181}
182
183void LatticeBeam3d :: computeGaussPoints()
184// Sets up the array of Gauss Points of the receiver.
185{
186 integrationRulesArray.resize(1);
187 integrationRulesArray [ 0 ].reset(new GaussIntegrationRule(1, this, 1, 3) );
188 integrationRulesArray [ 0 ]->SetUpPointsOnLine(1, _3dLattice);
189}
190
191
192
193double LatticeBeam3d :: giveArea() {
194 if ( geometryFlag == 0 ) {
196 }
197
198 return this->area;
199}
200
201
202bool
203LatticeBeam3d :: computeGtoLRotationMatrix(FloatMatrix &answer)
204{
205 FloatMatrix lcs;
206 int i, j;
207
208 answer.resize(12, 12);
209 answer.zero();
210
211 this->giveLocalCoordinateSystem(lcs);
212 for ( i = 1; i <= 3; i++ ) {
213 for ( j = 1; j <= 3; j++ ) {
214 answer.at(i, j) = lcs.at(i, j);
215 answer.at(i + 3, j + 3) = lcs.at(i, j);
216 answer.at(i + 6, j + 6) = lcs.at(i, j);
217 answer.at(i + 9, j + 9) = lcs.at(i, j);
218 }
219 }
220
221 return 1;
222}
223
224int
225LatticeBeam3d :: giveLocalCoordinateSystem(FloatMatrix &answer)
226{
227 if ( geometryFlag == 0 ) {
229 }
230
231 answer = this->localCoordinateSystem;
232
233 return 1;
234}
235
236
237double
238LatticeBeam3d :: computeVolumeAround(GaussPoint *aGaussPoint)
239{
240 if ( geometryFlag == 0 ) {
242 }
243
244 return this->area * this->length;
245}
246
247
248void
249LatticeBeam3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
250{
251 answer = {
252 D_u, D_v, D_w, R_u, R_v, R_w
253 };
254}
255
256void
257LatticeBeam3d :: initializeFrom(InputRecord &ir, int priority)
258{
259 ParameterManager &ppm = this->giveDomain()->elementPPM;
260 LatticeStructuralElement :: initializeFrom(ir, priority);
262}
263
264void
265LatticeBeam3d :: postInitialize()
266{
267 ParameterManager &ppm = this->giveDomain()->elementPPM;
268 this->LatticeStructuralElement :: postInitialize();
270}
271
272int
273LatticeBeam3d :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
274{
275 if ( geometryFlag == 0 ) {
277 }
278
279 answer.resize(3);
280 answer = this->globalCentroid;
281
282 return 1;
283}
284
285
286void
287LatticeBeam3d :: computeGeometryProperties()
288{
289 //coordinates of the two nodes
290 Node *nodeA, *nodeB;
291 FloatArray coordsA(3), coordsB(3);
292
293 nodeA = this->giveNode(1);
294 nodeB = this->giveNode(2);
295
296 for ( int i = 0; i < 3; i++ ) {
297 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1);
298 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1);
299 }
300
301 //Construct an initial temporary local coordinate system
302 FloatArray s(3), t(3);
303
304 //Calculate normal vector
305 this->normal.resize(3);
306 for ( int i = 0; i < 3; i++ ) {
307 this->normal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
308 }
309
310 this->length = sqrt(pow(normal.at(1), 2.) + pow(normal.at(2), 2.) + pow(normal.at(3), 2.) );
311
312 // Compute midpoint
313 this->midPoint.resize(3);
314 for ( int i = 0; i < 3; i++ ) {
315 this->midPoint.at(i + 1) = 0.5 * ( coordsB.at(i + 1) + coordsA.at(i + 1) );
316 }
317
318 for ( int i = 0; i < 3; i++ ) {
319 this->normal.at(i + 1) /= length;
320 }
321
322 this->globalCentroid = this->midPoint;
323
325
326 this->geometryFlag = 1;
327
328 return;
329}
330
331void
332LatticeBeam3d :: computeCrossSectionProperties() {
333 //Construct two perpendicular axis so that n is normal to the plane which they create
334 //Check, if one of the components of the normal-direction is zero
335 FloatArray s(3), t(3);
336 if ( this->normal.at(1) == 0 ) {
337 s.at(1) = 0.;
338 s.at(2) = this->normal.at(3);
339 s.at(3) = -this->normal.at(2);
340 } else if ( this->normal.at(2) == 0 ) {
341 s.at(1) = this->normal.at(3);
342 s.at(2) = 0.;
343 s.at(3) = -this->normal.at(1);
344 } else {
345 s.at(1) = this->normal.at(2);
346 s.at(2) = -this->normal.at(1);
347 s.at(3) = 0.;
348 }
349
350 s.normalize();
351
352 t.beVectorProductOf(this->normal, s);
353 t.normalize();
354
355 this->localCoordinateSystem.resize(3, 3);
356 this->localCoordinateSystem.zero();
357
358 //Set up rotation matrix
359 for ( int i = 1; i <= 3; i++ ) {
360 this->localCoordinateSystem.at(1, i) = this->normal.at(i);
361 this->localCoordinateSystem.at(2, i) = s.at(i);
362 this->localCoordinateSystem.at(3, i) = t.at(i);
363 }
364}
365
366void
367LatticeBeam3d :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
368{
369 answer = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->giveLatticeStress3d(strain, gp, tStep);
370}
371
372
373void
374LatticeBeam3d :: giveInternalForcesVector(FloatArray &answer,
375 TimeStep *tStep, int useUpdatedGpRecord)
376{
377 FloatArray u;
378 FloatMatrix stiffness;
379
380 FloatArray strain(6), stress(6);
381
382 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
383
384 // This function can be quite costly to do inside the loops when one has many slave dofs.
385 this->computeVectorOf(VM_Total, tStep, u);
386 // subtract initial displacements, if defined
387 if ( initialDisplacements ) {
389 }
390
391 // printf("Check what u is.\n");
392 // u.printYourself();
393
394 this->computeStiffnessMatrix(stiffness, ElasticStiffness, tStep);
395
396 // zero answer will resize accordingly when adding first contribution
397 answer.clear();
398 answer.beProductOf(stiffness, u);
399 // printf("answer.at(1) = %e, answer.at(7) = %e\n", answer.at(1), answer.at(7));
400
401 double area = pow(this->diameter / 2., 2.) * myPi;
402 //Apply yield limit to axial component
403 strain.zero();
404 strain.at(1) = ( u.at(7) - u.at(1) ) / this->giveLength();
405 this->computeStressVector(stress, strain, gp, tStep);
406 answer.at(1) = -stress.at(1) * area;
407 answer.at(7) = -answer.at(1);
408 // printf("answerNew.at(1) = %e, answerNew.at(7) = %e\n", answer.at(1), answer.at(7));
409
410 // if inactive update state, but no contribution to global system
411 if ( !this->isActivated(tStep) ) {
412 answer.zero();
413 return;
414 }
415 return;
416}
417
418#ifdef __OOFEG
419
420void
421LatticeBeam3d :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
422{
423 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
424
425 if ( mode == OGC_rawGeometry ) {
426 this->drawRawGeometry(gc, tStep);
427 } else if ( mode == OGC_deformedGeometry ) {
428 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
429 } else if ( mode == OGC_eigenVectorGeometry ) {
430 this->drawDeformedGeometry(gc, tStep, EigenVector);
431 } else if ( mode == OGC_scalarPlot ) {
432 this->drawScalar(gc, tStep);
433 } else if ( mode == OGC_elemSpecial ) {
434 this->drawSpecial(gc, tStep);
435 } else {
436 OOFEM_ERROR("unsupported mode");
437 }
438}
439
440
441
442void LatticeBeam3d :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
443{
444 GraphicObj *go;
445
446 WCRec p [ 2 ]; /* points */
447 if ( !gc.testElementGraphicActivity(this) ) {
448 return;
449 }
450
451 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
452 EASValsSetColor(gc.getElementColor() );
453 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
454
455 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
456 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
457 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
458 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
459 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
460 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
461
462 go = CreateLine3D(p);
463 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
464 EGAttachObject(go, ( EObjectP ) this);
465 EMAddGraphicsToModel(ESIModel(), go);
466}
467
468
469void LatticeBeam3d :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
470{
471 GraphicObj *go;
472
473 if ( !gc.testElementGraphicActivity(this) ) {
474 return;
475 }
476
477 double defScale = gc.getDefScale();
478
479 WCRec p [ 2 ]; /* points */
480
481 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
482 EASValsSetColor(gc.getDeformedElementColor() );
483 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
484
485 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
486 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
487 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
488
489 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
490 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
491 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
492
493 go = CreateLine3D(p);
494 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
495 EMAddGraphicsToModel(ESIModel(), go);
496}
497
498#endif
499} // end namespace oofem
#define REGISTER_Element(class)
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
virtual Material * giveMaterial()
Definition element.C:523
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition element.h:1077
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
Definition element.h:1078
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 beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void subtract(const FloatArray &src)
Definition floatarray.C:320
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
virtual double giveLength() override
FloatArray globalCentroid
FloatMatrix localCoordinateSystem
virtual int giveLocalCoordinateSystem(FloatMatrix &answer) override
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
virtual void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
virtual void computeGeometryProperties()
static ParamKey IPK_LatticeBeam3d_diameter
virtual void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
virtual void computeCrossSectionProperties()
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
#define OOFEM_ERROR(...)
Definition error.h:79
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_ELEMENT_ERROR_IFNOTSET(_pm, _componentnum, _paramkey)
#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