OOFEM 3.0
Loading...
Searching...
No Matches
libeam3d2.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 "timestep.h"
46#include "contextioerr.h"
47#include "mathfem.h"
48#include "classfactory.h"
49#include "engngm.h"
50#include "parametermanager.h"
51#include "paramkey.h"
52
53#ifdef __OOFEG
54 #include "oofeggraphiccontext.h"
55 #include "oofegutils.h"
56 #include "connectivitytable.h"
57#endif
58
59namespace oofem {
62
63LIBeam3d2::LIBeam3d2(int n, Domain *aDomain) : NLStructuralElement(n, aDomain), tc(), tempTc()
64{
66 referenceNode = 0;
67 length = 0.;
68
69 tempTcCounter = 0;
70}
71
72
73void
75// Returns the strain matrix of the receiver.
76// eeps = {\eps_x, \gamma_xz, \gamma_xy, \der{phi_x}{x}, \kappa_y, \kappa_z}^T
77{
78 double l, ksi, n1, n2, n1x, n2x;
79
80 if ( this->nlGeometry ) {
81 l = this->giveCurrentLength(domain->giveEngngModel()->giveCurrentStep() );
82 } else {
83 l = this->computeLength();
84 }
85
86 ksi = gp->giveNaturalCoordinate(1);
87
88 answer.resize(6, 12);
89 answer.zero();
90
91 n1 = 0.5 * ( 1 - ksi );
92 n2 = 0.5 * ( 1. + ksi );
93 n1x = -1.0 / l;
94 n2x = 1.0 / l;
95
96 answer.at(1, 1) = -1. / l;
97 answer.at(1, 7) = 1. / l;
98
99 answer.at(2, 3) = n1x;
100 answer.at(2, 5) = n1;
101 answer.at(2, 9) = n2x;
102 answer.at(2, 11) = n2;
103
104 answer.at(3, 2) = n1x;
105 answer.at(3, 6) = -n1;
106 answer.at(3, 8) = n2x;
107 answer.at(3, 12) = -n2;
108
109 answer.at(4, 4) = -1. / l;
110 answer.at(4, 10) = 1. / l;
111
112 answer.at(5, 5) = n1x;
113 answer.at(5, 11) = n2x;
114
115 answer.at(6, 6) = n1x;
116 answer.at(6, 12) = n2x;
117}
118
119void
121// Sets up the array of Gauss Points of the receiver.
122{
123 if ( integrationRulesArray.size() == 0 ) {
124 integrationRulesArray.resize(1);
125 integrationRulesArray [ 0 ] = std::make_unique< GaussIntegrationRule >(1, this, 1, 2);
127 }
128}
129
130
131
132void
134// Returns the lumped mass matrix of the receiver. This expression is
135// valid in both local and global axes.
136{
137 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
138 double density = this->giveStructuralCrossSection()->give('d', gp);
139 double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
140 answer.resize(12, 12);
141 answer.zero();
142 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass;
143 answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass;
144}
145
146
147void
149// Returns the displacement interpolation matrix {N} of the receiver, eva-
150// luated at gp.
151{
152 double ksi, n1, n2;
153
154 ksi = iLocCoord.at(1);
155 n1 = ( 1. - ksi ) * 0.5;
156 n2 = ( 1. + ksi ) * 0.5;
157
158 answer.resize(6, 12);
159 answer.zero();
160
161 answer.at(1, 1) = n1;
162 answer.at(1, 7) = n2;
163 answer.at(2, 2) = n1;
164 answer.at(2, 8) = n2;
165 answer.at(3, 3) = n1;
166 answer.at(3, 9) = n2;
167
168 answer.at(4, 4) = n1;
169 answer.at(4, 10) = n2;
170 answer.at(5, 5) = n1;
171 answer.at(5, 11) = n2;
172 answer.at(6, 6) = n1;
173 answer.at(6, 12) = n2;
174}
175
176
177void
178LIBeam3d2::computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
179// Returns the stiffness matrix of the receiver, expressed in the global
180// axes.
181{
182 StructuralElement::computeStiffnessMatrix(answer, rMode, tStep);
183}
184
185
186bool
188{
189 answer.resize(12, 12);
190 answer.zero();
191
192 if ( this->nlGeometry == 0 ) {
193 FloatMatrix lcs;
194
195 this->giveLocalCoordinateSystem(lcs);
196 for ( int i = 1; i <= 3; i++ ) {
197 for ( int j = 1; j <= 3; j++ ) {
198 answer.at(i, j) = lcs.at(i, j);
199 answer.at(i + 3, j + 3) = lcs.at(i, j);
200 answer.at(i + 6, j + 6) = lcs.at(i, j);
201 answer.at(i + 9, j + 9) = lcs.at(i, j);
202 }
203 }
204 } else {
205 this->updateTempTriad(domain->giveEngngModel()->giveCurrentStep() );
206
207 for ( int i = 1; i <= 3; i++ ) {
208 for ( int j = 1; j <= 3; j++ ) {
209 answer.at(i, j) = tempTc.at(j, i);
210 answer.at(i + 3, j + 3) = tempTc.at(j, i);
211 answer.at(i + 6, j + 6) = tempTc.at(j, i);
212 answer.at(i + 9, j + 9) = tempTc.at(j, i);
213 }
214 }
215 }
216
217 return 1;
218}
219
220
221double
223// Returns the length of the receiver. This method is valid only if 1
224// Gauss point is used.
225{
226 double weight = gp->giveWeight();
227 return weight * 0.5 * this->computeLength();
228}
229
230
231void
233{
234 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
235}
236
237
238int
240{
241 double ksi, n1, n2;
242
243 ksi = lcoords.at(1);
244 n1 = ( 1. - ksi ) * 0.5;
245 n2 = ( 1. + ksi ) * 0.5;
246
247 answer.resize(3);
248 answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 * this->giveNode(2)->giveCoordinate(1);
249 answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 * this->giveNode(2)->giveCoordinate(2);
250 answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 * this->giveNode(2)->giveCoordinate(3);
251
252 return 1;
253}
254
255
256void
257LIBeam3d2::computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
258{
259 answer = this->giveStructuralCrossSection()->give3dBeamStiffMtrx(rMode, gp, tStep);
260}
261
262void
264{
265 OOFEM_ERROR("computeConstitutiveMatrix_dPdF_At Not implemented");
266}
267
268
269void
271{
272 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Beam3d(strain, gp, tStep);
273}
274
275
276int
278{
279 return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 );
280}
281
282
283double
285// Returns the length of the receiver.
286{
287 double dx, dy, dz;
288 Node *nodeA, *nodeB;
289
290 if ( length == 0. ) {
291 nodeA = this->giveNode(1);
292 nodeB = this->giveNode(2);
293 dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
294 dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
295 dz = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
296 length = sqrt(dx * dx + dy * dy + dz * dz);
297 }
298
299 return length;
300}
301
302
303void
305{
306 // first call parent
308 ParameterManager &ppm = domain->elementPPM;
310}
311
312void
314{
315 ParameterManager &ppm = domain->elementPPM;
316
319 if ( referenceNode == 0 ) {
320 OOFEM_ERROR("wrong reference node specified.");
321 }
322 FloatMatrix lcs;
323 this->giveLocalCoordinateSystem(lcs);
324 this->tc.beTranspositionOf(lcs);
325}
326
327void
329{
330 /*
331 * provides dof mapping of local edge dofs (only nonzero are taken into account)
332 * to global element dofs
333 */
334 if ( iEdge != 1 ) {
335 OOFEM_ERROR("wrong edge number");
336 }
337
338
339 answer.resize(12);
340 for ( int i = 1; i <= 12; i++ ) {
341 answer.at(i) = i;
342 }
343}
344
345
346double
348{
349 if ( iEdge != 1 ) { // edge between nodes 1 2
350 OOFEM_ERROR("wrong egde number");
351 }
352
353 double weight = gp->giveWeight();
354 return 0.5 * this->computeLength() * weight;
355}
356
357
358int
360{
361 /*
362 * Returns transformation matrix from global coordinate system to local
363 * element coordinate system for element load vector components.
364 * If no transformation is necessary, answer is empty matrix (default);
365 */
366
367 // f(elemLocal) = T * f (global)
368
369 FloatMatrix lcs;
370
371 answer.resize(6, 6);
372 answer.zero();
373
374 this->giveLocalCoordinateSystem(lcs);
375 for ( int i = 1; i <= 3; i++ ) {
376 for ( int j = 1; j <= 3; j++ ) {
377 answer.at(i, j) = lcs.at(i, j);
378 answer.at(3 + i, 3 + j) = lcs.at(i, j);
379 }
380 }
381
382 return 1;
383}
384
385
386void
387LIBeam3d2::computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
388{
389 FloatArray lc(1);
390 NLStructuralElement::computeBodyLoadVectorAt(answer, load, tStep, mode);
391 answer.times(this->giveCrossSection()->give(CS_Area, lc, this) );
392}
393
394
395int
397{
398 // returns transformation matrix from
399 // edge local coordinate system
400 // to element local c.s
401 // (same as global c.s in this case)
402 //
403 // i.e. f(element local) = T * f(edge local)
404 //
405 answer.clear();
406 return 0;
407}
408
409
410int
412//
413// returns a unit vectors of local coordinate system at element
414// stored rowwise (mainly used by some materials with ortho and anisotrophy)
415//
416{
417 FloatArray lx(3), ly(3), lz(3), help(3);
418 double length = this->computeLength();
419 Node *nodeA, *nodeB, *refNode;
420
421 answer.resize(3, 3);
422 answer.zero();
423 nodeA = this->giveNode(1);
424 nodeB = this->giveNode(2);
425 refNode = this->giveDomain()->giveNode(this->referenceNode);
426
427 for ( int i = 1; i <= 3; i++ ) {
428 lx.at(i) = ( nodeB->giveCoordinate(i) - nodeA->giveCoordinate(i) ) / length;
429 help.at(i) = ( refNode->giveCoordinate(i) - nodeA->giveCoordinate(i) );
430 }
431
432 lz.beVectorProductOf(lx, help);
433 lz.normalize();
434 ly.beVectorProductOf(lz, lx);
435 ly.normalize();
436
437 for ( int i = 1; i <= 3; i++ ) {
438 answer.at(1, i) = lx.at(i);
439 answer.at(2, i) = ly.at(i);
440 answer.at(3, i) = lz.at(i);
441 }
442
443 return 1;
444}
445
446
447void
449{
450 // test if not previously done
451 if ( tStep->giveSolutionStateCounter() == tempTcCounter ) {
452 return;
453 }
454
455 FloatArray u, centreSpin(3);
456 FloatMatrix dR(3, 3);
457
458 // ask element's displacement increments
459 this->computeVectorOf(VM_Incremental, tStep, u);
460
461 // interpolate spin at the centre
462 centreSpin.at(1) = 0.5 * ( u.at(4) + u.at(10) );
463 centreSpin.at(2) = 0.5 * ( u.at(5) + u.at(11) );
464 centreSpin.at(3) = 0.5 * ( u.at(6) + u.at(12) );
465
466 // compute rotation matrix from centreSpin pseudovector
467 this->computeRotMtrx(dR, centreSpin);
468
469 // update triad
470 tempTc.beProductOf(dR, tc);
471 // remember timestamp
473}
474
475
476void
478{
479 FloatMatrix S(3, 3), SS(3, 3);
480 double psiSize;
481
482 if ( psi.giveSize() != 3 ) {
483 OOFEM_ERROR("psi param size mismatch");
484 }
485
486 answer.resize(3, 3);
487 answer.zero();
488
489 psiSize = psi.computeNorm();
490 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 1.;
491
492 if ( psiSize <= 1.e-40 ) {
493 return;
494 }
495
496 this->computeSMtrx(S, psi);
497 SS.beProductOf(S, S);
498 S.times(sin(psiSize) / psiSize);
499 SS.times( ( 1. - cos(psiSize) ) / ( psiSize * psiSize ) );
500
501 answer.add(S);
502 answer.add(SS);
503}
504
505
506void
508{
509 if ( vec.giveSize() != 3 ) {
510 OOFEM_ERROR("vec param size mismatch");
511 }
512
513 answer.resize(3, 3);
514
515 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 0.;
516 answer.at(1, 2) = -vec.at(3);
517 answer.at(1, 3) = vec.at(2);
518 answer.at(2, 1) = vec.at(3);
519 answer.at(2, 3) = -vec.at(1);
520 answer.at(3, 1) = -vec.at(2);
521 answer.at(3, 2) = vec.at(1);
522}
523
524
525void
527{
528 FloatArray ui, PrevEpsilon;
529 FloatMatrix b;
530
531 this->computeVectorOf(VM_Incremental, tStep, ui);
532
533 this->computeBmatrixAt(gp, b);
534 // increment of strains
535 answer.beProductOf(b, ui);
536
537 if ( this->nlGeometry ) {
538 // add increment to previous total value
539 PrevEpsilon = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
540 if ( PrevEpsilon.giveSize() ) {
541 answer.add(PrevEpsilon);
542 }
543 }
544}
545
546
547double
549// Returns the length of the receiver.
550{
551 double dx, dy, dz;
552 FloatArray u;
553
554 this->computeVectorOf(VM_Total, tStep, u);
555
556 dx = ( this->giveNode(2)->giveCoordinate(1) + u.at(7) ) -
557 ( this->giveNode(1)->giveCoordinate(1) + u.at(1) );
558 dy = ( this->giveNode(2)->giveCoordinate(2) + u.at(8) ) -
559 ( this->giveNode(1)->giveCoordinate(2) + u.at(2) );
560 dz = ( this->giveNode(2)->giveCoordinate(3) + u.at(9) ) -
561 ( this->giveNode(1)->giveCoordinate(3) + u.at(3) );
562
563 return sqrt(dx * dx + dy * dy + dz * dz);
564}
565
566
568// Updates the receiver at end of step.
569{
571
572 // update triad
573 this->updateTempTriad(tStep);
574 tc = tempTc;
575
576 // update curvature
577 //FloatArray curv;
578 //this->computeTempCurv (curv, tStep);
579 //kappa = curv;
580}
581
582
583void
585// initializes receiver to new time step or can be used
586// if current time step must be restarted
587{
589 tempTc = tc;
590}
591
592
593
595{
598
599 if ( ( iores = tc.storeYourself(stream) ) != CIO_OK ) {
600 THROW_CIOERR(iores);
601 }
602}
603
604
606{
609
610 if ( ( iores = tc.restoreYourself(stream) ) != CIO_OK ) {
611 THROW_CIOERR(iores);
612 }
613}
614
615
616void
618 GaussPoint *slaveGp, TimeStep *tStep)
619{
620 double layerYCoord, layerZCoord;
621
622 layerZCoord = slaveGp->giveNaturalCoordinate(2);
623 layerYCoord = slaveGp->giveNaturalCoordinate(1);
624
625 answer.resize(3); // {Exx,GMzx,GMxy}
626
627 answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(5) * layerZCoord - masterGpStrain.at(6) * layerYCoord;
628 answer.at(2) = masterGpStrain.at(2);
629 answer.at(3) = masterGpStrain.at(3);
630}
631
632
633Interface *
635{
636 if ( interface == FiberedCrossSectionInterfaceType ) {
637 return static_cast< FiberedCrossSectionInterface * >( this );
638 }
639
640 return NULL;
641}
642
643
644#ifdef __OOFEG
645void
647{
648 GraphicObj *go;
649
650 if ( !gc.testElementGraphicActivity(this) ) {
651 return;
652 }
653
654 // if (!go) { // create new one
655 WCRec p[ 2 ]; /* poin */
656 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
657 EASValsSetColor(gc.getElementColor() );
658 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
659 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
660 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
661 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
662 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
663 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
664 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
665 go = CreateLine3D(p);
666 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
667 EGAttachObject(go, ( EObjectP ) this);
668 EMAddGraphicsToModel(ESIModel(), go);
669}
670
671
672void
674{
675 GraphicObj *go;
676
677 if ( !gc.testElementGraphicActivity(this) ) {
678 return;
679 }
680
681 double defScale = gc.getDefScale();
682 // if (!go) { // create new one
683 WCRec p[ 2 ]; /* poin */
684 char const *colors[] = {
685 "red", "green", "blue"
686 };
687
688 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
689 EASValsSetColor(gc.getDeformedElementColor() );
690 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
691 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
692 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
693 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
694
695 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
696 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
697 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
698 go = CreateLine3D(p);
699 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
700 EMAddGraphicsToModel(ESIModel(), go);
701
702 // draw centre triad
703 int i, succ;
704 double coeff = this->computeLength() / 3.;
705 //this->updateTempTriad (tStep);
706
707 p [ 0 ].x = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale) );
708 p [ 0 ].y = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale) );
709 p [ 0 ].z = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale) );
710
711 // draw t1
712 for ( i = 1; i <= 3; i++ ) {
713 p [ 1 ].x = p [ 0 ].x + coeff * tc.at(1, i);
714 p [ 1 ].y = p [ 0 ].y + coeff * tc.at(2, i);
715 p [ 1 ].z = p [ 0 ].z + coeff * tc.at(3, i);
716
717 EASValsSetColor(ColorGetPixelFromString(const_cast< char * >( colors [ i - 1 ] ), & succ) );
718
719 go = CreateLine3D(p);
720 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
721 EMAddGraphicsToModel(ESIModel(), go);
722 }
723}
724
725
726void
728{
729 WCRec p[ 2 ];
730 GraphicObj *go;
731 FloatArray v;
732 double defScale;
733
734 if ( !gc.testElementGraphicActivity(this) ) {
735 return;
736 }
737
738 EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
739 if ( gc.getInternalVarsDefGeoFlag() ) {
740 // use deformed geometry
741 defScale = gc.getDefScale();
742 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
743 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
744 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
745 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
746 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
747 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
748 } else {
749 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
750 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
751 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
752 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
753 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
754 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
755 }
756
757 GaussPoint *gp;
758 double s;
759 gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
760 if ( giveIPValue(v, gp, gc.giveIntVarType(), tStep) == 0 ) {
761 return;
762 }
763
764 s = v.at(1);
765 gc.updateFringeTableMinMax(& s, 1);
766
767 go = CreateLine3D(p);
768 EASValsSetColor(gc.getElementColor() );
769 EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
770 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
771 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
772 EMAddGraphicsToModel(ESIModel(), go);
773
774 EASValsSetMType(FILLED_CIRCLE_MARKER);
775 go = CreateMarkerWD3D(p, v.at(1) );
776 EGWithMaskChangeAttributes(LAYER_MASK | FILL_MASK | MTYPE_MASK, go);
777 EMAddGraphicsToModel(ESIModel(), go);
778}
779
780#endif
781} // 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 n)
Definition domain.h:398
Node * giveNode(int i) const
Definition element.h:629
virtual void updateYourself(TimeStep *tStep)
Definition element.C:824
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
void postInitialize() override
Performs post initialization steps.
Definition element.C:797
void saveContext(DataStream &stream, ContextMode mode) override
Definition element.C:923
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
virtual void initForNewStep()
Definition element.C:879
CrossSection * giveCrossSection()
Definition element.C:534
void restoreContext(DataStream &stream, ContextMode mode) override
Definition element.C:999
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
double computeNorm() const
Definition floatarray.C:861
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 beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void add(const FloatArray &src)
Definition floatarray.C:218
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 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
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
StateCounterType tempTcCounter
Time stamp of temporary centre triad.
Definition libeam3d2.h:72
bool computeGtoLRotationMatrix(FloatMatrix &answer) override
Definition libeam3d2.C:187
double giveCurrentLength(TimeStep *tStep)
Definition libeam3d2.C:548
double computeLength() override
Definition libeam3d2.C:284
void initForNewStep() override
Definition libeam3d2.C:584
void giveEdgeDofMapping(IntArray &answer, int iEdge) const override
Definition libeam3d2.C:328
void updateYourself(TimeStep *tStep) override
Definition libeam3d2.C:567
void computeRotMtrx(FloatMatrix &answer, FloatArray &psi)
Definition libeam3d2.C:477
void initializeFrom(InputRecord &ir, int prio) override
Definition libeam3d2.C:304
void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override
Definition libeam3d2.C:673
void drawScalar(oofegGraphicContext &gc, TimeStep *tStep) override
Definition libeam3d2.C:727
static ParamKey IPK_LIBeam3d2_refnode
Definition libeam3d2.h:73
double computeVolumeAround(GaussPoint *gp) override
Definition libeam3d2.C:222
int computeLoadGToLRotationMtrx(FloatMatrix &answer) override
Definition libeam3d2.C:359
void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) override
Definition libeam3d2.C:387
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition libeam3d2.C:74
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override
Definition libeam3d2.C:178
void updateTempTriad(TimeStep *tStep)
Definition libeam3d2.C:448
int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) override
Definition libeam3d2.C:239
void computeConstitutiveMatrix_dPdF_At(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition libeam3d2.C:263
void saveContext(DataStream &stream, ContextMode mode) override
Definition libeam3d2.C:594
double computeEdgeVolumeAround(GaussPoint *gp, int iEdge) override
Definition libeam3d2.C:347
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition libeam3d2.C:133
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition libeam3d2.C:257
void postInitialize() override
Performs post initialization steps.
Definition libeam3d2.C:313
Interface * giveInterface(InterfaceType it) override
Definition libeam3d2.C:634
int testElementExtension(ElementExtension ext) override
Definition libeam3d2.C:277
int giveLocalCoordinateSystem(FloatMatrix &answer) override
Definition libeam3d2.C:411
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition libeam3d2.C:232
void restoreContext(DataStream &stream, ContextMode mode) override
Definition libeam3d2.C:605
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition libeam3d2.C:270
void computeSMtrx(FloatMatrix &answer, FloatArray &vec)
Definition libeam3d2.C:507
void computeGaussPoints() override
Definition libeam3d2.C:120
void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep) override
Definition libeam3d2.C:646
FloatMatrix tc
Last equilibrium triad at the centre.
Definition libeam3d2.h:68
LIBeam3d2(int n, Domain *d)
Definition libeam3d2.C:63
FloatMatrix tempTc
Temporary triad at the centre.
Definition libeam3d2.h:70
int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) override
Definition libeam3d2.C:396
void FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *slaveGp, TimeStep *tStep) override
Definition libeam3d2.C:617
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override
Definition libeam3d2.C:526
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition libeam3d2.C:148
void initializeFrom(InputRecord &ir, int priority) override
NLStructuralElement(int n, Domain *d)
int nlGeometry
Flag indicating if geometrical nonlinearities apply.
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Definition node.C:234
virtual FloatMatrixF< 6, 6 > give3dBeamStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatArrayF< 6 > giveGeneralizedStress_Beam3d(const FloatArrayF< 6 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
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.
StateCounterType giveSolutionStateCounter()
Definition timestep.h:211
#define THROW_CIOERR(e)
#define OOFEM_ERROR(...)
Definition error.h:79
#define S(p)
Definition mdm.C:469
long ContextMode
Definition contextmode.h:43
@ CS_Area
Area.
@ FiberedCrossSectionInterfaceType
@ Element_EdgeLoadSupport
Element extension for edge loads.
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_VARPLOT_PATTERN_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