OOFEM 3.0
Loading...
Searching...
No Matches
libeam2d.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 "fei2dlinelin.h"
38#include "node.h"
39#include "material.h"
40#include "crosssection.h"
41#include "gausspoint.h"
43#include "floatmatrix.h"
44#include "intarray.h"
45#include "floatarray.h"
46#include "mathfem.h"
47#include "classfactory.h"
48#include "parametermanager.h"
49#include "paramkey.h"
50
51namespace oofem {
55
56
57// Set up interpolation coordinates
58FEI2dLineLin LIBeam2d :: interpolationXZ(1, 3);
59FEI2dLineLin LIBeam2d :: interpolationXY(1, 2);
60
61LIBeam2d :: LIBeam2d(int n, Domain *aDomain) : StructuralElement(n, aDomain), LayeredCrossSectionInterface()
62{
64 length = 0.;
65 pitch = 10.; // a dummy value
66}
67
68
69FEInterpolation *LIBeam2d :: giveInterpolation() const { if (xy) {return & interpolationXY;} else {return & interpolationXZ;} }
70
71
73LIBeam2d :: giveInterface(InterfaceType interface)
74{
75 if ( interface == LayeredCrossSectionInterfaceType ) {
76 return static_cast< LayeredCrossSectionInterface * >(this);
77 }
78
79 return NULL;
80}
81
82
83void
84LIBeam2d :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
85// Returns the strain matrix of the receiver.
86{
87 double l, ksi, n1x, n4x, n3xx, n6xx, n2xxx, n3xxx, n5xxx, n6xxx;
88
89 l = this->computeLength();
90 ksi = gp->giveNaturalCoordinate(1);
91 n1x = -1.0 / l;
92 n4x = 1.0 / l;
93 n3xx = -1.0 / l;
94 n6xx = 1.0 / l;
95 n2xxx = -1.0 / l;
96 n3xxx = 0.5 * ( 1 - ksi );
97 n5xxx = 1.0 / l;
98 n6xxx = 0.5 * ( 1. + ksi );
99
100 answer.resize(3, 6);
101 answer.zero();
102
103 answer.at(1, 1) = n1x;
104 answer.at(1, 4) = n4x;
105 answer.at(2, 3) = n3xx;
106 answer.at(2, 6) = n6xx;
107 answer.at(3, 2) = n2xxx;
108 answer.at(3, 3) = n3xxx;
109 answer.at(3, 5) = n5xxx;
110 answer.at(3, 6) = n6xxx;
111}
112
113
114void
115LIBeam2d :: computeGaussPoints()
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);
121 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], 1, this);
122 }
123}
124
125
126void
127LIBeam2d :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
128{
129 answer = this->giveStructuralCrossSection()->give2dBeamStiffMtrx(rMode, gp, tStep);
130}
131
132
133void
134LIBeam2d :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
135{
136 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Beam2d(strain, gp, tStep);
137}
138
139
140int
141LIBeam2d :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
142{
143 double ksi, n1, n2;
144
145 ksi = lcoords.at(1);
146 n1 = ( 1. - ksi ) * 0.5;
147 n2 = ( 1. + ksi ) * 0.5;
148
149 answer.resize(3);
150 answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
151 if (xy) {
152 answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 *this->giveNode(2)->giveCoordinate(2);
153 } else {
154 answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
155 }
156
157 return 1;
158}
159
160
161void
162LIBeam2d :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
163// Returns the lumped mass matrix of the receiver. This expression is
164// valid in both local and global axes.
165{
166 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
167 double density = this->giveStructuralCrossSection()->give('d', gp);
168 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
169 answer.resize(6, 6);
170 answer.zero();
171 answer.at(1, 1) = halfMass;
172 answer.at(2, 2) = halfMass;
173 answer.at(4, 4) = halfMass;
174 answer.at(5, 5) = halfMass;
175}
176
177
178void
179LIBeam2d :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
180// Returns the displacement interpolation matrix {N} of the receiver, eva-
181// luated at gp.
182{
183 double ksi, n1, n2;
184
185 ksi = iLocCoord.at(1);
186 n1 = ( 1. - ksi ) * 0.5;
187 n2 = ( 1. + ksi ) * 0.5;
188
189 answer.resize(3, 6);
190 answer.zero();
191
192 answer.at(1, 1) = n1;
193 answer.at(1, 4) = n2;
194 answer.at(2, 2) = n1;
195 answer.at(2, 5) = n2;
196 answer.at(3, 3) = n1;
197 answer.at(3, 6) = n2;
198}
199
200
201void
202LIBeam2d :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
203// Returns the stiffness matrix of the receiver, expressed in the global
204// axes.
205{
206 StructuralElement :: computeStiffnessMatrix(answer, rMode, tStep);
207}
208
209
210bool
211LIBeam2d :: computeGtoLRotationMatrix(FloatMatrix &answer)
212{
213 double sine, cosine;
214
215 sine = sin( this->givePitch() );
216 cosine = cos(pitch);
217
218 answer.resize(6, 6);
219 answer.zero();
220
221 answer.at(1, 1) = cosine;
222 answer.at(1, 2) = sine;
223 answer.at(2, 1) = -sine;
224 answer.at(2, 2) = cosine;
225 answer.at(3, 3) = 1.;
226 answer.at(4, 4) = cosine;
227 answer.at(4, 5) = sine;
228 answer.at(5, 4) = -sine;
229 answer.at(5, 5) = cosine;
230 answer.at(6, 6) = 1.;
231
232 return true;
233}
234
235
236double
237LIBeam2d :: computeVolumeAround(GaussPoint *gp)
238// Returns the length of the receiver. This method is valid only if 1
239// Gauss point is used.
240{
241 double weight = gp->giveWeight();
242 return weight * 0.5 * this->computeLength();
243}
244
245
246void
247LIBeam2d :: computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
248//
249// returns full 3d strain vector of given layer (whose z-coordinate from center-line is
250// stored in slaveGp) for given tStep
251//
252{
253 double layerZeta, layerZCoord, top, bottom;
254
255 top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
256 bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
257 layerZeta = slaveGp->giveNaturalCoordinate(3);
258 layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
259
260 answer.resize(2); // {Exx,GMzx}
261
262 answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(2) * layerZCoord;
263 answer.at(2) = masterGpStrain.at(3);
264}
265
266
267int
268LIBeam2d :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
269{
270 if ( type == IST_BeamForceMomentTensor ) {
271 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
272 return 1;
273 } else if ( type == IST_BeamStrainCurvatureTensor ) {
274 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
275 return 1;
276 } else {
277 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
278 }
279}
280
281
282void
283LIBeam2d :: giveDofManDofIDMask(int inode, IntArray &answer) const
284{
285 if (xy) {
286 answer = {D_u, D_v, R_w};
287 } else {
288 answer = {D_u, D_w, R_v};
289 }
290}
291
292
293double
294LIBeam2d :: computeLength()
295// Returns the length of the receiver.
296{
297 double dx, dy;
298 Node *nodeA, *nodeB;
299
300 if ( length == 0. ) {
301 nodeA = this->giveNode(1);
302 nodeB = this->giveNode(2);
303 dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
304 if (xy) {
305 dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
306 } else {
307 dy = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
308 }
309 length = sqrt(dx * dx + dy * dy);
310 }
311
312 return length;
313}
314
315
316double
317LIBeam2d :: givePitch()
318// Returns the pitch of the receiver.
319{
320 double xA, xB, yA, yB;
321 Node *nodeA, *nodeB;
322
323 if ( pitch == 10. ) { // 10. : dummy initialization value
324 nodeA = this->giveNode(1);
325 nodeB = this->giveNode(2);
326 xA = nodeA->giveCoordinate(1);
327 xB = nodeB->giveCoordinate(1);
328 if (xy) {
329 yA = nodeA->giveCoordinate(2);
330 yB = nodeB->giveCoordinate(2);
331 } else {
332 yA = nodeA->giveCoordinate(3);
333 yB = nodeB->giveCoordinate(3);
334 }
335 pitch = atan2(yB - yA, xB - xA);
336 }
337
338 return pitch;
339}
340
341
342void
343LIBeam2d :: initializeFrom(InputRecord &ir, int priority)
344{
345 StructuralElement :: initializeFrom(ir, priority);
346 ParameterManager &ppm = domain->elementPPM;
347 PM_CHECK_FLAG_AND_REPORT(ppm, ir, this->number, IPK_LIBeam2d_XY, priority, xy) ;
348 PM_CHECK_FLAG_AND_REPORT(ppm, ir, this->number, IPK_LIBeam2d_XZ, priority, xz) ;
349}
350
351void
352LIBeam2d :: initializeFinish()
353{
354 StructuralElement :: initializeFinish();
355 // ParameterManager &ppm = domain->elementPPM;
356 xy ? (xz = false) : (xz = true);
357}
358
359
360void
361LIBeam2d :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
362{
363 /*
364 * provides dof mapping of local edge dofs (only nonzero are taken into account)
365 * to global element dofs
366 */
367
368 if ( iEdge != 1 ) {
369 OOFEM_ERROR("wrong edge number");
370 }
371
372 answer.resize(6);
373 answer.at(1) = 1;
374 answer.at(2) = 2;
375 answer.at(3) = 3;
376 answer.at(4) = 4;
377 answer.at(5) = 5;
378 answer.at(6) = 6;
379}
380
381
382double
383LIBeam2d :: computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
384{
385 if ( iEdge != 1 ) { // edge between nodes 1 2
386 OOFEM_ERROR("wrong egde number");
387 }
388
389 double weight = gp->giveWeight();
390 return 0.5 * this->computeLength() * weight;
391}
392
393
394void
395LIBeam2d :: computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
396{
397 FloatArray lc(1);
398 StructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
399 answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
400}
401
402
403int
404LIBeam2d :: computeLoadGToLRotationMtrx(FloatMatrix &answer)
405{
406 /*
407 * Returns transformation matrix from global coordinate system to local
408 * element coordinate system for element load vector components.
409 * If no transformation is necessary, answer is empty matrix (default);
410 */
411
412 // f(elemLocal) = T * f (global)
413
414 double sine, cosine;
415
416 answer.resize(3, 3);
417 answer.zero();
418
419 sine = sin( this->givePitch() );
420 cosine = cos(pitch);
421
422 answer.at(1, 1) = cosine;
423 answer.at(1, 2) = -sine;
424 answer.at(2, 1) = sine;
425 answer.at(2, 2) = cosine;
426 answer.at(3, 3) = 1.0;
427
428 return 1;
429}
430
431
432int
433LIBeam2d :: computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
434{
435 // returns transformation matrix from
436 // edge local coordinate system
437 // to element local c.s
438 // (same as global c.s in this case)
439 //
440 // i.e. f(element local) = T * f(edge local)
441 //
442 answer.clear();
443 return 0;
444}
445} // 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 * 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
double length
Definition libeam2d.h:64
double computeLength() override
Definition libeam2d.C:294
static ParamKey IPK_LIBeam2d_XZ
Definition libeam2d.h:60
static FEI2dLineLin interpolationXZ
Interpolation.
Definition libeam2d.h:57
static FEI2dLineLin interpolationXY
Definition libeam2d.h:58
double givePitch()
Definition libeam2d.C:317
static ParamKey IPK_LIBeam2d_XY
Definition libeam2d.h:61
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_BottomZCoord
Bottom z coordinate.
@ CS_Area
Area.
@ CS_TopZCoord
Top z coordinate.
@ LayeredCrossSectionInterfaceType
#define PM_CHECK_FLAG_AND_REPORT(_pm, _ir, _componentnum, _paramkey, _prio, _flag)

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