OOFEM 3.0
Loading...
Searching...
No Matches
latticebeam3dboundary.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 program is free software; you can redistribute it and/or modify
21 * it under the terms of the GNU General Public License as published by
22 * the Free Software Foundation; either version 2 of the License, or
23 * (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
28 * GNU General Public License for more details.
29 *
30 * You should have received a copy of the GNU General Public License
31 * along with this program; if not, write to the Free Software
32 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
33 */
34
35#include "domain.h"
39#include "node.h"
40#include "material.h"
41#include "gausspoint.h"
43#include "floatmatrix.h"
44#include "intarray.h"
45#include "floatarray.h"
46#include "mathfem.h"
48#include "classfactory.h"
50#include "contextioerr.h"
51#include "datastream.h"
53#include "latticebeam3d.h"
54#include "parametermanager.h"
55#include "paramkey.h"
56
57#ifdef __OOFEG
58 #include "oofeggraphiccontext.h"
59 #include "dof.h"
60
61#endif
62
63namespace oofem {
66
67LatticeBeam3dBoundary :: LatticeBeam3dBoundary(int n, Domain *aDomain) : LatticeBeam3d(n, aDomain), location(2)
68{
70 geometryFlag = 0;
71}
72
73LatticeBeam3dBoundary :: ~LatticeBeam3dBoundary()
74{}
75
76
77void
78LatticeBeam3dBoundary :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
79 TimeStep *tStep)
80// Computes numerically the stiffness matrix of the receiver.
81{
82 // double dV;
83 FloatMatrix d, bi, bj, dbj, dij, bjt;
84 FloatMatrix t(12, 18), tt;
85 FloatMatrix answerTemp(12, 12), ttk(18, 12);
86 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
87 answerTemp.zero();
88 t.zero();
89
90 //Call the parent and compute stiffness
91 LatticeBeam3d :: computeStiffnessMatrix(answerTemp, rMode, tStep);
92
94 answer.zero();
95
96 for ( int m = 1; m <= 12; m++ ) {
97 for ( int k = 1; k <= 12; k++ ) {
98 answer.at(m, k) = answerTemp.at(m, k);
99 }
100 }
101
102 if ( matStiffSymmFlag ) {
103 answer.symmetrized();
104 }
105
106 // Rotate to global system add the projections and rotate back to local system.
107
108 FloatMatrix R;
109 if ( this->giveRotationMatrix(R) ) {
110 answer.rotatedWith(R);
111 }
112
113 for ( int m = 1; m <= 12; m++ ) {
114 for ( int k = 1; k <= 12; k++ ) {
115 answerTemp.at(m, k) = answer.at(m, k);
116 }
117 }
118
119 IntArray projectionComponentNodeOne(3);
120 projectionComponentNodeOne.zero();
121 if ( location.at(1) != 0 ) {
122 giveSwitches(projectionComponentNodeOne, location.at(1) );
123 }
124
125 IntArray projectionComponentNodeTwo(3);
126 projectionComponentNodeTwo.zero();
127 if ( location.at(2) != 0 ) {
128 giveSwitches(projectionComponentNodeTwo, location.at(2) );
129 }
130
131 for ( int k = 1; k <= 12; k++ ) {
132 t.at(k, k) = 1.;
133 }
134 t.at(1, 13) = projectionComponentNodeOne.at(1);
135 t.at(2, 14) = projectionComponentNodeOne.at(2);
136 t.at(3, 15) = projectionComponentNodeOne.at(3);
137
138
139 t.at(2, 16) = projectionComponentNodeOne.at(3);
140 t.at(1, 17) = projectionComponentNodeOne.at(3);
141 t.at(1, 18) = projectionComponentNodeOne.at(2);
142
143
144 t.at(7, 13) = projectionComponentNodeTwo.at(1);
145 t.at(8, 14) = projectionComponentNodeTwo.at(2);
146 t.at(9, 15) = projectionComponentNodeTwo.at(3);
147
148 t.at(8, 16) = projectionComponentNodeTwo.at(3);
149 t.at(7, 17) = projectionComponentNodeTwo.at(3);
150 t.at(7, 18) = projectionComponentNodeTwo.at(2);
151
152
153 tt.beTranspositionOf(t);
154
155 ttk.beProductOf(tt, answerTemp);
156 answer.beProductOf(ttk, t);
157
158 //Rotate back to local system
159 FloatMatrix Rtranspose;
160 Rtranspose.beTranspositionOf(R);
161 answer.rotatedWith(Rtranspose);
162
163 return;
164}
165
166
167
168double
169LatticeBeam3dBoundary :: computeVolumeAround(GaussPoint *aGaussPoint)
170{
171 if ( geometryFlag == 0 ) {
173 }
174
175 return this->area * this->length;
176}
177
178void
179LatticeBeam3dBoundary :: giveVTKCoordinates(int nodeNumber, FloatArray &coords) {
180 coords.resize(3);
181 coords.zero();
182 Node *node;
183
184 FloatArray specimenDimension(3);
185 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
186 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
187 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
188
189 IntArray projectionComponent(3);
190 projectionComponent.zero();
191
192 if ( nodeNumber == 1 ) {
193 node = this->giveNode(1);
194 if ( location.at(1) != 0 ) {
195 giveSwitches(projectionComponent, location.at(1) );
196 }
197 } else if ( nodeNumber == 2 ) {
198 node = this->giveNode(2);
199 if ( location.at(2) != 0 ) {
200 giveSwitches(projectionComponent, location.at(2) );
201 }
202 } else {
203 OOFEM_ERROR("wrong element used in the vtk module");
204 }
205
206 for ( int i = 0; i < 3; i++ ) {
207 coords.at(i + 1) = node->giveCoordinate(i + 1) + projectionComponent.at(i + 1) * specimenDimension.at(i + 1);
208 }
209 return;
210}
211
212bool
213LatticeBeam3dBoundary :: computeGtoLRotationMatrix(FloatMatrix &answer)
214{
215 FloatMatrix lcs;
216 int i, j;
217
218 answer.resize(18, 18);
219 answer.zero();
220
221 this->giveLocalCoordinateSystem(lcs);
222 for ( i = 1; i <= 3; i++ ) {
223 for ( j = 1; j <= 3; j++ ) {
224 answer.at(i, j) = lcs.at(i, j);
225 answer.at(i + 3, j + 3) = lcs.at(i, j);
226 answer.at(i + 6, j + 6) = lcs.at(i, j);
227 answer.at(i + 9, j + 9) = lcs.at(i, j);
228 }
229 }
230 answer.at(13, 13) = 1.;
231 answer.at(14, 14) = 1.;
232 answer.at(15, 15) = 1.;
233 answer.at(16, 16) = 1.;
234 answer.at(17, 17) = 1.;
235 answer.at(18, 18) = 1.;
236
237 return 1;
238}
239
240int
241LatticeBeam3dBoundary :: giveLocalCoordinateSystem(FloatMatrix &answer)
242{
243 if ( geometryFlag == 0 ) {
245 }
246
247 answer = this->localCoordinateSystem;
248
249 return 1;
250}
251
252
253
254void
255LatticeBeam3dBoundary :: giveDofManDofIDMask(int inode, IntArray &answer) const
256{
257 answer = {
258 D_u, D_v, D_w, R_u, R_v, R_w
259 };
260}
261
262void
263LatticeBeam3dBoundary :: initializeFrom(InputRecord &ir, int priority)
264{
265 ParameterManager &ppm = this->giveDomain()->elementPPM;
266 LatticeBeam3d :: initializeFrom(ir, priority);
268}
269
270void
271LatticeBeam3dBoundary :: postInitialize()
272{
273 ParameterManager &ppm = this->giveDomain()->elementPPM;
274 LatticeBeam3d :: postInitialize();
276}
277
278void
279LatticeBeam3dBoundary :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
280{
281 FloatMatrix stiffness;
282 FloatArray u;
283
284 this->computeVectorOf(VM_Total, tStep, u);
285
286 // subtract initial displacements, if defined
287 if ( initialDisplacements ) {
289 }
290
291 // //Compute projection vector
292 // IntArray projectionComponentNodeOne(3);
293 // projectionComponentNodeOne.zero();
294 // if ( location.at(1) != 0 ) {
295 // giveSwitches( projectionComponentNodeOne, location.at(1) );
296 // }
297
298 // IntArray projectionComponentNodeTwo(3);
299 // projectionComponentNodeTwo.zero();
300
301 // if ( location.at(2) != 0 ) {
302 // giveSwitches( projectionComponentNodeTwo, location.at(2) );
303 // }
304
305 // //Rotate displacements into global coordinate system
306 // FloatMatrix rotationMatrix;
307 // if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
308 // u.rotatedWith(rotationMatrix, 't');
309 // }
310
311 // //General expressions for the corrected displacements. Rotations at nodes are not influenced. Only translations
312 // u.at(1) = u.at(1) + projectionComponentNodeOne.at(1) * u.at(13) + projectionComponentNodeOne.at(2) * u.at(18) + projectionComponentNodeOne.at(3) * u.at(17);
313 // u.at(2) = u.at(2) + projectionComponentNodeOne.at(2) * u.at(14) + projectionComponentNodeOne.at(3) * u.at(16);
314 // u.at(3) = u.at(3) + projectionComponentNodeOne.at(3) * u.at(15);
315
316 // u.at(7) = u.at(7) + projectionComponentNodeTwo.at(1) * u.at(13) + projectionComponentNodeTwo.at(2) * u.at(18) + projectionComponentNodeTwo.at(3) * u.at(17);
317 // u.at(8) = u.at(8) + projectionComponentNodeTwo.at(2) * u.at(14) + projectionComponentNodeTwo.at(3) * u.at(16);
318 // u.at(9) = u.at(9) + projectionComponentNodeTwo.at(3) * u.at(15);
319
320
321 // if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
322 // u.rotatedWith(rotationMatrix, 'n');
323 // }
324
325 this->computeStiffnessMatrix(stiffness, ElasticStiffness, tStep);
326
327 // zero answer will resize accordingly when adding first contribution
328 answer.clear();
329 answer.beProductOf(stiffness, u);
330
331 // if inactive update state, but no contribution to global system
332 if ( !this->isActivated(tStep) ) {
333 answer.zero();
334 return;
335 }
336 return;
337}
338
339
340void
341LatticeBeam3dBoundary :: giveSwitches(IntArray &answer, int location) {
342 int counter = 1;
343 for ( int x = -1; x < 2; x++ ) {
344 for ( int y = -1; y < 2; y++ ) {
345 for ( int z = -1; z < 2; z++ ) {
346 if ( !( z == 0 && y == 0 && x == 0 ) ) {
347 if ( counter == location ) {
348 answer(0) = x;
349 answer(1) = y;
350 answer(2) = z;
351 }
352 counter++;
353 }
354 }
355 }
356 }
357 return;
358}
359
360
361void
362LatticeBeam3dBoundary :: computeGeometryProperties()
363{
364 Node *nodeA, *nodeB;
365
366 FloatArray coordsA(3);
367 FloatArray coordsB(3);
368
369 nodeA = this->giveNode(1);
370 nodeB = this->giveNode(2);
371
372 FloatArray specimenDimension(3);
373 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
374 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
375 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
376
377 IntArray projectionComponentNodeOne(3);
378 projectionComponentNodeOne.zero();
379 if ( location.at(1) != 0 ) {
380 giveSwitches(projectionComponentNodeOne, location.at(1) );
381 }
382
383 IntArray projectionComponentNodeTwo(3);
384 projectionComponentNodeTwo.zero();
385 if ( location.at(2) != 0 ) {
386 giveSwitches(projectionComponentNodeTwo, location.at(2) );
387 }
388
389 for ( int i = 0; i < 3; i++ ) {
390 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1) + projectionComponentNodeOne.at(i + 1) * specimenDimension.at(i + 1);
391 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1) + projectionComponentNodeTwo.at(i + 1) * specimenDimension.at(i + 1);
392 }
393
394 //Calculate normal vector
395 FloatArray s(3), t(3);
396 this->midPoint.resize(3);
397
398 //Calculate normal vector
399 this->normal.resize(3);
400 for ( int i = 0; i < 3; i++ ) {
401 this->normal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
402 }
403
404 this->length = sqrt(pow(this->normal.at(1), 2.) + pow(this->normal.at(2), 2.) + pow(this->normal.at(3), 2.) );
405
406 for ( int i = 0; i < 3; i++ ) {
407 this->normal.at(i + 1) /= length;
408 }
409
410 // Compute midpoint
411 this->midPoint.resize(3);
412 for ( int i = 0; i < 3; i++ ) {
413 this->midPoint.at(i + 1) = 0.5 * ( coordsB.at(i + 1) + coordsA.at(i + 1) );
414 }
415
416 this->globalCentroid = this->midPoint;
417
419
420 //Set geometry flag to 1 so that this is done only once
421 this->geometryFlag = 1;
422
423 return;
424}
425
426
427void
428LatticeBeam3dBoundary :: saveContext(DataStream &stream, ContextMode mode)
429{
430 LatticeBeam3d :: saveContext(stream, mode);
431
433
434 if ( ( mode & CM_Definition ) ) {
435 if ( ( iores = location.storeYourself(stream) ) != CIO_OK ) {
436 THROW_CIOERR(iores);
437 }
438 }
439}
440
441void
442LatticeBeam3dBoundary :: restoreContext(DataStream &stream, ContextMode mode)
443{
444 LatticeBeam3d :: restoreContext(stream, mode);
445
447
448 if ( mode & CM_Definition ) {
449 if ( ( iores = this->location.restoreYourself(stream) ) != CIO_OK ) {
450 THROW_CIOERR(iores);
451 }
452 }
453}
454
455
456int
457LatticeBeam3dBoundary :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
458{
459 if ( geometryFlag == 0 ) {
461 }
462
463 answer.resize(3);
464 answer = this->globalCentroid;
465
466 return 1;
467}
468
469
470
471#ifdef __OOFEG
472
473void
474LatticeBeam3dBoundary :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
475{
476 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
477
478 if ( mode == OGC_rawGeometry ) {
479 this->drawRawGeometry(gc, tStep);
480 } else if ( mode == OGC_deformedGeometry ) {
481 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
482 } else if ( mode == OGC_elemSpecial ) {
483 this->drawSpecial(gc, tStep);
484 } else {
485 OOFEM_ERROR("drawYourself : unsupported mode");
486 }
487}
488
489
490void LatticeBeam3dBoundary :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
491{
492 GraphicObj *go;
493
494 WCRec p [ 2 ]; /* points */
495 if ( !gc.testElementGraphicActivity(this) ) {
496 return;
497 }
498
499 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
500 EASValsSetColor(gc.getActiveCrackColor() );
501 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
502
503
504 FloatArray specimenDimension(3);
505 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
506 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
507 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
508
509 //Compute projection vector
510 IntArray projectionComponentNodeOne(3);
511 projectionComponentNodeOne.zero();
512 if ( location.at(1) != 0 ) {
513 giveSwitches(projectionComponentNodeOne, location.at(1) );
514 }
515
516 IntArray projectionComponentNodeTwo(3);
517 projectionComponentNodeTwo.zero();
518
519 if ( location.at(2) != 0 ) {
520 giveSwitches(projectionComponentNodeTwo, location.at(2) );
521 }
522
523 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
524 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
525 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
526 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
527 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
528 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
529
530 go = CreateLine3D(p);
531 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
532 EGAttachObject(go, ( EObjectP ) this);
533 EMAddGraphicsToModel(ESIModel(), go);
534}
535
536
537void LatticeBeam3dBoundary :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
538{
539 GraphicObj *go;
540
541 if ( !gc.testElementGraphicActivity(this) ) {
542 return;
543 }
544
545 double defScale = gc.getDefScale();
546
547 WCRec p [ 2 ]; /* points */
548
549 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
550 EASValsSetColor(gc.getDeformedElementColor() );
551 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
552
553 FloatArray specimenDimension(3);
554 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
555 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
556 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
557
558 FloatArray dispOne(6);
559 dispOne.at(1) = this->giveNode(1)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
560 dispOne.at(2) = this->giveNode(1)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
561 dispOne.at(3) = this->giveNode(1)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
562 dispOne.at(4) = this->giveNode(1)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
563 dispOne.at(5) = this->giveNode(1)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
564 dispOne.at(6) = this->giveNode(1)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
565
566 FloatArray dispTwo(6);
567 dispTwo.at(1) = this->giveNode(2)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
568 dispTwo.at(2) = this->giveNode(2)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
569 dispTwo.at(3) = this->giveNode(2)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
570 dispTwo.at(4) = this->giveNode(2)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
571 dispTwo.at(5) = this->giveNode(2)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
572 dispTwo.at(6) = this->giveNode(2)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
573
574 FloatArray dispThree(6);
575 dispThree.at(1) = this->giveNode(3)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
576 dispThree.at(2) = this->giveNode(3)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
577 dispThree.at(3) = this->giveNode(3)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
578 dispThree.at(4) = this->giveNode(3)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
579 dispThree.at(5) = this->giveNode(3)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
580 dispThree.at(6) = this->giveNode(3)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
581
582 IntArray projectionComponentNodeOne(3);
583 projectionComponentNodeOne.zero();
584 if ( location.at(1) != 0 ) {
585 giveSwitches(projectionComponentNodeOne, location.at(1) );
586 }
587
588 IntArray projectionComponentNodeTwo(3);
589 projectionComponentNodeTwo.zero();
590 if ( location.at(2) != 0 ) {
591 giveSwitches(projectionComponentNodeTwo, location.at(2) );
592 }
593
594
595 //Modify dispOne and dispTwo
596 dispOne.at(1) = dispOne.at(1) + projectionComponentNodeOne.at(1) * dispThree.at(1) + projectionComponentNodeOne.at(2) * dispThree.at(4) + projectionComponentNodeOne.at(3) * dispThree.at(5);
597 dispOne.at(2) = dispOne.at(2) + projectionComponentNodeOne.at(2) * dispThree.at(2) + projectionComponentNodeOne.at(3) * dispThree.at(6);
598 dispOne.at(3) = dispOne.at(3) + projectionComponentNodeOne.at(3) * dispThree.at(3);
599
600
601 dispTwo.at(1) = dispTwo.at(1) + projectionComponentNodeTwo.at(1) * dispThree.at(1) + projectionComponentNodeTwo.at(2) * dispThree.at(4) + projectionComponentNodeTwo.at(3) * dispThree.at(5);
602 dispTwo.at(2) = dispTwo.at(2) + projectionComponentNodeTwo.at(2) * dispThree.at(2) + projectionComponentNodeTwo.at(3) * dispThree.at(6);
603 dispTwo.at(3) = dispTwo.at(3) + projectionComponentNodeTwo.at(3) * dispThree.at(3);
604
605 double x1, y1, z1, x2, y2, z2;
606 x1 = this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
607 y1 = this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
608 z1 = this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
609
610 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
611 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
612 z2 = this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
613
614 p [ 0 ].x = ( FPNum ) x1 + defScale * dispOne.at(1);
615 p [ 0 ].y = ( FPNum ) y1 + defScale * dispOne.at(2);
616 p [ 0 ].z = ( FPNum ) z1 + defScale * dispOne.at(3);
617
618 p [ 1 ].x = ( FPNum ) x2 + defScale * dispTwo.at(1);
619 p [ 1 ].y = ( FPNum ) y2 + defScale * dispTwo.at(2);
620 p [ 1 ].z = ( FPNum ) z2 + defScale * dispTwo.at(3);
621
622 go = CreateLine3D(p);
623 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
624 EMAddGraphicsToModel(ESIModel(), go);
625}
626
627#endif
628} // 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
virtual bool giveRotationMatrix(FloatMatrix &answer)
Definition element.C:298
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
CrossSection * giveCrossSection()
Definition element.C:534
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
Definition element.h:1078
Domain * giveDomain() const
Definition femcmpnn.h:97
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 zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void subtract(const FloatArray &src)
Definition floatarray.C:320
void rotatedWith(const FloatMatrix &r, char mode='n')
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
void zero()
Sets all component to zero.
Definition intarray.C:52
int & at(std::size_t i)
Definition intarray.h:104
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override
virtual void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
void giveSwitches(IntArray &answer, int location)
virtual void computeGeometryProperties() override
static ParamKey IPK_LatticeBeam3dBoundary_location
int giveLocalCoordinateSystem(FloatMatrix &answer) override
virtual void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
LatticeBeam3d(int n, Domain *)
FloatArray globalCentroid
FloatMatrix localCoordinateSystem
virtual void computeCrossSectionProperties()
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
#define THROW_CIOERR(e)
#define CM_Definition
Definition contextmode.h:47
#define OOFEM_ERROR(...)
Definition error.h:79
long ContextMode
Definition contextmode.h:43
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