OOFEM 3.0
Loading...
Searching...
No Matches
beam2d.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 "fei2dlinehermite.h"
39#include "node.h"
40#include "material.h"
41#include "crosssection.h"
42#include "gausspoint.h"
44#include "floatmatrix.h"
45#include "intarray.h"
46#include "floatarray.h"
47#include "engngm.h"
48#include "bodyload.h"
49#include "boundaryload.h"
50#include "mathfem.h"
51#include "classfactory.h"
53#include "masterdof.h"
54#include "parametermanager.h"
55#include "paramkey.h"
56
57#ifdef __OOFEG
58 #include "oofeggraphiccontext.h"
59#endif
60
61namespace oofem {
64
65// Set up interpolation coordinates
66FEI2dLineLin Beam2d :: interp_geom(1, 3);
67FEI2dLineHermite Beam2d :: interp_beam(1, 3);
68
69Beam2d :: Beam2d(int n, Domain *aDomain) : BeamBaseElement(n, aDomain), LayeredCrossSectionInterface()
70{
72
73 kappa = -1; // set kappa to undef value (should be always > 0.)
74 length = 0.;
75 pitch = 10.; // a dummy value
77
78 ghostNodes [ 0 ] = ghostNodes [ 1 ] = NULL;
80}
81
82
83Beam2d :: ~Beam2d()
84{
85 delete ghostNodes [ 0 ];
86 delete ghostNodes [ 1 ];
87}
88
89
90FEInterpolation *Beam2d :: giveInterpolation() const { return & interp_geom; }
91
92
94Beam2d :: giveInterface(InterfaceType interface)
95{
96 if ( interface == LayeredCrossSectionInterfaceType ) {
97 return static_cast< LayeredCrossSectionInterface * >(this);
98 }
99
100 return NULL;
101}
102
103
104void
105Beam2d :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
106{
107 double l, ksi, kappa, c1;
108 TimeStep *tStep = domain->giveEngngModel()->giveCurrentStep();
109
110 l = this->computeLength();
111 ksi = 0.5 + 0.5 * gp->giveNaturalCoordinate(1);
112 kappa = this->giveKappaCoeff(tStep);
113 c1 = 1. + 2. * kappa;
114
115 answer.resize(3, 6);
116 answer.zero();
117
118 answer.at(1, 1) = -1. / l;
119 answer.at(1, 4) = 1. / l;
120 answer.at(2, 2) = ( 6. - 12. * ksi ) / ( l * l * c1 );
121 answer.at(2, 3) = ( -2. * ( 2. + kappa ) + 6. * ksi ) / ( l * c1 );
122 answer.at(2, 5) = ( -6. + 12. * ksi ) / ( l * l * c1 );
123 answer.at(2, 6) = ( -2. * ( 1. - kappa ) + 6. * ksi ) / ( l * c1 );
124 answer.at(3, 2) = ( -2. * kappa ) / ( l * c1 );
125 answer.at(3, 3) = kappa / ( c1 );
126 answer.at(3, 5) = 2. * kappa / ( l * c1 );
127 answer.at(3, 6) = kappa / ( c1 );
128}
129
130
131void
132Beam2d :: computeGaussPoints()
133{
134 if ( integrationRulesArray.size() == 0 ) {
135 // the gauss point is used only when methods from crosssection and/or material
136 // classes are requested
137 integrationRulesArray.resize(1);
138 integrationRulesArray [ 0 ] = std :: make_unique< GaussIntegrationRule >(1, this, 1, 3);
139 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], this->numberOfGaussPoints, this);
140 }
141}
142
143
144void
145Beam2d :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
146// Returns the displacement interpolation matrix {N} of the receiver, eva-
147// luated at gp. Used for numerical calculation of consistent mass
148// matrix. Must contain only interpolation for displacement terms,
149// not for any rotations. (Inertia forces do not work on rotations).
150// r = {u1,w1,fi_y1,u2,w2,fi_y2}^T
151{
152 double l, ksi, ksi2, ksi3, kappa, c1;
153 TimeStep *tStep = this->domain->giveEngngModel()->giveCurrentStep();
154
155 l = this->computeLength();
156 ksi = 0.5 + 0.5 * iLocCoord.at(1);
157 kappa = this->giveKappaCoeff(tStep);
158 c1 = 1. + 2. * kappa;
159 ksi2 = ksi * ksi;
160 ksi3 = ksi2 * ksi;
161
162 answer.resize(3, 6);
163 answer.zero();
164
165 answer.at(1, 1) = 1. - ksi;
166 answer.at(1, 4) = ksi;
167 answer.at(2, 2) = ( c1 - 2. * kappa * ksi - 3. * ksi2 + 2. * ksi3 ) / c1;
168 answer.at(2, 3) = l * ( -( 1. + kappa ) * ksi + ( 2. + kappa ) * ksi2 - ksi3 ) / c1;
169 answer.at(2, 5) = ( 2. * kappa * ksi + 3. * ksi2 - 2. * ksi3 ) / c1;
170 answer.at(2, 6) = l * ( kappa * ksi + ( 1. - kappa ) * ksi2 - ksi3 ) / c1;
171 answer.at(3, 2) = ( 6. * ksi - 6. * ksi2 ) / ( l * c1 );
172 answer.at(3, 3) = ( c1 - 2. * ( 2. + kappa ) * ksi + 3. * ksi2 ) / c1;
173 answer.at(3, 5) = ( -6. * ksi + 6. * ksi2 ) / ( l * c1 );
174 answer.at(3, 6) = ( -2. * ( 1. - kappa ) * ksi + 3. * ksi2 ) / c1;
175}
176
177
178void
179Beam2d :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
180{
181 double l = this->computeLength();
182 FloatMatrix B, d, DB;
183 answer.clear();
184 for ( auto &gp : *this->giveDefaultIntegrationRulePtr() ) {
185 this->computeBmatrixAt(gp, B);
186 this->computeConstitutiveMatrixAt(d, rMode, gp, tStep);
187 double dV = gp->giveWeight() * 0.5 * l;
188 DB.beProductOf(d, B);
189 answer.plusProductSymmUpper(B, DB, dV);
190 }
191 answer.symmetrized();
192}
193
194
195void
196Beam2d :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
197{
198 answer = this->giveStructuralCrossSection()->give2dBeamStiffMtrx(rMode, gp, tStep);
199}
200
201
202void
203Beam2d :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
204{
205 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Beam2d(strain, gp, tStep);
206}
207
208
209bool
210Beam2d :: computeGtoLRotationMatrix(FloatMatrix &answer)
211{
212 double sine, cosine;
213
214 int ndofs = computeNumberOfGlobalDofs();
215 answer.resize(ndofs, ndofs);
216 answer.zero();
217
218 sine = sin( this->givePitch() );
219 cosine = cos(pitch);
220 answer.at(1, 1) = cosine;
221 answer.at(1, 2) = sine;
222 answer.at(2, 1) = -sine;
223 answer.at(2, 2) = cosine;
224 answer.at(3, 3) = 1.;
225 answer.at(4, 4) = cosine;
226 answer.at(4, 5) = sine;
227 answer.at(5, 4) = -sine;
228 answer.at(5, 5) = cosine;
229 answer.at(6, 6) = 1.;
230
231 for ( int i = 7; i <= ndofs; i++ ) {
232 answer.at(i, i) = 1.0;
233 }
234
235 if ( this->hasDofs2Condense() ) {
236 int condensedDofCounter = 0;
237 DofIDItem dofids[] = {
238 D_u, D_w, R_v
239 };
240 FloatMatrix l2p(6, ndofs); // local -> primary
241 l2p.zero();
242 // loop over nodes
243 for ( int inode = 0; inode < 2; inode++ ) {
244 // loop over DOFs
245 for ( int idof = 0; idof < 3; idof++ ) {
246 int eq = inode * 3 + idof + 1;
247 if ( ghostNodes [ inode ] ) {
248 if ( ghostNodes [ inode ]->hasDofID(dofids [ idof ]) ) {
249 condensedDofCounter++;
250 l2p.at(eq, 6 + condensedDofCounter) = 1.0;
251 continue;
252 }
253 }
254 l2p.at(eq, eq) = 1.0;
255 }
256 }
257
258 FloatMatrix g2l(answer);
259 answer.beProductOf(l2p, g2l);
260 }
261
262 return true;
263}
264
265
266double
267Beam2d :: computeVolumeAround(GaussPoint *gp)
268{
269 return 0.5 * this->computeLength() * gp->giveWeight();
270}
271
272
273void
274Beam2d :: computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
275//
276// returns full 3d strain vector of given layer (whose z-coordinate from center-line is
277// stored in slaveGp) for given tStep
278//
279{
280 double layerZeta, layerZCoord, top, bottom;
281
282 top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
283 bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
284 layerZeta = slaveGp->giveNaturalCoordinate(3);
285 layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
286
287 answer.resize(2); // {Exx,GMzx}
288
289 answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(2) * layerZCoord;
290 answer.at(2) = masterGpStrain.at(3);
291}
292
293
294void
295Beam2d :: giveDofManDofIDMask(int inode, IntArray &answer) const
296{
297 answer = {
298 D_u, D_w, R_v
299 };
300}
301
302
303double
304Beam2d :: computeLength()
305{
306 double dx, dy;
307 Node *nodeA, *nodeB;
308
309 if ( length == 0. ) {
310 nodeA = this->giveNode(1);
311 nodeB = this->giveNode(2);
312 dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
313 dy = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
314 length = sqrt(dx * dx + dy * dy);
315 }
316
317 return length;
318}
319
320
321double
322Beam2d :: givePitch()
323{
324 double xA, xB, yA, yB;
325 Node *nodeA, *nodeB;
326
327 if ( pitch == 10. ) { // 10. : dummy initialization value
328 nodeA = this->giveNode(1);
329 nodeB = this->giveNode(2);
330 xA = nodeA->giveCoordinate(1);
331 xB = nodeB->giveCoordinate(1);
332 yA = nodeA->giveCoordinate(3);
333 yB = nodeB->giveCoordinate(3);
334 pitch = atan2(yB - yA, xB - xA);
335 }
336
337 return pitch;
338}
339
340
341double
342Beam2d :: giveKappaCoeff(TimeStep *tStep)
343{
344 // kappa = (6*E*I)/(k*G*A*l^2)
345
346 if ( kappa < 0. ) {
347 FloatMatrix d;
348 double l = this->computeLength();
349
350 this->computeConstitutiveMatrixAt(d, ElasticStiffness, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
351 kappa = 6. * d.at(2, 2) / ( d.at(3, 3) * l * l );
352 }
353
354 return kappa;
355}
356
357
358int
359Beam2d :: giveLocalCoordinateSystem(FloatMatrix &answer)
360//
361// returns a unit vectors of local coordinate system at element
362// stored rowwise (mainly used by some materials with ortho and anisotrophy)
363//
364{
365 double sine, cosine;
366
367 answer.resize(3, 3);
368 answer.zero();
369
370 sine = sin( this->givePitch() );
371 cosine = cos(pitch);
372
373 answer.at(1, 1) = cosine;
374 answer.at(1, 2) = sine;
375 answer.at(2, 1) = -sine;
376 answer.at(2, 2) = cosine;
377 answer.at(3, 3) = 1.0;
378
379 return 1;
380}
381
382
383void
384Beam2d :: initializeFrom(InputRecord &ir, int priority)
385{
386 // first call parent
387 BeamBaseElement :: initializeFrom(ir, priority);
388 ParameterManager &ppm = giveDomain()->elementPPM;
390}
391
392void
393Beam2d :: postInitialize()
394{
395 BeamBaseElement :: postInitialize();
396 ParameterManager &ppm = giveDomain()->elementPPM;
397 if (ppm.checkIfSet(this->number, IPK_Beam2d_dofsToCondense.getIndex())) {
398 auto _val = ppm.getTempParam(this->number, IPK_Beam2d_dofsToCondense.getIndex());
399 IntArray val (std::get<IntArray>(*_val));
400
401 if ( val.giveSize() >= 6 ) {
402 throw ComponentInputException(IPK_Beam2d_dofsToCondense.getName(), ComponentInputException::ComponentType::ctElement, this->number, "wrong input data for condensed dofs");
403 }
404
405 DofIDItem mask[] = {
406 D_u, D_w, R_v
407 };
408 this->numberOfCondensedDofs = val.giveSize();
409 for ( int i = 1; i <= val.giveSize(); i++ ) {
410 if ( val.at(i) <= 3 ) {
411 if ( ghostNodes [ 0 ] == NULL ) {
412 ghostNodes [ 0 ] = new ElementDofManager(1, giveDomain(), this);
413 }
414 ghostNodes [ 0 ]->appendDof( new MasterDof(ghostNodes [ 0 ], mask [ val.at(i) - 1 ]) );
415 } else {
416 if ( ghostNodes [ 1 ] == NULL ) {
417 ghostNodes [ 1 ] = new ElementDofManager(1, giveDomain(), this);
418 }
419 ghostNodes [ 1 ]->appendDof( new MasterDof(ghostNodes [ 1 ], mask [ val.at(i) - 4 ]) );
420 }
421 }
422 }
423}
424
425void
426Beam2d :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
427{
428 BeamBaseElement :: giveInternalForcesVector(answer, tStep, useUpdatedGpRecord);
429}
430
431
432void
433Beam2d :: giveEndForcesVector(FloatArray &answer, TimeStep *tStep)
434{
435 // stress equivalent vector in nodes (vector of internal forces)
436 FloatArray load;
437
438 this->giveInternalForcesVector(answer, tStep, false);
439
440 // subtract exact end forces due to nonnodal loading
441 this->computeLocalForceLoadVector(load, tStep, VM_Total);
442 if ( load.isNotEmpty() ) {
443 answer.subtract(load);
444 }
445}
446
447
448void
449Beam2d :: computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
450{
451 answer.clear();
452
453 if ( edge != 1 ) {
454 OOFEM_ERROR("Beam2D only has 1 edge (the midline) that supports loads. Attempted to apply load to edge %d", edge);
455 }
456
457 if ( type != ExternalForcesVector ) {
458 return;
459 }
460
461 double l = this->computeLength();
462 FloatArray coords, t;
463 FloatMatrix N, T;
464
465 answer.clear();
466 for ( auto &gp : *this->giveDefaultIntegrationRulePtr() ) {
467 const FloatArray &lcoords = gp->giveNaturalCoordinates();
468 this->computeNmatrixAt(lcoords, N);
469 if ( load->giveFormulationType() == Load :: FT_Entity ) {
470 load->computeValues(t, tStep, lcoords, { D_u, D_w, R_v }, mode);
471 } else {
472 this->computeGlobalCoordinates(coords, lcoords);
473 load->computeValues(t, tStep, coords, { D_u, D_w, R_v }, mode);
474 }
475
476 if ( load->giveCoordSystMode() == Load :: CST_Global ) {
477 if ( this->computeLoadGToLRotationMtrx(T) ) {
478 t.rotatedWith(T, 'n');
479 }
480 }
481
482 double dl = gp->giveWeight() * 0.5 * l;
483 answer.plusProduct(N, t, dl);
484 }
485
486 if ( global ) {
487 // Loads from sets expects global c.s.
489 answer.rotatedWith(T, 't');
490 }
491}
492
493
494void
495Beam2d :: computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
496{
497 FloatArray lc(1);
498 BeamBaseElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
499 answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
500}
501
502
503int
504Beam2d :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
505{
506 if ( type == IST_BeamForceMomentTensor ) {
507 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
508 return 1;
509 } else if ( type == IST_BeamStrainCurvatureTensor ) {
510 answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
511 return 1;
512 } else if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
513 // Order in generalized strain is:
514 // {\eps_x, \gamma_xz, \kappa_y}
515 const FloatArray &help = type == IST_ShellForceTensor ?
518
519 answer.resize(6);
520 answer.at(1) = help.at(1); // nx
521 answer.at(2) = 0.0; // ny
522 answer.at(3) = 0.0; // nz
523 answer.at(4) = 0.0; // vyz
524 answer.at(5) = help.at(2); // vxz
525 answer.at(6) = 0.0; // vxy
526 return 1;
527 } else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
528 const FloatArray &help = type == IST_ShellMomentTensor ?
531 answer.resize(6);
532 answer.at(1) = 0.0; // mx
533 answer.at(2) = 0.0; // my
534 answer.at(3) = 0.0; // mz
535 answer.at(4) = 0.0; // myz
536 answer.at(5) = 0.0; // mxz
537 answer.at(6) = help.at(3); // mxy
538 return 1;
539 } else {
540 return BeamBaseElement :: giveIPValue(answer, gp, type, tStep);
541 }
542}
543
544
545void
546Beam2d :: printOutputAt(FILE *File, TimeStep *tStep)
547{
548 FloatArray rl, Fl;
549
550 fprintf( File, "beam element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
551
552 // ask for global element displacement vector
553 this->computeVectorOf(VM_Total, tStep, rl);
554
555 // ask for global element end forces vector
556 this->giveEndForcesVector(Fl, tStep);
557
558 fprintf(File, " local displacements ");
559 for ( auto &val : rl ) {
560 fprintf(File, " %.4e", val);
561 }
562
563 fprintf(File, "\n local end forces ");
564 for ( auto &val : Fl ) {
565 fprintf(File, " %.4e", val);
566 }
567
568 fprintf(File, "\n");
569
570 for ( auto &iRule : integrationRulesArray ) {
571 iRule->printOutputAt(File, tStep);
572 }
573}
574
575
576void
577Beam2d :: computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity)
578{
579 // computes mass matrix of the receiver
580
581 FloatMatrix stiff;
582 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
583
584 /*
585 * StructuralElement::computeMassMatrix(answer, tStep);
586 * answer.times(this->giveCrossSection()->give('A'));
587 */
588 double l = this->computeLength();
589 double kappa = this->giveKappaCoeff(tStep);
590 double kappa2 = kappa * kappa;
591
592 double density = this->giveStructuralCrossSection()->give('d', gp); // constant density assumed
593 if ( ipDensity != NULL ) {
594 // Override density if desired
595 density = * ipDensity;
596 }
597
598 double area = this->giveCrossSection()->give(CS_Area, gp); // constant area assumed
599 double c2 = ( area * density ) / ( ( 1. + 2. * kappa ) * ( 1. + 2. * kappa ) );
600 double c1 = ( area * density );
601
602 answer.resize(6, 6);
603 answer.zero();
604
605 answer.at(1, 1) = c1 * l / 3.0;
606 answer.at(1, 4) = c1 * l / 6.0;
607 answer.at(2, 2) = c2 * l * ( 13. / 35. + 7. * kappa / 5. + 4. * kappa2 / 3. );
608 answer.at(2, 3) = -c2 * l * l * ( 11. / 210. + kappa * 11. / 60. + kappa2 / 6. );
609 answer.at(2, 5) = c2 * l * ( 9. / 70. + kappa * 3. / 5. + kappa2 * 2. / 3. );
610 answer.at(2, 6) = c2 * l * l * ( 13. / 420. + kappa * 3. / 20. + kappa2 / 6. );
611 answer.at(3, 3) = c2 * l * l * l * ( 1. / 105. + kappa / 30. + kappa2 / 30. );
612 answer.at(3, 5) = -c2 * l * l * ( 13. / 420. + kappa * 3. / 20. + kappa2 / 6. );
613 answer.at(3, 6) = -c2 * l * l * l * ( 1. / 140. + kappa / 30. + kappa2 / 30. );
614
615 answer.at(4, 4) = c1 * l / 3.0;
616 answer.at(5, 5) = c2 * l * ( 13. / 35. + kappa * 7. / 5. + kappa2 * 4. / 3. );
617 answer.at(5, 6) = c2 * l * l * ( 11. / 210. + kappa * 11. / 60. + kappa2 / 6. );
618 answer.at(6, 6) = c2 * l * l * l * ( 1. / 105. + kappa / 30. + kappa2 / 30. );
619
620 answer.symmetrized();
621
622 mass = l * area * density;
623}
624
625
626void
627Beam2d :: computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
628{
629 // computes initial stress matrix of receiver (or geometric stiffness matrix)
630
631 FloatMatrix stiff;
632 FloatArray endForces;
633
634 double l = this->computeLength();
635 double kappa = this->giveKappaCoeff(tStep);
636 double kappa2 = kappa * kappa;
637 double N;
638
639 answer.resize(6, 6);
640 answer.zero();
641
642 answer.at(2, 2) = 4. * kappa2 + 4. * kappa + 6. / 5.;
643 answer.at(2, 3) = -l / 10.;
644 answer.at(2, 5) = -4. * kappa2 - 4. * kappa - 6. / 5.;
645 answer.at(2, 6) = -l / 10.;
646 answer.at(3, 3) = l * l * ( kappa2 / 3. + kappa / 3. + 2. / 15. );
647 answer.at(3, 5) = l / 10.;
648 answer.at(3, 6) = -l * l * ( kappa2 / 3. + kappa / 3. + 1. / 30. );
649
650 answer.at(5, 5) = 4. * kappa2 + 4. * kappa + 6. / 5.;
651 answer.at(5, 6) = l / 10.;
652 answer.at(6, 6) = l * l * ( kappa2 / 3. + kappa / 3. + 2. / 15. );
653
654 answer.at(1, 1) = min( fabs( answer.at(2, 2) ), fabs( answer.at(3, 3) ) ) / 1000.;
655 answer.at(1, 4) = -answer.at(1, 1);
656 answer.at(4, 4) = answer.at(1, 1);
657
658 answer.symmetrized();
659 // ask end forces in g.c.s
660 this->giveEndForcesVector(endForces, tStep);
661
662 N = ( -endForces.at(1) + endForces.at(4) ) / 2.;
663 answer.times( N / ( l * ( 1. + 2. * kappa ) * ( 1. + 2. * kappa ) ) );
664
665 //answer.beLumpedOf (mass);
666}
667
668void
669Beam2d :: computeLumpedInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
670{
671 // computes initial stress matrix of receiver (or geometric stiffness matrix)
672
673 FloatMatrix stiff;
674 FloatArray endForces;
675
676 double l = this->computeLength();
677 double N;
678
679 answer.resize(6, 6);
680 answer.zero();
681
682 answer.at(2, 2) = 1.;
683 answer.at(2, 5) =-1.;
684 answer.at(5, 2) =-1.;
685 answer.at(5, 5) = 1.;
686
687 // ask end forces in g.c.s
688 this->giveEndForcesVector(endForces, tStep);
689 N = ( -endForces.at(1) + endForces.at(4) ) / 2.;
690 answer.times( N / l);
691}
692
693
694
695#ifdef __OOFEG
696void Beam2d :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
697{
698 GraphicObj *go;
699
700 if ( !gc.testElementGraphicActivity(this) ) {
701 return;
702 }
703
704 // if (!go) { // create new one
705 WCRec p [ 2 ]; /* poin */
706 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
707 EASValsSetColor( gc.getElementColor() );
708 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
709 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
710 p [ 0 ].y = 0.;
711 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
712 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
713 p [ 1 ].y = 0.;
714 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
715 go = CreateLine3D(p);
716 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
717 EGAttachObject(go, ( EObjectP ) this);
718 EMAddGraphicsToModel(ESIModel(), go);
719}
720
721
722void Beam2d :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
723{
724 GraphicObj *go;
725
726 if ( !gc.testElementGraphicActivity(this) ) {
727 return;
728 }
729
730 double defScale = gc.getDefScale();
731 // if (!go) { // create new one
732 WCRec p [ 2 ]; /* poin */
733 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
734 EASValsSetColor( gc.getDeformedElementColor() );
735 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
736 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
737 p [ 0 ].y = 0.;
738 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
739
740 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
741 p [ 1 ].y = 0.;
742 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
743 go = CreateLine3D(p);
744 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
745 EMAddGraphicsToModel(ESIModel(), go);
746}
747#endif
748} // end namespace oofem
#define N(a, b)
#define REGISTER_Element(class)
double giveKappaCoeff(TimeStep *tStep)
Definition beam2d.C:342
int computeNumberOfGlobalDofs() override
Definition beam2d.h:102
DofManager * ghostNodes[2]
Definition beam2d.h:72
double computeLength() override
Definition beam2d.C:304
double length
Definition beam2d.h:65
void giveEndForcesVector(FloatArray &answer, TimeStep *tStep)
Definition beam2d.C:433
bool computeGtoLRotationMatrix(FloatMatrix &answer) override
Definition beam2d.C:210
double givePitch()
Definition beam2d.C:322
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override
Definition beam2d.C:426
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &) override
Definition beam2d.C:145
static FEI2dLineLin interp_geom
Definition beam2d.h:76
double kappa
Definition beam2d.h:65
static ParamKey IPK_Beam2d_dofsToCondense
[optional] DOFs to condense
Definition beam2d.h:79
double pitch
Definition beam2d.h:65
int numberOfCondensedDofs
number of condensed DOFs
Definition beam2d.h:74
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition beam2d.C:196
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
Definition beam2d.C:105
bool hasDofs2Condense()
Definition beam2d.h:174
virtual void computeLocalForceLoadVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
BeamBaseElement(int n, Domain *d)
CoordSystType giveCoordSystMode() override
double giveCoordinate(int i) const
Definition dofmanager.h:383
Node * giveNode(int i) const
Definition element.h:629
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Definition element.C:1225
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
int numberOfGaussPoints
Definition element.h:175
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
int giveLabel() const
Definition element.h:1125
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
int giveNumber() const
Definition femcmpnn.h:104
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Definition floatarray.C:288
void rotatedWith(FloatMatrix &r, char mode)
Definition floatarray.C:814
void subtract(const FloatArray &src)
Definition floatarray.C:320
void times(double s)
Definition floatarray.C:834
bool isNotEmpty() const
Returns true if receiver is not empty.
Definition floatarray.h:263
void times(double f)
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
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 beProductOf(const FloatMatrix &a, const FloatMatrix &b)
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
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
virtual void computeValues(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, const IntArray &dofids, ValueModeType mode)
Definition load.C:59
virtual FormulationType giveFormulationType()
Definition load.h:170
bool checkIfSet(size_t componentIndex, size_t paramIndex)
std::optional< paramValue > getTempParam(size_t componentIndex, size_t paramIndex) const
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
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.
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ 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
#define PM_UPDATE_TEMP_PARAMETER(_type, _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