OOFEM 3.0
Loading...
Searching...
No Matches
libeam3d.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 "crosssection.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "intarray.h"
44#include "floatarray.h"
45#include "mathfem.h"
46#include "classfactory.h"
47#include "parametermanager.h"
48#include "paramkey.h"
49
50#ifdef __OOFEG
51 #include "oofeggraphiccontext.h"
52#endif
53
54namespace oofem {
57
58LIBeam3d :: LIBeam3d(int n, Domain *aDomain) : StructuralElement(n, aDomain)
59 // Constructor.
60{
62 referenceNode = 0;
63 length = 0.;
64}
65
66
67void
68LIBeam3d :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
69// Returns the strain matrix of the receiver.
70// eeps = {\eps_x, \gamma_xz, \gamma_xy, \der{phi_x}{x}, \kappa_y, \kappa_z}^T
71{
72 double l, ksi, n1, n2, n1x, n2x;
73
74 l = this->computeLength();
75 ksi = gp->giveNaturalCoordinate(1);
76
77 answer.resize(6, 12);
78 answer.zero();
79
80 n1 = 0.5 * ( 1 - ksi );
81 n2 = 0.5 * ( 1. + ksi );
82 n1x = -1.0 / l;
83 n2x = 1.0 / l;
84
85 answer.at(1, 1) = -1. / l;
86 answer.at(1, 7) = 1. / l;
87
88 answer.at(2, 3) = n1x;
89 answer.at(2, 5) = n1;
90 answer.at(2, 9) = n2x;
91 answer.at(2, 11) = n2;
92
93 answer.at(3, 2) = n1x;
94 answer.at(3, 6) = -n1;
95 answer.at(3, 8) = n2x;
96 answer.at(3, 12) = -n2;
97
98 answer.at(4, 4) = -1. / l;
99 answer.at(4, 10) = 1. / l;
100
101 answer.at(5, 5) = n1x;
102 answer.at(5, 11) = n2x;
103
104 answer.at(6, 6) = n1x;
105 answer.at(6, 12) = n2x;
106}
107
108
109void
110LIBeam3d :: computeGaussPoints()
111// Sets up the array of Gauss Points of the receiver.
112{
113 if ( integrationRulesArray.size() == 0 ) {
114 integrationRulesArray.resize( 1 );
115 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 2);
116 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], 1, this);
117 }
118}
119
120
121void
122LIBeam3d :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
123{
124 answer = this->giveStructuralCrossSection()->give3dBeamStiffMtrx(rMode, gp, tStep);
125}
126
127
128void
129LIBeam3d :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
130{
131 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Beam3d(strain, gp, tStep);
132}
133
134
135void
136LIBeam3d :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
137// Returns the lumped mass matrix of the receiver. This expression is
138// valid in both local and global axes.
139{
140 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
141 double density = this->giveStructuralCrossSection()->give('d', gp);
142 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
143 answer.resize(12, 12);
144 answer.zero();
145 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass;
146 answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass;
147}
148
149
150void
151LIBeam3d :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
152// Returns the displacement interpolation matrix {N} of the receiver, eva-
153// luated at gp.
154{
155 double ksi, n1, n2;
156
157 ksi = iLocCoord.at(1);
158 n1 = ( 1. - ksi ) * 0.5;
159 n2 = ( 1. + ksi ) * 0.5;
160
161 answer.resize(6, 12);
162 answer.zero();
163
164 answer.at(1, 1) = n1;
165 answer.at(1, 7) = n2;
166 answer.at(2, 2) = n1;
167 answer.at(2, 8) = n2;
168 answer.at(3, 3) = n1;
169 answer.at(3, 9) = n2;
170
171 answer.at(4, 4) = n1;
172 answer.at(4, 10) = n2;
173 answer.at(5, 5) = n1;
174 answer.at(5, 11) = n2;
175 answer.at(6, 6) = n1;
176 answer.at(6, 12) = n2;
177}
178
179void
180LIBeam3d :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
181// Returns the stiffness matrix of the receiver, expressed in the global
182// axes.
183{
184 StructuralElement :: computeStiffnessMatrix(answer, rMode, tStep);
185}
186
187
188bool
189LIBeam3d :: computeGtoLRotationMatrix(FloatMatrix &answer)
190{
191 FloatMatrix lcs;
192
193 answer.resize(12, 12);
194 answer.zero();
195
196 this->giveLocalCoordinateSystem(lcs);
197 for ( int i = 1; i <= 3; i++ ) {
198 for ( int j = 1; j <= 3; j++ ) {
199 answer.at(i, j) = lcs.at(i, j);
200 answer.at(i + 3, j + 3) = lcs.at(i, j);
201 answer.at(i + 6, j + 6) = lcs.at(i, j);
202 answer.at(i + 9, j + 9) = lcs.at(i, j);
203 }
204 }
205 return true;
206}
207
208
209double
210LIBeam3d :: computeVolumeAround(GaussPoint *gp)
211// Returns the length of the receiver. This method is valid only if 1
212// Gauss point is used.
213{
214 double weight = gp->giveWeight();
215 return weight * 0.5 * this->computeLength();
216}
217
218
219int
220LIBeam3d :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
221{
224 if ( type == IST_BeamForceMomentTensor ) {
225 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
226 return 1;
227 } else if ( type == IST_BeamStrainCurvatureTensor ) {
228 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
229 return 1;
230 } else {
231 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
232 }
233}
234
235
236void
237LIBeam3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
238{
239 answer = {D_u, D_v, D_w, R_u, R_v, R_w};
240}
241
242
243int
244LIBeam3d :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
245{
246 double ksi, n1, n2;
247
248 ksi = lcoords.at(1);
249 n1 = ( 1. - ksi ) * 0.5;
250 n2 = ( 1. + ksi ) * 0.5;
251
252 answer.resize(3);
253 answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
254 answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 *this->giveNode(2)->giveCoordinate(2);
255 answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
256
257 return 1;
258}
259
260
261int
262LIBeam3d :: testElementExtension(ElementExtension ext)
263{
264 return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 );
265}
266
267
268double
269LIBeam3d :: computeLength()
270// Returns the length of the receiver.
271{
272 double dx, dy, dz;
273 Node *nodeA, *nodeB;
274
275 if ( length == 0. ) {
276 nodeA = this->giveNode(1);
277 nodeB = this->giveNode(2);
278 dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
279 dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
280 dz = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
281 length = sqrt(dx * dx + dy * dy + dz * dz);
282 }
283
284 return length;
285}
286
287
288void
289LIBeam3d :: initializeFrom(InputRecord &ir, int priority)
290{
291 StructuralElement :: initializeFrom(ir, priority);
292 ParameterManager &ppm = domain->elementPPM;
294}
295
296void
297LIBeam3d :: postInitialize()
298{
299 ParameterManager &ppm = domain->elementPPM;
300
301 StructuralElement :: postInitialize();
303 if ( referenceNode == 0 ) {
304 OOFEM_ERROR("wrong reference node specified.");
305 }
306}
307
308void
309LIBeam3d :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
310{
311 /*
312 * provides dof mapping of local edge dofs (only nonzero are taken into account)
313 * to global element dofs
314 */
315 if ( iEdge != 1 ) {
316 OOFEM_ERROR("wrong edge number");
317 }
318
319
320 answer.resize(12);
321 for ( int i = 1; i <= 12; i++ ) {
322 answer.at(i) = i;
323 }
324}
325
326double
327LIBeam3d :: computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
328{
329 if ( iEdge != 1 ) { // edge between nodes 1 2
330 OOFEM_ERROR("wrong egde number");
331 }
332
333 double weight = gp->giveWeight();
334 return 0.5 * this->computeLength() * weight;
335}
336
337int
338LIBeam3d :: computeLoadGToLRotationMtrx(FloatMatrix &answer)
339{
340 /*
341 * Returns transformation matrix from global coordinate system to local
342 * element coordinate system for element load vector components.
343 * If no transformation is necessary, answer is empty matrix (default);
344 */
345
346 // f(elemLocal) = T * f (global)
347
348 FloatMatrix lcs;
349
350 answer.resize(6, 6);
351 answer.zero();
352
353 this->giveLocalCoordinateSystem(lcs);
354 for ( int i = 1; i <= 3; i++ ) {
355 for ( int j = 1; j <= 3; j++ ) {
356 answer.at(i, j) = lcs.at(i, j);
357 answer.at(3 + i, 3 + j) = lcs.at(i, j);
358 }
359 }
360 return 1;
361}
362
363
364int
365LIBeam3d :: computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
366{
367 // returns transformation matrix from
368 // edge local coordinate system
369 // to element local c.s
370 // (same as global c.s in this case)
371 //
372 // i.e. f(element local) = T * f(edge local)
373 //
374 answer.clear();
375 return 0;
376}
377
378void
379LIBeam3d :: computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
380{
381 FloatArray lc(1);
382 StructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
383 answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
384}
385
386Node* LIBeam3d :: giveReferenceNode(int refNode)
387//
388// returns pointer to the reference node in the general case that
389// the global number of reference node is not equal to its consecutive number
390// in the domain
391//
392{
393 int nnodes = this->giveDomain()->giveNumberOfDofManagers();
394
395 for ( int i = 1; i <= nnodes; i++ ) {
396 Node* referenceNode = this->giveDomain()->giveNode(i);
397 int globalNumber = referenceNode->giveGlobalNumber();
398 if ( globalNumber == refNode ) {
399 return referenceNode;
400 }
401 }
402
403 OOFEM_ERROR("Could not find the reference node. Check numbering.");
404}
405
406int
407LIBeam3d :: giveLocalCoordinateSystem(FloatMatrix &answer)
408//
409// returns a unit vectors of local coordinate system at element
410// stored rowwise (mainly used by some materials with ortho and anisotrophy)
411//
412{
413 FloatArray lx(3), ly(3), lz(3), help(3);
414 double length = this->computeLength();
415 Node *nodeA, *nodeB, *refNode;
416
417 answer.resize(3, 3);
418 answer.zero();
419 nodeA = this->giveNode(1);
420 nodeB = this->giveNode(2);
421 refNode = this->giveReferenceNode(referenceNode);
422
423 for ( int i = 1; i <= 3; i++ ) {
424 lx.at(i) = ( nodeB->giveCoordinate(i) - nodeA->giveCoordinate(i) ) / length;
425 help.at(i) = ( refNode->giveCoordinate(i) - nodeA->giveCoordinate(i) );
426 }
427
428 lz.beVectorProductOf(lx, help);
429 lz.normalize();
430 ly.beVectorProductOf(lz, lx);
431 ly.normalize();
432
433 for ( int i = 1; i <= 3; i++ ) {
434 answer.at(1, i) = lx.at(i);
435 answer.at(2, i) = ly.at(i);
436 answer.at(3, i) = lz.at(i);
437 }
438
439 return 1;
440}
441
442void
443LIBeam3d :: FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain,
444 GaussPoint *slaveGp, TimeStep *tStep)
445{
446 double layerYCoord, layerZCoord;
447
448 layerZCoord = slaveGp->giveNaturalCoordinate(2);
449 layerYCoord = slaveGp->giveNaturalCoordinate(1);
450
451 answer.resize(3); // {Exx,GMzx,GMxy}
452
453 answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(5) * layerZCoord - masterGpStrain.at(6) * layerYCoord;
454 answer.at(2) = masterGpStrain.at(2);
455 answer.at(3) = masterGpStrain.at(3);
456}
457
458void
459LIBeam3d :: NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node,
460 InternalStateType type, TimeStep *tStep)
461{
462 GaussPoint* gp = integrationRulesArray[0]->getIntegrationPoint(0);
463 this->giveIPValue(answer, gp, type, tStep);
464}
465
466Interface *
467LIBeam3d :: giveInterface(InterfaceType interface)
468{
469 if ( interface == FiberedCrossSectionInterfaceType ) {
470 return static_cast< FiberedCrossSectionInterface * >(this);
471 } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
472 return static_cast< NodalAveragingRecoveryModelInterface * >(this);
473 }
474
475 return NULL;
476}
477
478
479#ifdef __OOFEG
480void
481LIBeam3d :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
482{
483 GraphicObj *go;
484
485 if ( !gc.testElementGraphicActivity(this) ) {
486 return;
487 }
488
489 // if (!go) { // create new one
490 WCRec p [ 2 ]; /* poin */
491 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
492 EASValsSetColor( gc.getElementColor() );
493 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
494 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
495 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
496 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
497 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
498 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
499 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
500 go = CreateLine3D(p);
501 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
502 EGAttachObject(go, ( EObjectP ) this);
503 EMAddGraphicsToModel(ESIModel(), go);
504}
505
506
507void
508LIBeam3d :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
509{
510 GraphicObj *go;
511
512 if ( !gc.testElementGraphicActivity(this) ) {
513 return;
514 }
515
516 double defScale = gc.getDefScale();
517 // if (!go) { // create new one
518 WCRec p [ 2 ]; /* poin */
519 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
520 EASValsSetColor( gc.getDeformedElementColor() );
521 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
522 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
523 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
524 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
525
526 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
527 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
528 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
529 go = CreateLine3D(p);
530 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
531 EMAddGraphicsToModel(ESIModel(), go);
532}
533#endif
534} // 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
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
Domain * giveDomain() const
Definition femcmpnn.h:97
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
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 times(double s)
Definition floatarray.C:834
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 zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition gausspoint.h:136
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
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
static ParamKey IPK_LIBeam3d_refnode
Definition libeam3d.h:61
double computeLength() override
Definition libeam3d.C:269
int giveLocalCoordinateSystem(FloatMatrix &answer) override
Definition libeam3d.C:407
double length
Definition libeam3d.h:58
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition libeam3d.C:220
Node * giveReferenceNode(int refNode)
Definition libeam3d.C:386
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
StructuralElement(int n, Domain *d)
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver's strain vector.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
#define OOFEM_ERROR(...)
Definition error.h:79
@ CS_Area
Area.
@ NodalAveragingRecoveryModelInterfaceType
@ FiberedCrossSectionInterfaceType
@ Element_EdgeLoadSupport
Element extension for edge loads.
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