OOFEM 3.0
Loading...
Searching...
No Matches
lattice3dboundary.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"
36#include "lattice3dboundary.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
64Lattice3dBoundary :: Lattice3dBoundary(int n, Domain *aDomain) : Lattice3d(n, aDomain), location(2)
65{
67 geometryFlag = 0;
68}
69
70Lattice3dBoundary :: ~Lattice3dBoundary()
71{}
72
73
74void
75Lattice3dBoundary :: 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
188Lattice3dBoundary :: 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, 18), tt;
195 FloatMatrix answerTemp(12, 12), answerHelp, ttk(18, 12);
196 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
197 answerTemp.zero();
198 answerHelp.zero();
199 t.zero();
200
201
202
203 if ( geometryFlag == 0 ) {
205 }
206
207 double volume = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
208
209 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), bj);
210 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
211 for ( int i = 1; i <= 6; i++ ) {
212 d.at(i, i) *= volume;
213 }
214
215 dbj.beProductOf(d, bj);
216 bjt.beTranspositionOf(bj);
217 answerTemp.beProductOf(bjt, dbj);
218
220 answer.zero();
221
223 answerHelp.zero();
224
225 for ( int m = 1; m <= 12; m++ ) {
226 for ( int k = 1; k <= 12; k++ ) {
227 answer.at(m, k) = answerTemp.at(m, k);
228 }
229 }
230
231 if ( matStiffSymmFlag ) {
232 answer.symmetrized();
233 }
234
235 // Rotate to global system add the projections and rotate back to local system.
236
237 FloatMatrix R;
238 if ( this->giveRotationMatrix(R) ) {
239 answer.rotatedWith(R);
240 }
241
242 for ( int m = 1; m <= 12; m++ ) {
243 for ( int k = 1; k <= 12; k++ ) {
244 answerTemp.at(m, k) = answer.at(m, k);
245 }
246 }
247
248 IntArray projectionComponentNodeOne(3);
249 projectionComponentNodeOne.zero();
250 if ( location.at(1) != 0 ) {
251 giveSwitches(projectionComponentNodeOne, location.at(1) );
252 }
253
254 IntArray projectionComponentNodeTwo(3);
255 projectionComponentNodeTwo.zero();
256 if ( location.at(2) != 0 ) {
257 giveSwitches(projectionComponentNodeTwo, location.at(2) );
258 }
259
260 for ( int k = 1; k <= 12; k++ ) {
261 t.at(k, k) = 1.;
262 }
263 t.at(1, 13) = projectionComponentNodeOne.at(1);
264 t.at(2, 14) = projectionComponentNodeOne.at(2);
265 t.at(3, 15) = projectionComponentNodeOne.at(3);
266
267
268 t.at(2, 16) = projectionComponentNodeOne.at(3);
269 t.at(1, 17) = projectionComponentNodeOne.at(3);
270 t.at(1, 18) = projectionComponentNodeOne.at(2);
271
272
273 t.at(7, 13) = projectionComponentNodeTwo.at(1);
274 t.at(8, 14) = projectionComponentNodeTwo.at(2);
275 t.at(9, 15) = projectionComponentNodeTwo.at(3);
276
277 t.at(8, 16) = projectionComponentNodeTwo.at(3);
278 t.at(7, 17) = projectionComponentNodeTwo.at(3);
279 t.at(7, 18) = projectionComponentNodeTwo.at(2);
280
281
282 tt.beTranspositionOf(t);
283
284 ttk.beProductOf(tt, answerTemp);
285 answer.beProductOf(ttk, t);
286
287 //Rotate back to local system
288 FloatMatrix Rtranspose;
289 Rtranspose.beTranspositionOf(R);
290 answer.rotatedWith(Rtranspose);
291
292 return;
293}
294
295
296
297double
298Lattice3dBoundary :: computeVolumeAround(GaussPoint *aGaussPoint)
299{
300 if ( geometryFlag == 0 ) {
302 }
303
304 return this->area * this->length;
305}
306
307void
308Lattice3dBoundary :: recalculateCoordinates(int nodeNumber, FloatArray &coords) {
309 coords.resize(3);
310 coords.zero();
311 Node *node;
312
313 FloatArray specimenDimension(3);
314 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
315 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
316 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
317
318 IntArray projectionComponent(3);
319 projectionComponent.zero();
320
321 if ( nodeNumber == 1 ) {
322 node = this->giveNode(1);
323 if ( location.at(1) != 0 ) {
324 giveSwitches(projectionComponent, location.at(1) );
325 }
326 } else if ( nodeNumber == 2 ) {
327 node = this->giveNode(2);
328 if ( location.at(2) != 0 ) {
329 giveSwitches(projectionComponent, location.at(2) );
330 }
331 } else {
332 OOFEM_ERROR("wrong element used in the vtk module");
333 }
334
335 for ( int i = 0; i < 3; i++ ) {
336 coords.at(i + 1) = node->giveCoordinate(i + 1) + projectionComponent.at(i + 1) * specimenDimension.at(i + 1);
337 }
338 return;
339}
340
341void
342Lattice3dBoundary :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN)
343// Computes the vector containing the strains at the Gauss point gp of
344// the receiver, at time step stepN. The nature of these strains depends
345// on the element's type.
346{
347 FloatMatrix b;
348 FloatArray u;
349
350
351 //Compute strain vector
352 //Get the 18 components of the displacement vector of this element
353 this->computeVectorOf(VM_Total, stepN, u);
354 this->computeBmatrixAt(gp, b);
355
356
357 // subtract initial displacements, if defined
358 if ( initialDisplacements ) {
360 }
361
362 //Compute projection vector
363 IntArray projectionComponentNodeOne(3);
364 projectionComponentNodeOne.zero();
365 if ( location.at(1) != 0 ) {
366 giveSwitches(projectionComponentNodeOne, location.at(1) );
367 }
368
369 IntArray projectionComponentNodeTwo(3);
370 projectionComponentNodeTwo.zero();
371
372 if ( location.at(2) != 0 ) {
373 giveSwitches(projectionComponentNodeTwo, location.at(2) );
374 }
375
376 FloatMatrix rotationMatrix;
377 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
378 u.rotatedWith(rotationMatrix, 't');
379 }
380
381 //General expressions for the corrected displacements. Rotations at nodes are not influenced. Only translations
382 u.at(1) = u.at(1) + projectionComponentNodeOne.at(1) * u.at(13) + projectionComponentNodeOne.at(2) * u.at(18) + projectionComponentNodeOne.at(3) * u.at(17);
383 u.at(2) = u.at(2) + projectionComponentNodeOne.at(2) * u.at(14) + projectionComponentNodeOne.at(3) * u.at(16);
384 u.at(3) = u.at(3) + projectionComponentNodeOne.at(3) * u.at(15);
385
386 u.at(7) = u.at(7) + projectionComponentNodeTwo.at(1) * u.at(13) + projectionComponentNodeTwo.at(2) * u.at(18) + projectionComponentNodeTwo.at(3) * u.at(17);
387 u.at(8) = u.at(8) + projectionComponentNodeTwo.at(2) * u.at(14) + projectionComponentNodeTwo.at(3) * u.at(16);
388 u.at(9) = u.at(9) + projectionComponentNodeTwo.at(3) * u.at(15);
389
390
391 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
392 u.rotatedWith(rotationMatrix, 'n');
393 }
394 //define a temp displacement vector
395 FloatArray uTemp(12);
396
397 for ( int i = 1; i <= 12; i++ ) {
398 uTemp.at(i) = u.at(i);
399 }
400
401 answer.beProductOf(b, uTemp);
402}
403
404bool
405Lattice3dBoundary :: computeGtoLRotationMatrix(FloatMatrix &answer)
406{
407 FloatMatrix lcs;
408 int i, j;
409
410 answer.resize(18, 18);
411 answer.zero();
412
413 this->giveLocalCoordinateSystem(lcs);
414 for ( i = 1; i <= 3; i++ ) {
415 for ( j = 1; j <= 3; j++ ) {
416 answer.at(i, j) = lcs.at(i, j);
417 answer.at(i + 3, j + 3) = lcs.at(i, j);
418 answer.at(i + 6, j + 6) = lcs.at(i, j);
419 answer.at(i + 9, j + 9) = lcs.at(i, j);
420 }
421 }
422 answer.at(13, 13) = 1.;
423 answer.at(14, 14) = 1.;
424 answer.at(15, 15) = 1.;
425 answer.at(16, 16) = 1.;
426 answer.at(17, 17) = 1.;
427 answer.at(18, 18) = 1.;
428
429 return 1;
430}
431
432int
433Lattice3dBoundary :: giveLocalCoordinateSystem(FloatMatrix &answer)
434{
435 if ( geometryFlag == 0 ) {
437 }
438
439 answer = this->localCoordinateSystem;
440
441 return 1;
442}
443
444
445
446void
447Lattice3dBoundary :: giveDofManDofIDMask(int inode, IntArray &answer) const
448{
449 if ( inode == 3 ) {
450 answer = { E_xx, E_yy, E_zz, G_yz, G_xz, G_xy };
451 } else {
452 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
453 }
454}
455
456void
457Lattice3dBoundary :: initializeFrom(InputRecord &ir, int priority)
458{
459 ParameterManager &ppm = this->giveDomain()->elementPPM;
460 Lattice3d :: initializeFrom(ir, priority);
462}
463
464void
465Lattice3dBoundary :: postInitialize()
466{
467 ParameterManager &ppm = this->giveDomain()->elementPPM;
468 Lattice3d :: postInitialize();
470}
471
472void
473Lattice3dBoundary :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
474{
475 Material *mat = this->giveMaterial();
476
477 FloatMatrix b, bt, A, R, GNT;
478 FloatArray bs, TotalStressVector, u, strain;
479 double dV;
480
481 this->computeVectorOf(VM_Total, tStep, u);
482
483 // subtract initial displacements, if defined
484 if ( initialDisplacements ) {
486 }
487
488 answer.resize(18);
489 answer.zero();
490
491 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), b);
492
493 bt.beTranspositionOf(b);
494 if ( useUpdatedGpRecord == 1 ) {
495 TotalStressVector = ( ( LatticeMaterialStatus * ) mat->giveStatus(integrationRulesArray [ 0 ]->getIntegrationPoint(0) ) )
496 ->giveLatticeStress();
497 } else
498 if ( !this->isActivated(tStep) ) {
499 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(integrationRulesArray [ 0 ]->getIntegrationPoint(0)->giveMaterialMode() ) );
500 strain.zero();
501 }
502 this->computeStrainVector(strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
503
504 this->computeStressVector(TotalStressVector, strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
505
506 dV = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
507 bs.beProductOf(bt, TotalStressVector);
508 bs.times(dV);
509
510 for ( int m = 1; m <= 12; m++ ) {
511 answer.at(m) = bs.at(m);
512 }
513
514 if ( this->giveRotationMatrix(R) ) {
515 answer.rotatedWith(R, 't');
516 }
517
518 IntArray projectionComponentNodeOne(3);
519 projectionComponentNodeOne.zero();
520 if ( location.at(1) != 0 ) {
521 giveSwitches(projectionComponentNodeOne, location.at(1) );
522 }
523
524 IntArray projectionComponentNodeTwo(3);
525 projectionComponentNodeTwo.zero();
526
527 if ( location.at(2) != 0 ) {
528 giveSwitches(projectionComponentNodeTwo, location.at(2) );
529 }
530
531 //Normal stresses
532 answer.at(13) = projectionComponentNodeOne.at(1) * answer.at(1) + projectionComponentNodeTwo.at(1) * answer.at(7);
533 answer.at(14) = projectionComponentNodeOne.at(2) * answer.at(2) + projectionComponentNodeTwo.at(2) * answer.at(8);
534 answer.at(15) = projectionComponentNodeOne.at(3) * answer.at(3) + projectionComponentNodeTwo.at(3) * answer.at(9);
535
536 //Shear stresses
537 answer.at(16) = projectionComponentNodeOne.at(3) * answer.at(2) + projectionComponentNodeTwo.at(3) * answer.at(8);
538 answer.at(17) = projectionComponentNodeOne.at(3) * answer.at(1) + projectionComponentNodeTwo.at(3) * answer.at(7);
539 answer.at(18) = projectionComponentNodeOne.at(2) * answer.at(1) + projectionComponentNodeTwo.at(2) * answer.at(7);
540
541 //Rotate to local system
542 if ( this->giveRotationMatrix(R) ) {
543 answer.rotatedWith(R, 'n');
544 }
545
546 return;
547}
548
549
550void
551Lattice3dBoundary :: giveSwitches(IntArray &answer, int location) {
552 int counter = 1;
553 for ( int x = -1; x < 2; x++ ) {
554 for ( int y = -1; y < 2; y++ ) {
555 for ( int z = -1; z < 2; z++ ) {
556 if ( !( z == 0 && y == 0 && x == 0 ) ) {
557 if ( counter == location ) {
558 answer(0) = x;
559 answer(1) = y;
560 answer(2) = z;
561 }
562 counter++;
563 }
564 }
565 }
566 }
567 return;
568}
569
570
571void
572Lattice3dBoundary :: computeGeometryProperties()
573{
574 Node *nodeA, *nodeB;
575
576 FloatArray coordsA(3);
577 FloatArray coordsB(3);
578
579 nodeA = this->giveNode(1);
580 nodeB = this->giveNode(2);
581
582 FloatArray specimenDimension(3);
583 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
584 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
585 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
586
587 IntArray projectionComponentNodeOne(3);
588 projectionComponentNodeOne.zero();
589 if ( location.at(1) != 0 ) {
590 giveSwitches(projectionComponentNodeOne, location.at(1) );
591 }
592
593 IntArray projectionComponentNodeTwo(3);
594 projectionComponentNodeTwo.zero();
595 if ( location.at(2) != 0 ) {
596 giveSwitches(projectionComponentNodeTwo, location.at(2) );
597 }
598
599 for ( int i = 0; i < 3; i++ ) {
600 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1) + projectionComponentNodeOne.at(i + 1) * specimenDimension.at(i + 1);
601 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1) + projectionComponentNodeTwo.at(i + 1) * specimenDimension.at(i + 1);
602 }
603
604 //Calculate normal vector
605 FloatArray s(3), t(3);
606 this->midPoint.resize(3);
607
608 //Calculate normal vector
609 this->normal.resize(3);
610 for ( int i = 0; i < 3; i++ ) {
611 this->normal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
612 }
613
614 this->length = sqrt(pow(this->normal.at(1), 2.) + pow(this->normal.at(2), 2.) + pow(this->normal.at(3), 2.) );
615
616 for ( int i = 0; i < 3; i++ ) {
617 this->normal.at(i + 1) /= length;
618 }
619
620 // Compute midpoint
621 this->midPoint.resize(3);
622 for ( int i = 0; i < 3; i++ ) {
623 this->midPoint.at(i + 1) = 0.5 * ( coordsB.at(i + 1) + coordsA.at(i + 1) );
624 }
625
627
628 //Set geometry flag to 1 so that this is done only once
629 this->geometryFlag = 1;
630
631 return;
632}
633
634
635void
636Lattice3dBoundary :: saveContext(DataStream &stream, ContextMode mode)
637{
638 Lattice3d :: saveContext(stream, mode);
639
641
642 if ( ( mode & CM_Definition ) ) {
643 if ( ( iores = location.storeYourself(stream) ) != CIO_OK ) {
644 THROW_CIOERR(iores);
645 }
646 }
647}
648
649
650
651void
652Lattice3dBoundary :: restoreContext(DataStream &stream, ContextMode mode)
653{
654 Lattice3d :: restoreContext(stream, mode);
655
657
658 if ( mode & CM_Definition ) {
659 if ( ( iores = this->location.restoreYourself(stream) ) != CIO_OK ) {
660 THROW_CIOERR(iores);
661 }
662 }
663}
664
665
666
667#ifdef __OOFEG
668
669void
670Lattice3dBoundary :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
671{
672 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
673
674 if ( mode == OGC_rawGeometry ) {
675 this->drawRawGeometry(gc, tStep);
676 this->drawRawCrossSections(gc, tStep);
677 } else if ( mode == OGC_deformedGeometry ) {
678 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
679 } else if ( mode == OGC_elemSpecial ) {
680 this->drawSpecial(gc, tStep);
681 } else {
682 OOFEM_ERROR("drawYourself : unsupported mode");
683 }
684}
685
686
687
688
689void Lattice3dBoundary :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
690{
691 GraphicObj *go;
692
693 WCRec p [ 2 ]; /* points */
694 if ( !gc.testElementGraphicActivity(this) ) {
695 return;
696 }
697
698 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
699 EASValsSetColor(gc.getActiveCrackColor() );
700 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
701
702
703 FloatArray specimenDimension(3);
704 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
705 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
706 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
707
708 //Compute projection vector
709 IntArray projectionComponentNodeOne(3);
710 projectionComponentNodeOne.zero();
711 if ( location.at(1) != 0 ) {
712 giveSwitches(projectionComponentNodeOne, location.at(1) );
713 }
714
715 IntArray projectionComponentNodeTwo(3);
716 projectionComponentNodeTwo.zero();
717
718 if ( location.at(2) != 0 ) {
719 giveSwitches(projectionComponentNodeTwo, location.at(2) );
720 }
721
722 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
723 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
724 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
725 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
726 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
727 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
728
729 go = CreateLine3D(p);
730 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
731 EGAttachObject(go, ( EObjectP ) this);
732 EMAddGraphicsToModel(ESIModel(), go);
733}
734
735
736void
737Lattice3dBoundary :: drawRawCrossSections(oofegGraphicContext &gc, TimeStep *tStep)
738{
739 GraphicObj *go;
740
741 // if (!go) { // create new one
742 //Create as many points as we have polygon vertices
743 numberOfPolygonVertices = this->polygonCoords.giveSize() / 3.;
744 WCRec p [ numberOfPolygonVertices ]; /* poin */
745
746 if ( !gc.testElementGraphicActivity(this) ) {
747 return;
748 }
749
750 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
751 // EASValsSetColor(gc.getNodeColor());
752 EASValsSetLayer(OOFEG_RAW_CROSSSECTION_LAYER);
753
754 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
755 p [ i ].x = ( FPNum ) polygonCoords(3 * i);
756 p [ i ].y = ( FPNum ) polygonCoords(3 * i + 1);
757 p [ i ].z = ( FPNum ) polygonCoords(3 * i + 2);
758 }
759
760
761 WCRec pTemp [ 2 ]; /* points */
762 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
763 if ( i < numberOfPolygonVertices - 1 ) {
764 pTemp [ 0 ] = p [ i ];
765 pTemp [ 1 ] = p [ i + 1 ];
766 } else {
767 pTemp [ 0 ] = p [ i ];
768 pTemp [ 1 ] = p [ 0 ];
769 }
770
771 go = CreateLine3D(pTemp);
772 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
773 EGAttachObject(go, ( EObjectP ) this);
774 EMAddGraphicsToModel(ESIModel(), go);
775 }
776}
777
778
779void Lattice3dBoundary :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
780{
781 //That seems to be wrong. The strain field should be ordered exx, eyy, ezz, gyz, gzx, gyx
782 //Therefore, the x displacement should include 5th and 6th strain components.
783 GraphicObj *go;
784
785 if ( !gc.testElementGraphicActivity(this) ) {
786 return;
787 }
788
789 double defScale = gc.getDefScale();
790
791 WCRec p [ 2 ]; /* points */
792
793 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
794 EASValsSetColor(gc.getDeformedElementColor() );
795 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
796
797 FloatArray specimenDimension(3);
798 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
799 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
800 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
801
802 FloatArray dispOne(6);
803 dispOne.at(1) = this->giveNode(1)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
804 dispOne.at(2) = this->giveNode(1)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
805 dispOne.at(3) = this->giveNode(1)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
806 dispOne.at(4) = this->giveNode(1)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
807 dispOne.at(5) = this->giveNode(1)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
808 dispOne.at(6) = this->giveNode(1)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
809
810 FloatArray dispTwo(6);
811 dispTwo.at(1) = this->giveNode(2)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
812 dispTwo.at(2) = this->giveNode(2)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
813 dispTwo.at(3) = this->giveNode(2)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
814 dispTwo.at(4) = this->giveNode(2)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
815 dispTwo.at(5) = this->giveNode(2)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
816 dispTwo.at(6) = this->giveNode(2)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
817
818 FloatArray dispThree(6);
819 dispThree.at(1) = this->giveNode(3)->giveDofWithID(E_xx)->giveUnknown(VM_Total, tStep);
820 dispThree.at(2) = this->giveNode(3)->giveDofWithID(E_yy)->giveUnknown(VM_Total, tStep);
821 dispThree.at(3) = this->giveNode(3)->giveDofWithID(E_zz)->giveUnknown(VM_Total, tStep);
822 dispThree.at(4) = this->giveNode(3)->giveDofWithID(G_yz)->giveUnknown(VM_Total, tStep);
823 dispThree.at(5) = this->giveNode(3)->giveDofWithID(G_xz)->giveUnknown(VM_Total, tStep);
824 dispThree.at(6) = this->giveNode(3)->giveDofWithID(G_xy)->giveUnknown(VM_Total, tStep);
825
826 IntArray projectionComponentNodeOne(3);
827 projectionComponentNodeOne.zero();
828 if ( location.at(1) != 0 ) {
829 giveSwitches(projectionComponentNodeOne, location.at(1) );
830 }
831
832 IntArray projectionComponentNodeTwo(3);
833 projectionComponentNodeTwo.zero();
834 if ( location.at(2) != 0 ) {
835 giveSwitches(projectionComponentNodeTwo, location.at(2) );
836 }
837
838 //Modify dispOne and dispTwo
839 //Seems to be wrong. Should be
840 dispOne.at(1) = dispOne.at(1) + projectionComponentNodeOne.at(1) * dispThree.at(1) + projectionComponentNodeOne.at(3) * dispThree.at(5) + projectionComponentNodeOne.at(2) * dispThree.at(6);
841 dispOne.at(2) = dispOne.at(2) + projectionComponentNodeOne.at(2) * dispThree.at(2) + projectionComponentNodeOne.at(3) * dispThree.at(4);
842 dispOne.at(3) = dispOne.at(3) + projectionComponentNodeOne.at(3) * dispThree.at(3);
843
844
845 dispTwo.at(1) = dispTwo.at(1) + projectionComponentNodeTwo.at(1) * dispThree.at(1) + projectionComponentNodeTwo.at(3) * dispThree.at(5) + projectionComponentNodeTwo.at(2) * dispThree.at(6);
846 dispTwo.at(2) = dispTwo.at(2) + projectionComponentNodeTwo.at(2) * dispThree.at(2) + projectionComponentNodeTwo.at(3) * dispThree.at(4);
847 dispTwo.at(3) = dispTwo.at(3) + projectionComponentNodeTwo.at(3) * dispThree.at(3);
848
849 double x1, y1, z1, x2, y2, z2;
850 x1 = this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
851 y1 = this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
852 z1 = this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
853
854 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
855 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
856 z2 = this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
857
858 p [ 0 ].x = ( FPNum ) x1 + defScale * dispOne.at(1);
859 p [ 0 ].y = ( FPNum ) y1 + defScale * dispOne.at(2);
860 p [ 0 ].z = ( FPNum ) z1 + defScale * dispOne.at(3);
861
862 p [ 1 ].x = ( FPNum ) x2 + defScale * dispTwo.at(1);
863 p [ 1 ].y = ( FPNum ) y2 + defScale * dispTwo.at(2);
864 p [ 1 ].z = ( FPNum ) z2 + defScale * dispTwo.at(3);
865
866 go = CreateLine3D(p);
867 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
868 EMAddGraphicsToModel(ESIModel(), go);
869}
870
871#endif
872} // 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 drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
int computeNumberOfDofs() override
int giveLocalCoordinateSystem(FloatMatrix &answer) override
static ParamKey IPK_Lattice3dBoundary_location
double computeVolumeAround(GaussPoint *) override
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN) override
void drawRawCrossSections(oofegGraphicContext &, TimeStep *tStep)
void computeGeometryProperties() override
void giveSwitches(IntArray &answer, int location)
bool computeGtoLRotationMatrix(FloatMatrix &) override
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
FloatArray normal
Definition lattice3d.h:66
Lattice3d(int n, Domain *)
Definition lattice3d.C:70
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