OOFEM 3.0
Loading...
Searching...
No Matches
lattice3dboundarytruss.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"
38#include "node.h"
39#include "material.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "intarray.h"
44#include "floatarray.h"
45#include "mathfem.h"
47#include "classfactory.h"
49#include "contextioerr.h"
50#include "datastream.h"
51#include "crosssection.h"
52#include "dof.h"
53#include "parametermanager.h"
54#include "paramkey.h"
55
56#ifdef __OOFEG
57 #include "oofeggraphiccontext.h"
58#endif
59
60namespace oofem {
63
64Lattice3dBoundaryTruss :: Lattice3dBoundaryTruss(int n, Domain *aDomain) : Lattice3dBoundary(n, aDomain), location (2)
65{
67 geometryFlag = 0;
68}
69
70Lattice3dBoundaryTruss :: ~Lattice3dBoundaryTruss()
71{}
72
73
74void
75Lattice3dBoundaryTruss :: computeBmatrixAt(GaussPoint *aGaussPoint, FloatMatrix &answer, int li, int ui)
76// Returns the strain matrix of the receiver.
77{
78 if ( geometryFlag == 0 ) {
80 }
81
82 //Assemble Bmatrix (used to compute strains and rotations}
83 answer.resize(6, 12);
84 answer.zero();
85
86 //Normal displacement jump in x-direction
87 //First node
88 answer.at(1, 1) = -1.;
89 answer.at(1, 2) = 0.;
90 answer.at(1, 3) = 0.;
91 answer.at(1, 4) = 0.;
92 answer.at(1, 5) = -this->eccT;
93 answer.at(1, 6) = this->eccS;
94 //Second node
95 answer.at(1, 7) = 1.;
96 answer.at(1, 8) = 0.;
97 answer.at(1, 9) = 0.;
98 answer.at(1, 10) = 0.;
99 answer.at(1, 11) = this->eccT;
100 answer.at(1, 12) = -this->eccS;
101
102 //Shear displacement jump in y-plane
103 //first node
104 answer.at(2, 1) = 0.;
105 answer.at(2, 2) = -1.;
106 answer.at(2, 3) = 0.;
107 answer.at(2, 4) = this->eccT;
108 answer.at(2, 5) = 0;
109 answer.at(2, 6) = -this->length / 2.;
110 //Second node
111 answer.at(2, 7) = 0.;
112 answer.at(2, 8) = 1.;
113 answer.at(2, 9) = 0.;
114 answer.at(2, 10) = -this->eccT;
115 answer.at(2, 11) = 0;
116 answer.at(2, 12) = -this->length / 2.;
117
118 //Shear displacement jump in z-plane
119 //first node
120 answer.at(3, 1) = 0.;
121 answer.at(3, 2) = 0.;
122 answer.at(3, 3) = -1.;
123 answer.at(3, 4) = -this->eccS;
124 answer.at(3, 5) = this->length / 2.;
125 answer.at(3, 6) = 0.;
126 //Second node
127 answer.at(3, 7) = 0.;
128 answer.at(3, 8) = 0.;
129 answer.at(3, 9) = 1.;
130 answer.at(3, 10) = this->eccS;
131 answer.at(3, 11) = this->length / 2.;
132 answer.at(3, 12) = 0.;
133
134 //Rotation around x-axis
135 //First node
136 answer.at(4, 1) = 0.;
137 answer.at(4, 2) = 0;
138 answer.at(4, 3) = 0.;
139 answer.at(4, 4) = -sqrt(Ip / this->area);
140 answer.at(4, 5) = 0.;
141 answer.at(4, 6) = 0.;
142 //Second node
143 answer.at(4, 7) = 0.;
144 answer.at(4, 8) = 0.;
145 answer.at(4, 9) = 0.;
146 answer.at(4, 10) = sqrt(Ip / this->area);
147 answer.at(4, 11) = 0.;
148 answer.at(4, 12) = 0.;
149
150 //Rotation around y-axis
151 //First node
152 answer.at(5, 1) = 0.;
153 answer.at(5, 2) = 0.;
154 answer.at(5, 3) = 0.;
155 answer.at(5, 4) = 0.;
156 answer.at(5, 5) = -sqrt(I1 / this->area);
157 answer.at(5, 6) = 0.;
158 //Second node
159 answer.at(5, 7) = 0.;
160 answer.at(5, 8) = 0.;
161 answer.at(5, 9) = 0.;
162 answer.at(5, 10) = 0.;
163 answer.at(5, 11) = sqrt(I1 / this->area);
164 answer.at(5, 12) = 0.;
165
166 //Rotation around z-axis
167 //First node
168 answer.at(6, 1) = 0.;
169 answer.at(6, 2) = 0.;
170 answer.at(6, 3) = 0.;
171 answer.at(6, 4) = 0.;
172 answer.at(6, 5) = 0.;
173 answer.at(6, 6) = -sqrt(I2 / this->area);
174 //Second node
175 answer.at(6, 7) = 0.;
176 answer.at(6, 8) = 0.;
177 answer.at(6, 9) = 0.;
178 answer.at(6, 10) = 0.;
179 answer.at(6, 11) = 0.;
180 answer.at(6, 12) = sqrt(I2 / this->area);
181
182 answer.times(1. / this->length);
183
184 return;
185}
186
187void
188Lattice3dBoundaryTruss :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
189 TimeStep *tStep)
190// Computes numerically the stiffness matrix of the receiver.
191{
192 // double dV;
193 FloatMatrix d, bi, bj, dbj, dij, bjt;
194 FloatMatrix t(12, 13), tt;
195 FloatMatrix answerTemp(12, 12), answerHelp, ttk(13, 12);
196 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
197 answerTemp.zero();
198 answerHelp.zero();
199 t.zero();
200
201 if ( geometryFlag == 0 ) {
203 }
204
205 double volume = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
206
207 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), bj);
208 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
209 for ( int i = 1; i <= 6; i++ ) {
210 d.at(i, i) *= volume;
211 }
212
213 dbj.beProductOf(d, bj);
214 bjt.beTranspositionOf(bj);
215 answerTemp.beProductOf(bjt, dbj);
216
218 answer.zero();
219
221 answerHelp.zero();
222
223 for ( int m = 1; m <= 12; m++ ) {
224 for ( int k = 1; k <= 12; k++ ) {
225 answer.at(m, k) = answerTemp.at(m, k);
226 }
227 }
228
229 if ( matStiffSymmFlag ) {
230 answer.symmetrized();
231 }
232
233 // Rotate to global system add the projections and rotate back to local system.
234
235 FloatMatrix R;
236 if ( this->giveRotationMatrix(R) ) {
237 answer.rotatedWith(R);
238 }
239
240 for ( int m = 1; m <= 12; m++ ) {
241 for ( int k = 1; k <= 12; k++ ) {
242 answerTemp.at(m, k) = answer.at(m, k);
243 }
244 }
245
246 IntArray projectionComponentNodeOne(3);
247 projectionComponentNodeOne.zero();
248 if ( location.at(1) != 0 ) {
249 giveSwitches(projectionComponentNodeOne, location.at(1) );
250 }
251
252 IntArray projectionComponentNodeTwo(3);
253 projectionComponentNodeTwo.zero();
254 if ( location.at(2) != 0 ) {
255 giveSwitches(projectionComponentNodeTwo, location.at(2) );
256 }
257
258 for ( int k = 1; k <= 12; k++ ) {
259 t.at(k, k) = 1.;
260 }
261 t.at(1, 13) = projectionComponentNodeOne.at(1);
262 t.at(7, 13) = projectionComponentNodeTwo.at(1);
263
264 tt.beTranspositionOf(t);
265
266 ttk.beProductOf(tt, answerTemp);
267 answer.beProductOf(ttk, t);
268
269 //Rotate back to local system
270 FloatMatrix Rtranspose;
271 Rtranspose.beTranspositionOf(R);
272 answer.rotatedWith(Rtranspose);
273
274 return;
275}
276
277
278double
279Lattice3dBoundaryTruss :: computeVolumeAround(GaussPoint *aGaussPoint)
280{
281 if ( geometryFlag == 0 ) {
283 }
284
285 return this->area * this->length;
286}
287
288void
289Lattice3dBoundaryTruss :: recalculateCoordinates(int nodeNumber, FloatArray &coords) {
290 coords.resize(3);
291 coords.zero();
292 Node *node;
293
294 FloatArray specimenDimension(3);
295 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
296 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
297 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
298
299 IntArray projectionComponent(3);
300 projectionComponent.zero();
301
302 if ( nodeNumber == 1 ) {
303 node = this->giveNode(1);
304 if ( location.at(1) != 0 ) {
305 giveSwitches(projectionComponent, location.at(1) );
306 }
307 } else if ( nodeNumber == 2 ) {
308 node = this->giveNode(2);
309 if ( location.at(2) != 0 ) {
310 giveSwitches(projectionComponent, location.at(2) );
311 }
312 } else {
313 OOFEM_ERROR("wrong element used in the vtk module");
314 }
315
316 for ( int i = 0; i < 3; i++ ) {
317 coords.at(i + 1) = node->giveCoordinate(i + 1) + projectionComponent.at(i + 1) * specimenDimension.at(i + 1);
318 }
319 return;
320}
321
322void
323Lattice3dBoundaryTruss :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN)
324// Computes the vector containing the strains at the Gauss point gp of
325// the receiver, at time step stepN. The nature of these strains depends
326// on the element's type.
327{
328 FloatMatrix b;
329 FloatArray u;
330
331 //Compute strain vector
332 //Get the 13 components of the displacement vector of this element
333 this->computeVectorOf(VM_Total, stepN, u);
334 this->computeBmatrixAt(gp, b);
335
336 // subtract initial displacements, if defined
337 if ( initialDisplacements ) {
339 }
340
341 //Compute projection vector
342 IntArray projectionComponentNodeOne(3);
343 projectionComponentNodeOne.zero();
344 if ( location.at(1) != 0 ) {
345 giveSwitches(projectionComponentNodeOne, location.at(1) );
346 }
347
348 IntArray projectionComponentNodeTwo(3);
349 projectionComponentNodeTwo.zero();
350
351 if ( location.at(2) != 0 ) {
352 giveSwitches(projectionComponentNodeTwo, location.at(2) );
353 }
354
355 FloatMatrix rotationMatrix;
356 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
357 u.rotatedWith(rotationMatrix, 't');
358 }
359
360 //General expressions for the corrected displacements. Rotations at nodes are not influenced. Only translations
361 u.at(1) = u.at(1) + projectionComponentNodeOne.at(1) * u.at(13);
362
363 u.at(7) = u.at(7) + projectionComponentNodeTwo.at(1) * u.at(13);
364
365 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
366 u.rotatedWith(rotationMatrix, 'n');
367 }
368 //define a temp displacement vector
369 FloatArray uTemp(12);
370
371 for ( int i = 1; i <= 12; i++ ) {
372 uTemp.at(i) = u.at(i);
373 }
374
375 answer.beProductOf(b, uTemp);
376}
377
378bool
379Lattice3dBoundaryTruss :: computeGtoLRotationMatrix(FloatMatrix &answer)
380{
381 FloatMatrix lcs;
382 int i, j;
383
384 answer.resize(13, 13);
385 answer.zero();
386
387 this->giveLocalCoordinateSystem(lcs);
388 for ( i = 1; i <= 3; i++ ) {
389 for ( j = 1; j <= 3; j++ ) {
390 answer.at(i, j) = lcs.at(i, j);
391 answer.at(i + 3, j + 3) = lcs.at(i, j);
392 answer.at(i + 6, j + 6) = lcs.at(i, j);
393 answer.at(i + 9, j + 9) = lcs.at(i, j);
394 }
395 }
396 answer.at(13, 13) = 1.;
397
398 return 1;
399}
400
401int
402Lattice3dBoundaryTruss :: giveLocalCoordinateSystem(FloatMatrix &answer)
403{
404 if ( geometryFlag == 0 ) {
406 }
407
408 answer = this->localCoordinateSystem;
409
410 return 1;
411}
412
413
414
415void
416Lattice3dBoundaryTruss :: giveDofManDofIDMask(int inode, IntArray &answer) const
417{
418 if ( inode == 3 ) {
419 answer = { E_xx };
420 } else {
421 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
422 }
423}
424
425void
426Lattice3dBoundaryTruss :: initializeFrom(InputRecord &ir, int priority)
427{
428 ParameterManager &ppm = this->giveDomain()->elementPPM;
429 Lattice3d :: initializeFrom(ir, priority);
431}
432
433void
434Lattice3dBoundaryTruss :: postInitialize()
435{
436 ParameterManager &ppm = this->giveDomain()->elementPPM;
437 // first call parent
438 Lattice3d :: postInitialize();
440}
441
442
443void
444Lattice3dBoundaryTruss :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
445{
446 Material *mat = this->giveMaterial();
447
448 FloatMatrix b, bt, A, R, GNT;
449 FloatArray bs, TotalStressVector, u, strain;
450 double dV;
451
452 this->computeVectorOf(VM_Total, tStep, u);
453
454 // subtract initial displacements, if defined
455 if ( initialDisplacements ) {
457 }
458
459 answer.resize(13);
460 answer.zero();
461
462 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), b);
463
464 bt.beTranspositionOf(b);
465 if ( useUpdatedGpRecord == 1 ) {
466 TotalStressVector = ( ( LatticeMaterialStatus * ) mat->giveStatus(integrationRulesArray [ 0 ]->getIntegrationPoint(0) ) )
467 ->giveLatticeStress();
468 } else
469 if ( !this->isActivated(tStep) ) {
470 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(integrationRulesArray [ 0 ]->getIntegrationPoint(0)->giveMaterialMode() ) );
471 strain.zero();
472 }
473 this->computeStrainVector(strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
474
475 this->computeStressVector(TotalStressVector, strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
476
477 dV = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
478 bs.beProductOf(bt, TotalStressVector);
479 bs.times(dV);
480
481 for ( int m = 1; m <= 12; m++ ) {
482 answer.at(m) = bs.at(m);
483 }
484
485 if ( this->giveRotationMatrix(R) ) {
486 answer.rotatedWith(R, 't');
487 }
488
489 IntArray projectionComponentNodeOne(3);
490 projectionComponentNodeOne.zero();
491 if ( location.at(1) != 0 ) {
492 giveSwitches(projectionComponentNodeOne, location.at(1) );
493 }
494
495 IntArray projectionComponentNodeTwo(3);
496 projectionComponentNodeTwo.zero();
497
498 if ( location.at(2) != 0 ) {
499 giveSwitches(projectionComponentNodeTwo, location.at(2) );
500 }
501
502 //Normal stresses
503 answer.at(13) = projectionComponentNodeOne.at(1) * answer.at(1) + projectionComponentNodeTwo.at(1) * answer.at(7);
504
505 //Rotate to local system
506 if ( this->giveRotationMatrix(R) ) {
507 answer.rotatedWith(R, 'n');
508 }
509
510 return;
511}
512
513
514void
515Lattice3dBoundaryTruss :: giveSwitches(IntArray &answer, int location) {
516 int counter = 1;
517 for ( int x = -1; x < 2; x++ ) {
518 for ( int y = -1; y < 2; y++ ) {
519 for ( int z = -1; z < 2; z++ ) {
520 if ( !( z == 0 && y == 0 && x == 0 ) ) {
521 if ( counter == location ) {
522 answer(0) = x;
523 answer(1) = y;
524 answer(2) = z;
525 }
526 counter++;
527 }
528 }
529 }
530 }
531 return;
532}
533
534
535void
536Lattice3dBoundaryTruss :: computeGeometryProperties()
537{
538 Node *nodeA, *nodeB;
539
540 FloatArray coordsA(3);
541 FloatArray coordsB(3);
542
543 nodeA = this->giveNode(1);
544 nodeB = this->giveNode(2);
545
546 FloatArray specimenDimension(3);
547 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
548 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
549 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
550
551 IntArray projectionComponentNodeOne(3);
552 projectionComponentNodeOne.zero();
553 if ( location.at(1) != 0 ) {
554 giveSwitches(projectionComponentNodeOne, location.at(1) );
555 }
556
557 IntArray projectionComponentNodeTwo(3);
558 projectionComponentNodeTwo.zero();
559 if ( location.at(2) != 0 ) {
560 giveSwitches(projectionComponentNodeTwo, location.at(2) );
561 }
562
563 for ( int i = 0; i < 3; i++ ) {
564 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1) + projectionComponentNodeOne.at(i + 1) * specimenDimension.at(i + 1);
565 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1) + projectionComponentNodeTwo.at(i + 1) * specimenDimension.at(i + 1);
566 }
567
568 //Calculate normal vector
569 FloatArray s(3), t(3);
570 this->midPoint.resize(3);
571
572 //Calculate normal vector
573 this->normal.resize(3);
574 for ( int i = 0; i < 3; i++ ) {
575 this->normal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
576 }
577
578 this->length = sqrt(pow(this->normal.at(1), 2.) + pow(this->normal.at(2), 2.) + pow(this->normal.at(3), 2.) );
579
580 for ( int i = 0; i < 3; i++ ) {
581 this->normal.at(i + 1) /= length;
582 }
583
584 // Compute midpoint
585 this->midPoint.resize(3);
586 for ( int i = 0; i < 3; i++ ) {
587 this->midPoint.at(i + 1) = 0.5 * ( coordsB.at(i + 1) + coordsA.at(i + 1) );
588 }
589
591
592 //Set geometry flag to 1 so that this is done only once
593 this->geometryFlag = 1;
594
595 return;
596}
597
598
599void
600Lattice3dBoundaryTruss :: saveContext(DataStream &stream, ContextMode mode)
601{
602 Lattice3d :: saveContext(stream, mode);
603
605
606 if ( ( mode & CM_Definition ) ) {
607 if ( ( iores = location.storeYourself(stream) ) != CIO_OK ) {
608 THROW_CIOERR(iores);
609 }
610 }
611}
612
613
614
615void
616Lattice3dBoundaryTruss :: restoreContext(DataStream &stream, ContextMode mode)
617{
618 Lattice3d :: restoreContext(stream, mode);
619
621
622 if ( mode & CM_Definition ) {
623 if ( ( iores = this->location.restoreYourself(stream) ) != CIO_OK ) {
624 THROW_CIOERR(iores);
625 }
626 }
627}
628
629
630
631#ifdef __OOFEG
632
633void
634Lattice3dBoundaryTruss :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
635{
636 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
637
638 if ( mode == OGC_rawGeometry ) {
639 this->drawRawGeometry(gc, tStep);
640 this->drawRawCrossSections(gc, tStep);
641 } else if ( mode == OGC_deformedGeometry ) {
642 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
643 } else if ( mode == OGC_elemSpecial ) {
644 this->drawSpecial(gc, tStep);
645 } else {
646 OOFEM_ERROR("drawYourself : unsupported mode");
647 }
648}
649
650
651
652
653void Lattice3dBoundaryTruss :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
654{
655 GraphicObj *go;
656
657 WCRec p [ 2 ]; /* points */
658 if ( !gc.testElementGraphicActivity(this) ) {
659 return;
660 }
661
662 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
663 EASValsSetColor(gc.getActiveCrackColor() );
664 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
665
666
667 FloatArray specimenDimension(3);
668 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
669 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
670 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
671
672 //Compute projection vector
673 IntArray projectionComponentNodeOne(3);
674 projectionComponentNodeOne.zero();
675 if ( location.at(1) != 0 ) {
676 giveSwitches(projectionComponentNodeOne, location.at(1) );
677 }
678
679 IntArray projectionComponentNodeTwo(3);
680 projectionComponentNodeTwo.zero();
681
682 if ( location.at(2) != 0 ) {
683 giveSwitches(projectionComponentNodeTwo, location.at(2) );
684 }
685
686 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
687 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
688 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
689 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
690 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
691 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
692
693 go = CreateLine3D(p);
694 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
695 EGAttachObject(go, ( EObjectP ) this);
696 EMAddGraphicsToModel(ESIModel(), go);
697}
698
699
700void
701Lattice3dBoundaryTruss :: drawRawCrossSections(oofegGraphicContext &gc, TimeStep *tStep)
702{
703 GraphicObj *go;
704
705 // if (!go) { // create new one
706 //Create as many points as we have polygon vertices
707 numberOfPolygonVertices = this->polygonCoords.giveSize() / 3.;
708 WCRec p [ numberOfPolygonVertices ]; /* poin */
709
710 if ( !gc.testElementGraphicActivity(this) ) {
711 return;
712 }
713
714 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
715 // EASValsSetColor(gc.getNodeColor());
716 EASValsSetLayer(OOFEG_RAW_CROSSSECTION_LAYER);
717
718 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
719 p [ i ].x = ( FPNum ) polygonCoords(3 * i);
720 p [ i ].y = ( FPNum ) polygonCoords(3 * i + 1);
721 p [ i ].z = ( FPNum ) polygonCoords(3 * i + 2);
722 }
723
724
725 WCRec pTemp [ 2 ]; /* points */
726 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
727 if ( i < numberOfPolygonVertices - 1 ) {
728 pTemp [ 0 ] = p [ i ];
729 pTemp [ 1 ] = p [ i + 1 ];
730 } else {
731 pTemp [ 0 ] = p [ i ];
732 pTemp [ 1 ] = p [ 0 ];
733 }
734
735 go = CreateLine3D(pTemp);
736 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
737 EGAttachObject(go, ( EObjectP ) this);
738 EMAddGraphicsToModel(ESIModel(), go);
739 }
740}
741
742
743void Lattice3dBoundaryTruss :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
744{
745 //That seems to be wrong. The strain field should be ordered exx, eyy, ezz, gyz, gzx, gyx
746 //Therefore, the x displacement should include 5th and 6th strain components.
747 GraphicObj *go;
748
749 if ( !gc.testElementGraphicActivity(this) ) {
750 return;
751 }
752
753 double defScale = gc.getDefScale();
754
755 WCRec p [ 2 ]; /* points */
756
757 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
758 EASValsSetColor(gc.getDeformedElementColor() );
759 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
760
761 FloatArray specimenDimension(3);
762 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
763 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
764 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
765
766 FloatArray dispOne(6);
767 dispOne.at(1) = this->giveNode(1)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
768 dispOne.at(2) = this->giveNode(1)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
769 dispOne.at(3) = this->giveNode(1)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
770 dispOne.at(4) = this->giveNode(1)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
771 dispOne.at(5) = this->giveNode(1)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
772 dispOne.at(6) = this->giveNode(1)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
773
774 FloatArray dispTwo(6);
775 dispTwo.at(1) = this->giveNode(2)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
776 dispTwo.at(2) = this->giveNode(2)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
777 dispTwo.at(3) = this->giveNode(2)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
778 dispTwo.at(4) = this->giveNode(2)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
779 dispTwo.at(5) = this->giveNode(2)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
780 dispTwo.at(6) = this->giveNode(2)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
781
782 FloatArray dispThree(1);
783 dispThree.at(1) = this->giveNode(3)->giveDofWithID(E_xx)->giveUnknown(VM_Total, tStep);
784
785 IntArray projectionComponentNodeOne(3);
786 projectionComponentNodeOne.zero();
787 if ( location.at(1) != 0 ) {
788 giveSwitches(projectionComponentNodeOne, location.at(1) );
789 }
790
791 IntArray projectionComponentNodeTwo(3);
792 projectionComponentNodeTwo.zero();
793 if ( location.at(2) != 0 ) {
794 giveSwitches(projectionComponentNodeTwo, location.at(2) );
795 }
796
797 //Modify dispOne and dispTwo
798 dispOne.at(1) = dispOne.at(1) + projectionComponentNodeOne.at(1) * dispThree.at(1);
799 dispOne.at(2) = dispOne.at(2);
800 dispOne.at(3) = dispOne.at(3);
801
802
803 dispTwo.at(1) = dispTwo.at(1) + projectionComponentNodeTwo.at(1) * dispThree.at(1);
804 dispTwo.at(2) = dispTwo.at(2);
805 dispTwo.at(3) = dispTwo.at(3);
806
807 double x1, y1, z1, x2, y2, z2;
808 x1 = this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
809 y1 = this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
810 z1 = this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
811
812 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
813 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
814 z2 = this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
815
816 p [ 0 ].x = ( FPNum ) x1 + defScale * dispOne.at(1);
817 p [ 0 ].y = ( FPNum ) y1 + defScale * dispOne.at(2);
818 p [ 0 ].z = ( FPNum ) z1 + defScale * dispOne.at(3);
819
820 p [ 1 ].x = ( FPNum ) x2 + defScale * dispTwo.at(1);
821 p [ 1 ].y = ( FPNum ) y2 + defScale * dispTwo.at(2);
822 p [ 1 ].z = ( FPNum ) z2 + defScale * dispTwo.at(3);
823
824 go = CreateLine3D(p);
825 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
826 EMAddGraphicsToModel(ESIModel(), go);
827}
828
829#endif
830} // 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
virtual Material * giveMaterial()
Definition element.C:523
virtual MaterialMode giveMaterialMode()
Definition element.h:738
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
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 rotatedWith(FloatMatrix &r, char mode)
Definition floatarray.C:814
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void subtract(const FloatArray &src)
Definition floatarray.C:320
void times(double s)
Definition floatarray.C:834
void times(double f)
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 giveSwitches(IntArray &answer, int location)
bool computeGtoLRotationMatrix(FloatMatrix &) override
double computeVolumeAround(GaussPoint *) override
int giveLocalCoordinateSystem(FloatMatrix &answer) override
void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
static ParamKey IPK_Lattice3dBoundaryTruss_location
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN) override
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
void drawRawCrossSections(oofegGraphicContext &, TimeStep *tStep)
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
Lattice3dBoundary(int n, Domain *)
FloatArray normal
Definition lattice3d.h:66
FloatArray midPoint
Definition lattice3d.h:66
FloatArray polygonCoords
Definition lattice3d.h:62
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition lattice3d.C:268
virtual void computeCrossSectionProperties()
Definition lattice3d.C:478
FloatMatrix localCoordinateSystem
Definition lattice3d.h:64
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition lattice3d.C:262
int numberOfPolygonVertices
Definition lattice3d.h:63
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Definition material.C:206
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
#define OOFEG_RAW_CROSSSECTION_LAYER
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