OOFEM 3.0
Loading...
Searching...
No Matches
libeam2dnl.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
48#ifdef __OOFEG
49 #include "oofeggraphiccontext.h"
50#endif
51
52namespace oofem {
54
56{
58 length = 0.;
59 pitch = 10.; // a dummy value
60}
61
64{
65 if ( interface == LayeredCrossSectionInterfaceType ) {
66 return static_cast< LayeredCrossSectionInterface * >( this );
67 }
68
69 return NULL;
70}
71
72
73void
75// Returns the strain matrix of the receiver.
76{
77 double l, ksi, n1, n2, n1x, n2x, n3x;
78
79 l = this->computeLength();
80 ksi = gp->giveNaturalCoordinate(1);
81
82 n1 = 0.5 * ( 1.0 - ksi );
83 n2 = 0.5 * ( 1.0 + ksi );
84 n1x = -1.0 / l;
85 n2x = 1.0 / l;
86 n3x = -4.0 * ksi / l;
87
88 answer.resize(3, 6);
89 answer.zero();
90
91 // eps_x
92 answer.at(1, 1) = n1x;
93 answer.at(1, 4) = n2x;
94 // kappa
95 answer.at(2, 3) = n1x;
96 answer.at(2, 6) = n2x;
97 // gamma_xz
98 answer.at(3, 2) = n1x;
99 answer.at(3, 3) = n1 - n3x * l / 8.0;
100 answer.at(3, 5) = n2x;
101 answer.at(3, 6) = n2 + n3x * l / 8.0;
102}
103
104
105void
107{
108 //
109 // Returns nonlinear part of geometrical equations of the receiver at gp.
110 //
111 // Returns A matrix (see Bittnar & Sejnoha Num. Met. Mech. part II, chap 9)
112 double l, l8, ll88, ksi, n1x, n2x, n3x;
113
114 l = this->computeLength();
115 ksi = gp->giveNaturalCoordinate(1);
116
117 //n1 = 0.5*(1.0 - ksi);
118 //n2 = 0.5*(1.0 + ksi);
119 n1x = -1.0 / l;
120 n2x = 1.0 / l;
121 n3x = -4.0 * ksi / l;
122 l8 = l / 8.;
123 ll88 = l8 * l8;
124
125 answer.resize(6, 6);
126 answer.zero();
127 /*
128 * if (i==1){
129 * answer.at(1,1)= n1x*n1x;
130 * answer.at(1,4)= n1x*n2x;
131 * answer.at(2,2)= n1x*n1x;
132 * answer.at(2,3)=-n1x*n3x*l8;
133 * answer.at(2,5)= n1x*n2x;
134 * answer.at(2,6)= n1x*n3x*l8;
135 * answer.at(3,2)=-n3x*l8*n1x;
136 * answer.at(3,3)= n3x*n3x*ll88;
137 * answer.at(3,5)=-n3x*l8*n2x;
138 * answer.at(3,6)=-n3x*n3x*ll88;
139 *
140 * answer.at(4,1)= n2x*n1x;
141 * answer.at(4,4)= n2x*n2x;
142 * answer.at(5,2)= n2x*n1x;
143 * answer.at(5,3)=-n2x*n3x*l8;
144 * answer.at(5,5)= n2x*n2x;
145 * answer.at(5,6)= n2x*n3x*l8;
146 * answer.at(6,2)= n3x*l8*n1x;
147 * answer.at(6,3)=-n3x*n3x*ll88;
148 * answer.at(6,5)= n3x*l8*n2x;
149 * answer.at(6,2)= n3x*n3x*ll88;
150 *
151 * }
152 *
153 * if (i==2){
154 * answer.at(1,2)=n1x*n1x;
155 * answer.at(1,5)=n1x*n2x;
156 * answer.at(2,1)=n1x*n1x;
157 * answer.at(2,4)=n1x*n2x;
158 * answer.at(4,2)=n2x*n1x;
159 * answer.at(4,5)=n2x*n2x;
160 * answer.at(5,1)=n2x*n1x;
161 * answer.at(5,4)=n2x*n2x;
162 *
163 * }
164 *
165 * if (i==3){
166 * answer.at(1,3)=n1x*n1;
167 * answer.at(1,6)=n1x*n2;
168 * answer.at(3,1)=n1*n1x;
169 * answer.at(3,4)=n1*n2x;
170 * answer.at(4,3)=n2x*n1;
171 * answer.at(4,6)=n2x*n2;
172 * answer.at(6,1)=n2*n1x;
173 * answer.at(6,4)=n2*n2x;
174 * }
175 */
176
177 /* Working
178 *
179 * // large deflection
180 * // based on von Karman equation
181 *
182 * if (i==1) {
183 * answer.at(2,2)= n1x*n1x;
184 * answer.at(2,3)=-n1x*n3x*l8;
185 * answer.at(2,5)= n1x*n2x;
186 * answer.at(2,6)= n1x*n3x*l8;
187 * answer.at(3,2)=-n3x*l8*n1x;
188 * answer.at(3,3)= n3x*n3x*ll88;
189 * answer.at(3,5)=-n3x*l8*n2x;
190 * answer.at(3,6)=-n3x*n3x*ll88;
191 *
192 * answer.at(5,2)= n2x*n1x;
193 * answer.at(5,3)=-n2x*n3x*l8;
194 * answer.at(5,5)= n2x*n2x;
195 * answer.at(5,6)= n2x*n3x*l8;
196 * answer.at(6,2)= n3x*l8*n1x;
197 * answer.at(6,3)=-n3x*n3x*ll88;
198 * answer.at(6,5)= n3x*l8*n2x;
199 * answer.at(6,2)= n3x*n3x*ll88;
200 *
201 * }
202 */
203
204 if ( i == 1 ) {
205 // answer.at(1,1)= n1x*n1x;
206 // answer.at(1,4)= n1x*n2x;
207 answer.at(2, 2) = n1x * n1x;
208 answer.at(2, 3) = -n1x * n3x * l8;
209 answer.at(2, 5) = n1x * n2x;
210 answer.at(2, 6) = n1x * n3x * l8;
211 answer.at(3, 2) = -n3x * l8 * n1x;
212 answer.at(3, 3) = n3x * n3x * ll88;
213 answer.at(3, 5) = -n3x * l8 * n2x;
214 answer.at(3, 6) = -n3x * n3x * ll88;
215
216 // answer.at(4,1)= n2x*n1x;
217 // answer.at(4,4)= n2x*n2x;
218 answer.at(5, 2) = n2x * n1x;
219 answer.at(5, 3) = -n2x * n3x * l8;
220 answer.at(5, 5) = n2x * n2x;
221 answer.at(5, 6) = n2x * n3x * l8;
222 answer.at(6, 2) = n3x * l8 * n1x;
223 answer.at(6, 3) = -n3x * n3x * ll88;
224 answer.at(6, 5) = n3x * l8 * n2x;
225 answer.at(6, 2) = n3x * n3x * ll88;
226 }
227
228 /*
229 * if (i==2){
230 * answer.at(1,2)=n1x*n1x;
231 * answer.at(1,5)=n1x*n2x;
232 * answer.at(2,1)=n1x*n1x;
233 * answer.at(2,4)=n1x*n2x;
234 * answer.at(4,2)=n2x*n1x;
235 * answer.at(4,5)=n2x*n2x;
236 * answer.at(5,1)=n2x*n1x;
237 * answer.at(5,4)=n2x*n2x;
238 *
239 * }
240 *
241 * if (i==3){
242 * answer.at(1,3)=n1x*n1;
243 * answer.at(1,6)=n1x*n2;
244 * answer.at(3,1)=n1*n1x;
245 * answer.at(3,4)=n1*n2x;
246 * answer.at(4,3)=n2x*n1;
247 * answer.at(4,6)=n2x*n2;
248 * answer.at(6,1)=n2*n1x;
249 * answer.at(6,4)=n2*n2x;
250 * }
251 */
252}
253
254
255void
257{
258 double dV;
259 FloatArray stress;
260 FloatMatrix A;
261
262 answer.resize(6, 6);
263 answer.zero();
264
265 // assemble initial stress matrix
266 for ( GaussPoint *gp: * this->giveDefaultIntegrationRulePtr() ) {
267 dV = this->computeVolumeAround(gp);
268 stress = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
269 if ( stress.giveSize() ) {
270 for ( int j = 1; j <= stress.giveSize(); j++ ) {
271 // loop over each component of strain vector
272 this->computeNLBMatrixAt(A, gp, j);
273 if ( A.isNotEmpty() ) {
274 A.times(stress.at(j) * dV);
275 answer.add(A);
276 }
277 }
278 }
279 }
280}
281
282
284// Sets up the array of Gauss Points of the receiver.
285{
286 if ( integrationRulesArray.size() == 0 ) {
287 integrationRulesArray.resize(1);
288 integrationRulesArray [ 0 ] = std::make_unique< GaussIntegrationRule >(1, this, 1, 2);
290 }
291}
292
293
294void
296// Returns the lumped mass matrix of the receiver. This expression is
297// valid in both local and global axes.
298{
299 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
300 double density = this->giveStructuralCrossSection()->give('d', gp);
301 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
302 answer.resize(6, 6);
303 answer.zero();
304 answer.at(1, 1) = halfMass;
305 answer.at(2, 2) = halfMass;
306 answer.at(4, 4) = halfMass;
307 answer.at(5, 5) = halfMass;
308}
309
310
311void
313// Returns the displacement interpolation matrix {N} of the receiver, eva-
314// luated at gp.
315{
316 double ksi, n1, n2, n3;
317 double l = this->computeLength();
318
319 ksi = iLocCoord.at(1);
320 n1 = ( 1. - ksi ) * 0.5;
321 n2 = ( 1. + ksi ) * 0.5;
322 n3 = ( 1. - ksi * ksi );
323
324 answer.resize(3, 6);
325 answer.zero();
326
327 // us
328 answer.at(1, 1) = n1;
329 answer.at(1, 4) = n2;
330 // w
331 answer.at(2, 2) = n1;
332 answer.at(2, 3) = -n3 * l / 8.;
333 answer.at(2, 5) = n2;
334 answer.at(2, 6) = n3 * l / 8.;
335 // fi_y
336 answer.at(3, 3) = n1;
337 answer.at(3, 6) = n2;
338}
339
340
341bool
343{
344 double sine, cosine;
345
346 sine = sin(this->givePitch() );
347 cosine = cos(pitch);
348
349 answer.resize(6, 6);
350 answer.zero();
351
352 answer.at(1, 1) = cosine;
353 answer.at(1, 2) = sine;
354 answer.at(2, 1) = -sine;
355 answer.at(2, 2) = cosine;
356 answer.at(3, 3) = 1.;
357 answer.at(4, 4) = cosine;
358 answer.at(4, 5) = sine;
359 answer.at(5, 4) = -sine;
360 answer.at(5, 5) = cosine;
361 answer.at(6, 6) = 1.;
362
363 return true;
364}
365
366
368// Returns the length of the receiver. This method is valid only if 1
369// Gauss point is used.
370{
371 double weight = gp->giveWeight();
372 return weight * 0.5 * this->computeLength();
373}
374
375
376int
378{
381 if ( type == IST_BeamForceMomentTensor ) {
382 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
383 return 1;
384 } else if ( type == IST_BeamStrainCurvatureTensor ) {
385 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
386 return 1;
387 } else {
388 return StructuralElement::giveIPValue(answer, gp, type, tStep);
389 }
390}
391
392
393void
395{
396 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Beam2d(strain, gp, tStep);
397}
398
399
400void
401LIBeam2dNL::computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
402{
403 answer = this->giveStructuralCrossSection()->give2dBeamStiffMtrx(rMode, gp, tStep);
404}
405
406void
408{
409 OOFEM_ERROR("computeConstitutiveMatrix_dPdF_At Not implemented for libeam2dnl");
410}
411
412
413void
414LIBeam2dNL::computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
415//
416// returns full 3d strain vector of given layer (whose z-coordinate from center-line is
417// stored in slaveGp) for given tStep
418//
419{
420 double layerZeta, layerZCoord, top, bottom;
421
422 top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
423 bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
424 layerZeta = slaveGp->giveNaturalCoordinate(3);
425 layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
426
427 answer.resize(2); // {Exx,GMzx}
428
429 answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(2) * layerZCoord;
430 answer.at(2) = masterGpStrain.at(3);
431}
432
433
434void
436{
437 answer = { D_u, D_w, R_v };
438}
439
440
441int
443{
444 double ksi, n1, n2;
445
446 ksi = lcoords.at(1);
447 n1 = ( 1. - ksi ) * 0.5;
448 n2 = ( 1. + ksi ) * 0.5;
449
450 answer.resize(3);
451 answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 * this->giveNode(2)->giveCoordinate(1);
452 answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 * this->giveNode(2)->giveCoordinate(3);
453
454 return 1;
455}
456
457
459// Returns the length of the receiver.
460{
461 double dx, dy;
462 Node *nodeA, *nodeB;
463
464 if ( length == 0. ) {
465 nodeA = this->giveNode(1);
466 nodeB = this->giveNode(2);
467 dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
468 dy = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
469 length = sqrt(dx * dx + dy * dy);
470 }
471
472 return length;
473}
474
475
477// Returns the pitch of the receiver.
478{
479 double xA, xB, yA, yB;
480 Node *nodeA, *nodeB;
481
482 if ( pitch == 10. ) { // 10. : dummy initialization value
483 nodeA = this->giveNode(1);
484 nodeB = this->giveNode(2);
485 xA = nodeA->giveCoordinate(1);
486 xB = nodeB->giveCoordinate(1);
487 yA = nodeA->giveCoordinate(3);
488 yB = nodeB->giveCoordinate(3);
489 pitch = atan2(yB - yA, xB - xA);
490 }
491
492 return pitch;
493}
494
495
496void
498{
500}
501
502
503void
505{
506 /*
507 * provides dof mapping of local edge dofs (only nonzero are taken into account)
508 * to global element dofs
509 */
510
511 if ( iEdge != 1 ) {
512 OOFEM_ERROR("wrong edge number");
513 }
514
515
516 answer.resize(6);
517 answer.at(1) = 1;
518 answer.at(2) = 2;
519 answer.at(3) = 3;
520 answer.at(4) = 4;
521 answer.at(5) = 5;
522 answer.at(6) = 6;
523}
524
525double
527{
528 if ( iEdge != 1 ) { // edge between nodes 1 2
529 OOFEM_ERROR("wrong egde number");
530 }
531
532 double weight = gp->giveWeight();
533 return 0.5 * this->computeLength() * weight;
534}
535
536
537int
539{
540 /*
541 * Returns transformation matrix from global coordinate system to local
542 * element coordinate system for element load vector components.
543 * If no transformation is necessary, answer is empty matrix (default);
544 */
545
546 // f(elemLocal) = T * f (global)
547 double sine, cosine;
548
549 answer.resize(3, 3);
550 answer.zero();
551
552 sine = sin(this->givePitch() );
553 cosine = cos(pitch);
554
555 answer.at(1, 1) = cosine;
556 answer.at(1, 2) = -sine;
557 answer.at(2, 1) = sine;
558 answer.at(2, 2) = cosine;
559 answer.at(3, 3) = 1.0;
560
561 return 1;
562}
563
564int
566{
567 // returns transformation matrix from
568 // edge local coordinate system
569 // to element local c.s
570 // (same as global c.s in this case)
571 //
572 // i.e. f(element local) = T * f(edge local)
573 //
574 answer.clear();
575 return 0;
576}
577
578
579void
580LIBeam2dNL::computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
581{
582 FloatArray lc(1);
583 StructuralElement::computeBodyLoadVectorAt(answer, load, tStep, mode);
584 answer.times(this->giveCrossSection()->give(CS_Area, lc, this) );
585}
586
587
588#ifdef __OOFEG
590{
591 GraphicObj *go;
592
593 if ( !gc.testElementGraphicActivity(this) ) {
594 return;
595 }
596
597 // if (!go) { // create new one
598 WCRec p[ 2 ]; /* poin */
599 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
600 EASValsSetColor(gc.getElementColor() );
601 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
602 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
603 p [ 0 ].y = 0.;
604 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
605 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
606 p [ 1 ].y = 0.;
607 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
608 go = CreateLine3D(p);
609 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
610 EGAttachObject(go, ( EObjectP ) this);
611 EMAddGraphicsToModel(ESIModel(), go);
612}
613
614
616{
617 GraphicObj *go;
618 if ( !gc.testElementGraphicActivity(this) ) {
619 return;
620 }
621
622 double defScale = gc.getDefScale();
623 // if (!go) { // create new one
624 WCRec p[ 2 ]; /* poin */
625 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
626 EASValsSetColor(gc.getDeformedElementColor() );
627 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
628 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
629 p [ 0 ].y = 0.;
630 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
631
632 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
633 p [ 1 ].y = 0.;
634 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
635 go = CreateLine3D(p);
636 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
637 EMAddGraphicsToModel(ESIModel(), go);
638}
639#endif
640} // 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
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
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void times(double s)
Definition floatarray.C:834
void times(double f)
void add(const FloatMatrix &a)
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.
bool isNotEmpty() const
Tests for empty matrix.
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
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition libeam2dnl.C:401
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition libeam2dnl.C:394
void computeConstitutiveMatrix_dPdF_At(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition libeam2dnl.C:407
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition libeam2dnl.C:74
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition libeam2dnl.C:312
void computeNLBMatrixAt(FloatMatrix &answer, GaussPoint *gp, int)
Definition libeam2dnl.C:106
void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition libeam2dnl.C:256
double computeEdgeVolumeAround(GaussPoint *, int) override
Definition libeam2dnl.C:526
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition libeam2dnl.C:435
void computeGaussPoints() override
Definition libeam2dnl.C:283
int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) override
Definition libeam2dnl.C:442
bool computeGtoLRotationMatrix(FloatMatrix &answer) override
Definition libeam2dnl.C:342
void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override
Definition libeam2dnl.C:615
Interface * giveInterface(InterfaceType it) override
Definition libeam2dnl.C:63
void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep) override
Definition libeam2dnl.C:589
double computeLength() override
Definition libeam2dnl.C:458
void giveEdgeDofMapping(IntArray &answer, int) const override
Definition libeam2dnl.C:504
int computeLoadGToLRotationMtrx(FloatMatrix &answer) override
Definition libeam2dnl.C:538
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition libeam2dnl.C:295
LIBeam2dNL(int n, Domain *d)
Definition libeam2dnl.C:55
void initializeFrom(InputRecord &ir, int priority) override
Definition libeam2dnl.C:497
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition libeam2dnl.C:377
int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *) override
Definition libeam2dnl.C:565
void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep) override
Definition libeam2dnl.C:414
double computeVolumeAround(GaussPoint *gp) override
Definition libeam2dnl.C:367
void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) override
Definition libeam2dnl.C:580
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< 3, 3 > give2dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatArrayF< 3 > giveGeneralizedStress_Beam2d(const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
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
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