OOFEM 3.0
Loading...
Searching...
No Matches
lattice3d.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 "lattice3d.h"
38#include "node.h"
39#include "material.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "floatmatrixf.h"
44#include "intarray.h"
45#include "floatarray.h"
46#include "floatarrayf.h"
47#include "mathfem.h"
49#include "contextioerr.h"
50#include "datastream.h"
51#include "classfactory.h"
53#include "parametermanager.h"
54#include "paramkey.h"
55
56#ifdef __OOFEG
57 #include "oofeggraphiccontext.h"
59#endif
60
61namespace oofem {
63
69
70Lattice3d :: Lattice3d(int n, Domain *aDomain) : LatticeStructuralElement(n, aDomain)
71{
73 geometryFlag = 0;
74 this->eccS = 0.;
75 this->eccT = 0.;
76 minLength = 1.e-20;
77 couplingFlag = 0;
78
79}
80
81Lattice3d :: ~Lattice3d()
82{}
83
84
85void
86Lattice3d :: computeBmatrixAt(GaussPoint *aGaussPoint, FloatMatrix &answer, int li, int ui)
87// Returns the strain matrix of the receiver.
88{
89 if ( geometryFlag == 0 ) {
91 }
92
93 //Assemble Bmatrix (used to compute strains and rotations}
94 answer.resize(6, 12);
95 answer.zero();
96
97 //Normal displacement jump in x-direction
98 //First node
99 answer.at(1, 1) = -1.;
100 answer.at(1, 2) = 0.;
101 answer.at(1, 3) = 0.;
102 answer.at(1, 4) = 0.;
103 answer.at(1, 5) = -this->eccT;
104 answer.at(1, 6) = this->eccS;
105 //Second node
106 answer.at(1, 7) = 1.;
107 answer.at(1, 8) = 0.;
108 answer.at(1, 9) = 0.;
109 answer.at(1, 10) = 0.;
110 answer.at(1, 11) = this->eccT;
111 answer.at(1, 12) = -this->eccS;
112
113 //Shear displacement jump in y-plane
114 //first node
115 answer.at(2, 1) = 0.;
116 answer.at(2, 2) = -1.;
117 answer.at(2, 3) = 0.;
118 answer.at(2, 4) = this->eccT;
119 answer.at(2, 5) = 0;
120 answer.at(2, 6) = -this->length / 2.;
121 //Second node
122 answer.at(2, 7) = 0.;
123 answer.at(2, 8) = 1.;
124 answer.at(2, 9) = 0.;
125 answer.at(2, 10) = -this->eccT;
126 answer.at(2, 11) = 0;
127 answer.at(2, 12) = -this->length / 2.;
128
129 //Shear displacement jump in z-plane
130 //first node
131 answer.at(3, 1) = 0.;
132 answer.at(3, 2) = 0.;
133 answer.at(3, 3) = -1.;
134 answer.at(3, 4) = -this->eccS;
135 answer.at(3, 5) = this->length / 2.;
136 answer.at(3, 6) = 0.;
137 //Second node
138 answer.at(3, 7) = 0.;
139 answer.at(3, 8) = 0.;
140 answer.at(3, 9) = 1.;
141 answer.at(3, 10) = this->eccS;
142 answer.at(3, 11) = this->length / 2.;
143 answer.at(3, 12) = 0.;
144
145 //Rotation around x-axis
146 //First node
147 answer.at(4, 1) = 0.;
148 answer.at(4, 2) = 0;
149 answer.at(4, 3) = 0.;
150 answer.at(4, 4) = -sqrt(Ip / this->area);
151 answer.at(4, 5) = 0.;
152 answer.at(4, 6) = 0.;
153 //Second node
154 answer.at(4, 7) = 0.;
155 answer.at(4, 8) = 0.;
156 answer.at(4, 9) = 0.;
157 answer.at(4, 10) = sqrt(Ip / this->area);
158 answer.at(4, 11) = 0.;
159 answer.at(4, 12) = 0.;
160
161 //Rotation around y-axis
162 //First node
163 answer.at(5, 1) = 0.;
164 answer.at(5, 2) = 0.;
165 answer.at(5, 3) = 0.;
166 answer.at(5, 4) = 0.;
167 answer.at(5, 5) = -sqrt(I1 / this->area);
168 answer.at(5, 6) = 0.;
169 //Second node
170 answer.at(5, 7) = 0.;
171 answer.at(5, 8) = 0.;
172 answer.at(5, 9) = 0.;
173 answer.at(5, 10) = 0.;
174 answer.at(5, 11) = sqrt(I1 / this->area);
175 answer.at(5, 12) = 0.;
176
177 //Rotation around z-axis
178 //First node
179 answer.at(6, 1) = 0.;
180 answer.at(6, 2) = 0.;
181 answer.at(6, 3) = 0.;
182 answer.at(6, 4) = 0.;
183 answer.at(6, 5) = 0.;
184 answer.at(6, 6) = -sqrt(I2 / this->area);
185 //Second node
186 answer.at(6, 7) = 0.;
187 answer.at(6, 8) = 0.;
188 answer.at(6, 9) = 0.;
189 answer.at(6, 10) = 0.;
190 answer.at(6, 11) = 0.;
191 answer.at(6, 12) = sqrt(I2 / this->area);
192
193 answer.times(1. / this->length);
194
195 return;
196}
197
198void
199Lattice3d :: giveGPCoordinates(FloatArray &coords)
200{
201 if ( geometryFlag == 0 ) {
203 }
204 coords.resize(3);
205 coords = this->globalCentroid;
206 return;
207}
208
209double
210Lattice3d :: giveLength()
211{
212 if ( geometryFlag == 0 ) {
214 }
215
216 return this->length;
217}
218
219
220int
221Lattice3d :: giveCrackFlag()
222{
223 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
224 LatticeMaterialStatus *status = static_cast< LatticeMaterialStatus * >( gp->giveMaterialStatus() );
225 int crackFlag = 0;
226 crackFlag = status->giveCrackFlag();
227
228 return crackFlag;
229}
230
231
232double
233Lattice3d :: giveCrackWidth()
234{
235 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
236 LatticeMaterialStatus *status = static_cast< LatticeMaterialStatus * >( gp->giveMaterialStatus() );
237 double crackWidth = 0;
238 crackWidth = status->giveCrackWidth();
239
240 return crackWidth;
241}
242
243void
244Lattice3d :: givePlasticStrain(FloatArray &plasticStrain)
245{
246 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
247 LatticeMaterialStatus *status = static_cast< LatticeMaterialStatus * >( this->giveMaterial()->giveStatus(gp) );
248 plasticStrain = status->givePlasticLatticeStrain();
249 return;
250}
251
252void
253Lattice3d :: giveOldPlasticStrain(FloatArray &plasticStrain)
254{
255 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
256 LatticeMaterialStatus *status = static_cast< LatticeMaterialStatus * >( this->giveMaterial()->giveStatus(gp) );
257 plasticStrain = status->giveOldPlasticLatticeStrain();
258 return;
259}
260
261void
262Lattice3d :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
263{
264 answer = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->give3dStiffnessMatrix(rMode, gp, tStep);
265}
266
267void
268Lattice3d :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
269{
270 answer = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->giveLatticeStress3d(strain, gp, tStep);
271}
272
273void
274Lattice3d :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
275 TimeStep *tStep)
276// Computes numerically the stiffness matrix of the receiver.
277{
278 FloatMatrix d, bi, bj, bjt, dbj, dij;
279
280 answer.resize(12, 12);
281 answer.zero();
282 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), bj);
283 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
284
285 double volume = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
286
287 for ( int i = 1; i <= 6; i++ ) {
288 d.at(i, i) *= volume;
289 }
290
291 dbj.beProductOf(d, bj);
292 bjt.beTranspositionOf(bj);
293 answer.beProductOf(bjt, dbj);
294
295 return;
296}
297
298void Lattice3d :: computeGaussPoints()
299// Sets up the array of Gauss Points of the receiver.
300{
301 integrationRulesArray.resize(1);
302 integrationRulesArray [ 0 ].reset(new GaussIntegrationRule(1, this, 1, 3) );
303 integrationRulesArray [ 0 ]->SetUpPointsOnLine(1, _3dLattice);
304}
305
306
307
308double Lattice3d :: giveArea() {
309 if ( geometryFlag == 0 ) {
311 }
312
313 return this->area;
314}
315
316
317bool
318Lattice3d :: computeGtoLRotationMatrix(FloatMatrix &answer)
319{
320 FloatMatrix lcs;
321 int i, j;
322
323 answer.resize(12, 12);
324 answer.zero();
325
326 this->giveLocalCoordinateSystem(lcs);
327 for ( i = 1; i <= 3; i++ ) {
328 for ( j = 1; j <= 3; j++ ) {
329 answer.at(i, j) = lcs.at(i, j);
330 answer.at(i + 3, j + 3) = lcs.at(i, j);
331 answer.at(i + 6, j + 6) = lcs.at(i, j);
332 answer.at(i + 9, j + 9) = lcs.at(i, j);
333 }
334 }
335
336 return 1;
337}
338
339
340double
341Lattice3d :: giveNormalStress()
342{
343 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
344 LatticeMaterialStatus *status = static_cast< LatticeMaterialStatus * >( gp->giveMaterialStatus() );
345 double normalStress = 0;
346 normalStress = status->giveNormalLatticeStress();
347
348 return normalStress;
349}
350
351
352int
353Lattice3d :: giveLocalCoordinateSystem(FloatMatrix &answer)
354{
355 if ( geometryFlag == 0 ) {
357 }
358
359 answer = this->localCoordinateSystem;
360
361 return 1;
362}
363
364
365double
366Lattice3d :: computeVolumeAround(GaussPoint *aGaussPoint)
367{
368 if ( geometryFlag == 0 ) {
370 }
371
372 return this->area * this->length;
373}
374
375void
376Lattice3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
377{
378 answer = {
379 D_u, D_v, D_w, R_u, R_v, R_w
380 };
381}
382
383void
384Lattice3d :: initializeFrom(InputRecord &ir, int priority)
385{
386 ParameterManager &ppm = this->giveDomain()->elementPPM;
387 LatticeStructuralElement :: initializeFrom(ir, priority);
388
389 PM_UPDATE_PARAMETER(minLength, ppm, ir, this->number, IPK_Lattice3d_mlength, priority);
394}
395
396void
397Lattice3d :: postInitialize()
398{
399 this->LatticeStructuralElement :: postInitialize();
400 numberOfPolygonVertices = (int) (polygonCoords.giveSize() / 3.);
401}
402
403int
404Lattice3d :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
405{
406 if ( geometryFlag == 0 ) {
408 }
409
410 answer.resize(3);
411 answer = this->globalCentroid;
412
413 return 1;
414}
415
416
417void
418Lattice3d :: computeGeometryProperties()
419{
420 //coordinates of the two nodes
421 Node *nodeA, *nodeB;
422 FloatArray coordsA(3), coordsB(3);
423
424 nodeA = this->giveNode(1);
425 nodeB = this->giveNode(2);
426
427 for ( int i = 0; i < 3; i++ ) {
428 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1);
429 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1);
430 }
431
432 //Construct an initial temporary local coordinate system
433 FloatArray s(3), t(3);
434
435 //Calculate normal vector
436 this->normal.resize(3);
437 for ( int i = 0; i < 3; i++ ) {
438 this->normal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
439 }
440
441 this->length = sqrt(pow(normal.at(1), 2.) + pow(normal.at(2), 2.) + pow(normal.at(3), 2.) );
442
443 // Compute midpoint
444 this->midPoint.resize(3);
445 for ( int i = 0; i < 3; i++ ) {
446 this->midPoint.at(i + 1) = 0.5 * ( coordsB.at(i + 1) + coordsA.at(i + 1) );
447 }
448
449 for ( int i = 0; i < 3; i++ ) {
450 this->normal.at(i + 1) /= length;
451 }
452
454
455 this->geometryFlag = 1;
456
457 return;
458}
459
460
461void
462Lattice3d :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
463// Returns the lumped mass matrix of the receiver. This expression is
464// valid in both local and global axes.
465{
466 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
467 double density = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->give('d', gp);
468 double halfMass = density * computeVolumeAround(gp) / 2.;
469 answer.resize(12, 12);
470 answer.zero();
471 answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass;
472 answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass;
473}
474
475
476
477void
478Lattice3d :: computeCrossSectionProperties() {
479 if ( this->numberOfPolygonVertices < 3 ) {
480 printf("Too small number of polygon vertices\n");
481 return;
482 }
483
484 //Construct two perpendicular axis so that n is normal to the plane which they create
485 //Check, if one of the components of the normal-direction is zero
486 FloatArray s(3), t(3);
487 if ( this->normal.at(1) == 0 ) {
488 s.at(1) = 0.;
489 s.at(2) = this->normal.at(3);
490 s.at(3) = -this->normal.at(2);
491 } else if ( this->normal.at(2) == 0 ) {
492 s.at(1) = this->normal.at(3);
493 s.at(2) = 0.;
494 s.at(3) = -this->normal.at(1);
495 } else {
496 s.at(1) = this->normal.at(2);
497 s.at(2) = -this->normal.at(1);
498 s.at(3) = 0.;
499 }
500
501 s.normalize();
502
503 t.beVectorProductOf(this->normal, s);
504 t.normalize();
505
506 //Set up rotation matrix
507 FloatMatrix lcs(3, 3);
508
509 for ( int i = 1; i <= 3; i++ ) {
510 lcs.at(1, i) = this->normal.at(i);
511 lcs.at(2, i) = s.at(i);
512 lcs.at(3, i) = t.at(i);
513 }
514
515
516 //Calculate the local coordinates of the polygon vertices
517 FloatArray help(3), test(3);
519 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
520 for ( int n = 0; n < 3; n++ ) {
521 help(n) = polygonCoords(3 * k + n);
522 }
523
524 test.beProductOf(lcs, help);
525 for ( int n = 0; n < 3; n++ ) {
526 lpc(3 * k + n) = test(n);
527 }
528 }
529
530 this->area = 0.;
531
532 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
533 if ( k < numberOfPolygonVertices - 1 ) {
534 this->area += lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2);
535 } else { //Back to zero for n+1
536 this->area += lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2);
537 }
538 }
539
540 this->area *= 0.5;
541
542 FloatArray tempCoords(3 * numberOfPolygonVertices);
543 if ( this->area < 0 ) { //Set area to a positive value and rearrange the coordinate entries
544 this->area *= -1.;
545
546 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
547 for ( int m = 0; m < 3; m++ ) {
548 tempCoords.at(3 * k + m + 1) = polygonCoords.at(3 * ( numberOfPolygonVertices - k - 1 ) + m + 1);
549 }
550 }
551
552 polygonCoords = tempCoords;
553
554 // Calculate again local co-ordinate system for different order
555 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
556 for ( int n = 0; n < 3; n++ ) {
557 help(n) = polygonCoords(3 * k + n);
558 }
559
560 test.beProductOf(lcs, help);
561 for ( int n = 0; n < 3; n++ ) {
562 lpc(3 * k + n) = test(n);
563 }
564 }
565 }
566
567 if ( this->area < pow(minLength, 2.) ) {
568 this->area = pow(minLength, 2.);
569 }
570
571 //Calculate centroids
572 centroid.resize(3);
573 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
574 if ( k < numberOfPolygonVertices - 1 ) {
575 centroid.at(2) += ( lpc(3 * k + 1) + lpc(3 * ( k + 1 ) + 1) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
576 centroid.at(3) += ( lpc(3 * k + 2) + lpc(3 * ( k + 1 ) + 2) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
577 } else { //Back to zero for n+1
578 centroid.at(2) += ( lpc(3 * k + 1) + lpc(1) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
579 centroid.at(3) += ( lpc(3 * k + 2) + lpc(2) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
580 }
581 }
582
583 centroid.times(1. / ( 6. * this->area ) );
584
585 centroid.at(1) = lpc.at(1); //The first component of all lpcs should be the same
586
587 //Shift coordinates to centroid
588 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
589 for ( int l = 0; l < 3; l++ ) {
590 lpc(3 * k + l) -= centroid(l);
591 }
592 }
593
594 //Compute second moments of area.
595 //This is for the temporary coordinate system
596 double Ixx = 0.;
597 double Iyy = 0.;
598 double Ixy = 0.;
599 double a;
600
601 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
602 if ( k < numberOfPolygonVertices - 1 ) {
603 a = lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2);
604
605 Ixx += ( ( pow(lpc(3 * k + 2), 2.) + lpc(3 * k + 2) * lpc(3 * ( k + 1 ) + 2) + pow(lpc(3 * ( k + 1 ) + 2), 2.) ) * a ) / 12.;
606
607 Iyy += ( ( pow(lpc(3 * k + 1), 2.) + lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 1) + pow(lpc(3 * ( k + 1 ) + 1), 2.) ) * a ) / 12.;
608
609 Ixy += ( ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) + 2. * lpc(3 * k + 1) * lpc(3 * k + 2) +
610 2 * lpc(3 * ( k + 1 ) + 1) * lpc(3 * ( k + 1 ) + 2) + lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) ) * a ) / 24.;
611 } else { //Back to zero for n+1
612 a = lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2);
613
614 Ixx += ( ( pow(lpc(3 * k + 2), 2.) + lpc(3 * k + 2) * lpc(2) + pow(lpc(2), 2.) ) * a ) / 12.;
615
616 Iyy += ( ( pow(lpc(3 * k + 1), 2.) + lpc(3 * k + 1) * lpc(1) + pow(lpc(1), 2.) ) * a ) / 12.;
617
618 Ixy += ( ( lpc(3 * k + 1) * lpc(2) + 2. * lpc(3 * k + 1) * lpc(3 * k + 2) +
619 2 * lpc(1) * lpc(2) + lpc(1) * lpc(3 * k + 2) ) * a ) / 24.;
620 }
621 }
622
623 //Compute main axis of the cross-section
624 double angleChange = 0.;
625 double sum = fabs(Ixx + Iyy);
626 double pi = 3.14159265;
627 if ( ( fabs(Ixx - Iyy) / sum > 1.e-6 ) && fabs(Ixy) / sum > 1.e-6 ) {
628 angleChange = 0.5 * atan(-2 * Ixy / ( Ixx - Iyy ) );
629 } else if ( ( fabs(Ixx - Iyy) / sum < 1.e-6 ) && fabs(Ixy) / sum > 1.e-6 ) {
630 angleChange = pi / 4.;
631 }
632
633 if ( Iyy > Ixx ) {
634 angleChange = angleChange + pi / 2.;
635 }
636
637 //Moment of inertias saved in the element
638
639 this->I1 = ( Ixx + Iyy ) / 2. + sqrt(pow( ( Ixx - Iyy ) / 2., 2. ) + pow(Ixy, 2.) );
640 this->I2 = ( Ixx + Iyy ) / 2. - sqrt(pow( ( Ixx - Iyy ) / 2., 2. ) + pow(Ixy, 2.) );
641
642 this->Ip = I1 + I2;
643
644 //Rotation around normal axis by angleChange
645
646 FloatMatrix rotationChange(3, 3);
647 rotationChange.zero();
648
649 rotationChange.at(1, 1) = 1.;
650 rotationChange.at(2, 2) = cos(angleChange);
651 rotationChange.at(2, 3) = -sin(angleChange);
652
653 rotationChange.at(3, 2) = sin(angleChange);
654 rotationChange.at(3, 3) = cos(angleChange);
655
656 this->localCoordinateSystem.beProductOf(rotationChange, lcs);
657
658 //Calculate the polygon vertices in the new coordinate system
659 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
660 for ( int n = 0; n < 3; n++ ) {
661 help(n) = polygonCoords(3 * k + n);
662 }
663
664 test.beProductOf(this->localCoordinateSystem, help);
665 for ( int n = 0; n < 3; n++ ) {
666 lpc(3 * k + n) = test(n);
667 }
668 }
669
670 //Calculate centroid again in local coordinate system
671 centroid.zero();
672 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
673 if ( k < numberOfPolygonVertices - 1 ) {
674 centroid.at(2) += ( lpc(3 * k + 1) + lpc(3 * ( k + 1 ) + 1) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
675 centroid.at(3) += ( lpc(3 * k + 2) + lpc(3 * ( k + 1 ) + 2) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
676 } else { //Back to zero for n+1
677 centroid.at(2) += ( lpc(3 * k + 1) + lpc(1) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
678 centroid.at(3) += ( lpc(3 * k + 2) + lpc(2) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
679 }
680 }
681
682 centroid.times(1. / ( 6. * this->area ) );
683
684 centroid.at(1) = lpc.at(1); //The first component of all lpcs should be the same
685
686 FloatArray midPointLocal(3);
687 midPointLocal.beProductOf(this->localCoordinateSystem, midPoint);
688
689 //eccentricities stored in the element
690 this->eccS = centroid.at(2) - midPointLocal.at(2);
691 this->eccT = centroid.at(3) - midPointLocal.at(3);
692
693 FloatMatrix transposeLCS;
694 transposeLCS.beTranspositionOf(this->localCoordinateSystem);
695
696 this->globalCentroid.beProductOf(transposeLCS, centroid);
697
698 return;
699}
700
701
702void
703Lattice3d :: saveContext(DataStream &stream, ContextMode mode)
704{
705 LatticeStructuralElement :: saveContext(stream, mode);
706
708
709 if ( ( mode & CM_Definition ) ) {
710 if ( ( iores = polygonCoords.storeYourself(stream) ) != CIO_OK ) {
711 THROW_CIOERR(iores);
712 }
713
714 if ( !stream.write(couplingFlag) ) {
716 }
717
718 if ( ( iores = couplingNumbers.storeYourself(stream) ) != CIO_OK ) {
719 THROW_CIOERR(iores);
720 }
721 }
722}
723
724
725void
726Lattice3d :: restoreContext(DataStream &stream, ContextMode mode)
727{
728 LatticeStructuralElement :: restoreContext(stream, mode);
729
731
732 if ( mode & CM_Definition ) {
733 if ( ( iores = this->polygonCoords.restoreYourself(stream) ) != CIO_OK ) {
734 THROW_CIOERR(iores);
735 }
736
737 if ( !stream.read(this->couplingFlag) ) {
739 }
740
741 if ( ( iores = this->couplingNumbers.restoreYourself(stream) ) != CIO_OK ) {
742 THROW_CIOERR(iores);
743 }
744 }
745}
746
747#ifdef __OOFEG
748
749void
750Lattice3d :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
751{
752 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
753
754 if ( mode == OGC_rawGeometry ) {
755 this->drawRawGeometry(gc, tStep);
756 this->drawRawCrossSections(gc, tStep);
757 } else if ( mode == OGC_deformedGeometry ) {
758 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
759 } else if ( mode == OGC_eigenVectorGeometry ) {
760 this->drawDeformedGeometry(gc, tStep, EigenVector);
761 } else if ( mode == OGC_scalarPlot ) {
762 this->drawScalar(gc, tStep);
763 } else if ( mode == OGC_elemSpecial ) {
764 this->drawSpecial(gc, tStep);
765 } else {
766 OOFEM_ERROR("unsupported mode");
767 }
768}
769
770
771
772
773void Lattice3d :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
774{
775 GraphicObj *go;
776
777 WCRec p [ 2 ]; /* points */
778 if ( !gc.testElementGraphicActivity(this) ) {
779 return;
780 }
781
782 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
783 EASValsSetColor(gc.getElementColor() );
784 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
785
786 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
787 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
788 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
789 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
790 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
791 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
792
793 go = CreateLine3D(p);
794 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
795 EGAttachObject(go, ( EObjectP ) this);
796 EMAddGraphicsToModel(ESIModel(), go);
797}
798
799
800void
801Lattice3d :: drawRawCrossSections(oofegGraphicContext &gc, TimeStep *tStep)
802{
803 GraphicObj *go;
804
805 //Create as many points as we have polygon vertices
806 numberOfPolygonVertices = this->polygonCoords.giveSize() / 3.;
807 WCRec p [ numberOfPolygonVertices ]; /* poin */
808
809 if ( !gc.testElementGraphicActivity(this) ) {
810 return;
811 }
812
813 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
814
815 EASValsSetLayer(OOFEG_RAW_CROSSSECTION_LAYER);
816
817 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
818 p [ i ].x = ( FPNum ) polygonCoords(3 * i);
819 p [ i ].y = ( FPNum ) polygonCoords(3 * i + 1);
820 p [ i ].z = ( FPNum ) polygonCoords(3 * i + 2);
821 }
822
823
824 WCRec pTemp [ 2 ]; /* points */
825 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
826 if ( i < numberOfPolygonVertices - 1 ) {
827 pTemp [ 0 ] = p [ i ];
828 pTemp [ 1 ] = p [ i + 1 ];
829 } else {
830 pTemp [ 0 ] = p [ i ];
831 pTemp [ 1 ] = p [ 0 ];
832 }
833
834 go = CreateLine3D(pTemp);
835 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
836 EGAttachObject(go, ( EObjectP ) this);
837 EMAddGraphicsToModel(ESIModel(), go);
838 }
839}
840
841
842
843void Lattice3d :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
844{
845 GraphicObj *go;
846
847 if ( !gc.testElementGraphicActivity(this) ) {
848 return;
849 }
850
851 double defScale = gc.getDefScale();
852
853 WCRec p [ 2 ]; /* points */
854
855 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
856 EASValsSetColor(gc.getDeformedElementColor() );
857 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
858
859 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
860 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
861 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
862
863 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
864 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
865 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
866
867 go = CreateLine3D(p);
868 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
869 EMAddGraphicsToModel(ESIModel(), go);
870}
871
872#endif
873} // end namespace oofem
#define REGISTER_Element(class)
virtual int read(int *data, std::size_t count)=0
Reads count integer values into array pointed by data.
virtual int write(const int *data, std::size_t count)=0
Writes count integer values from array pointed by data.
double giveCoordinate(int i) const
Definition dofmanager.h:383
Node * giveNode(int i) const
Definition element.h:629
virtual Material * giveMaterial()
Definition element.C:523
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition element.h:1077
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
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 beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void times(double f)
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
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
void drawRawCrossSections(oofegGraphicContext &, TimeStep *tStep)
Definition lattice3d.C:801
FloatArray normal
Definition lattice3d.h:66
FloatArray midPoint
Definition lattice3d.h:66
static ParamKey IPK_Lattice3d_mlength
Definition lattice3d.h:76
IntArray couplingNumbers
Definition lattice3d.h:69
FloatArray polygonCoords
Definition lattice3d.h:62
FloatArray centroid
Definition lattice3d.h:66
void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
Definition lattice3d.C:843
static ParamKey IPK_Lattice3d_pressures
Definition lattice3d.h:75
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
Definition lattice3d.C:86
static ParamKey IPK_Lattice3d_couplingnumber
Definition lattice3d.h:74
double computeVolumeAround(GaussPoint *) override
Definition lattice3d.C:366
static ParamKey IPK_Lattice3d_polycoords
Definition lattice3d.h:72
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
Definition lattice3d.C:773
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
FloatArray pressures
Definition lattice3d.h:70
int giveLocalCoordinateSystem(FloatMatrix &answer) override
Definition lattice3d.C:353
virtual void computeGeometryProperties()
Definition lattice3d.C:418
FloatArray globalCentroid
Definition lattice3d.h:66
static ParamKey IPK_Lattice3d_couplingflag
Definition lattice3d.h:73
int numberOfPolygonVertices
Definition lattice3d.h:63
const FloatArrayF< 6 > & giveOldPlasticLatticeStrain() const
Returns plastic lattice strain.
virtual double giveCrackWidth() const
virtual int giveCrackFlag() const
const FloatArrayF< 6 > & givePlasticLatticeStrain() const
Returns plastic lattice strain.
double giveNormalLatticeStress() const
Gives the last equilibrated normal stress.
#define THROW_CIOERR(e)
#define CM_Definition
Definition contextmode.h:47
#define OOFEM_ERROR(...)
Definition error.h:79
long ContextMode
Definition contextmode.h:43
double sum(const FloatArray &x)
@ CIO_IOERR
General IO error.
#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_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