OOFEM 3.0
Loading...
Searching...
No Matches
latticelink3dboundary.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
64LatticeLink3dBoundary :: LatticeLink3dBoundary(int n, Domain *aDomain) : LatticeLink3d(n, aDomain), location(2)
65{
67 geometryFlag = 0;
68}
69
70LatticeLink3dBoundary :: ~LatticeLink3dBoundary()
71{}
72
73
74void
75LatticeLink3dBoundary :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
76 TimeStep *tStep)
77// Computes numerically the stiffness matrix of the receiver.
78{
79 FloatMatrix d, bi, bj, dbj, dij, bjt;
80 FloatMatrix t(12, 18), tt;
81 FloatMatrix answerTemp(12, 12), answerHelp, ttk(18, 12);
82 bool matStiffSymmFlag = this->giveCrossSection()->isCharacteristicMtrxSymmetric(rMode);
83 answerTemp.zero();
84 answerHelp.zero();
85 t.zero();
86
87 if ( geometryFlag == 0 ) {
89 }
90
91
92 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), bj);
93 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
94
95 dbj.beProductOf(d, bj);
96 bjt.beTranspositionOf(bj);
97 answerTemp.beProductOf(bjt, dbj);
98
100 answer.zero();
101
103 answerHelp.zero();
104
105 for ( int m = 1; m <= 12; m++ ) {
106 for ( int k = 1; k <= 12; k++ ) {
107 answer.at(m, k) = answerTemp.at(m, k);
108 }
109 }
110
111 if ( matStiffSymmFlag ) {
112 answer.symmetrized();
113 }
114
115 // Rotate to global system add the projections and rotate back to local system.
116
117 FloatMatrix R;
118 if ( this->giveRotationMatrix(R) ) {
119 answer.rotatedWith(R);
120 }
121
122 for ( int m = 1; m <= 12; m++ ) {
123 for ( int k = 1; k <= 12; k++ ) {
124 answerTemp.at(m, k) = answer.at(m, k);
125 }
126 }
127
128 IntArray projectionComponentNodeOne(3);
129 projectionComponentNodeOne.zero();
130 if ( location.at(1) != 0 ) {
131 giveSwitches(projectionComponentNodeOne, location.at(1) );
132 }
133
134 IntArray projectionComponentNodeTwo(3);
135 projectionComponentNodeTwo.zero();
136 if ( location.at(2) != 0 ) {
137 giveSwitches(projectionComponentNodeTwo, location.at(2) );
138 }
139
140 for ( int k = 1; k <= 12; k++ ) {
141 t.at(k, k) = 1.;
142 }
143 t.at(1, 13) = projectionComponentNodeOne.at(1);
144 t.at(2, 14) = projectionComponentNodeOne.at(2);
145 t.at(3, 15) = projectionComponentNodeOne.at(3);
146
147
148 t.at(2, 16) = projectionComponentNodeOne.at(3);
149 t.at(1, 17) = projectionComponentNodeOne.at(3);
150 t.at(1, 18) = projectionComponentNodeOne.at(2);
151
152
153 t.at(7, 13) = projectionComponentNodeTwo.at(1);
154 t.at(8, 14) = projectionComponentNodeTwo.at(2);
155 t.at(9, 15) = projectionComponentNodeTwo.at(3);
156
157 t.at(8, 16) = projectionComponentNodeTwo.at(3);
158 t.at(7, 17) = projectionComponentNodeTwo.at(3);
159 t.at(7, 18) = projectionComponentNodeTwo.at(2);
160
161
162 tt.beTranspositionOf(t);
163
164 ttk.beProductOf(tt, answerTemp);
165 answer.beProductOf(ttk, t);
166
167 //Rotate back to local system
168 FloatMatrix Rtranspose;
169 Rtranspose.beTranspositionOf(R);
170 answer.rotatedWith(Rtranspose);
171
172 return;
173}
174
175void
176LatticeLink3dBoundary :: recalculateCoordinates(int nodeNumber, FloatArray &coords) {
177 coords.resize(3);
178 coords.zero();
179 Node *node;
180
181 FloatArray specimenDimension(3);
182 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
183 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
184 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
185
186 IntArray projectionComponent(3);
187 projectionComponent.zero();
188
189 if ( nodeNumber == 1 ) {
190 node = this->giveNode(1);
191 if ( location.at(1) != 0 ) {
192 giveSwitches(projectionComponent, location.at(1) );
193 }
194 } else if ( nodeNumber == 2 ) {
195 node = this->giveNode(2);
196 if ( location.at(2) != 0 ) {
197 giveSwitches(projectionComponent, location.at(2) );
198 }
199 } else {
200 OOFEM_ERROR("wrong element used in the vtk module");
201 }
202
203 for ( int i = 0; i < 3; i++ ) {
204 coords.at(i + 1) = node->giveCoordinate(i + 1) + projectionComponent.at(i + 1) * specimenDimension.at(i + 1);
205 }
206 return;
207}
208
209void
210LatticeLink3dBoundary :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN)
211// Computes the vector containing the strains at the Gauss point gp of
212// the receiver, at time step stepN. The nature of these strains depends
213// on the element's type.
214{
215 FloatMatrix b;
216 FloatArray u;
217
218
219 //Compute strain vector
220 //Get the 18 components of the displacement vector of this element
221 this->computeVectorOf(VM_Total, stepN, u);
222 this->computeBmatrixAt(gp, b);
223
224
225 // subtract initial displacements, if defined
226 if ( initialDisplacements ) {
228 }
229
230 //Compute projection vector
231 IntArray projectionComponentNodeOne(3);
232 projectionComponentNodeOne.zero();
233 if ( location.at(1) != 0 ) {
234 giveSwitches(projectionComponentNodeOne, location.at(1) );
235 }
236
237 IntArray projectionComponentNodeTwo(3);
238 projectionComponentNodeTwo.zero();
239
240 if ( location.at(2) != 0 ) {
241 giveSwitches(projectionComponentNodeTwo, location.at(2) );
242 }
243
244 FloatMatrix rotationMatrix;
245 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
246 u.rotatedWith(rotationMatrix, 't');
247 }
248
249 //General expressions for the corrected displacements. Rotations at nodes are not influenced. Only translations
250 u.at(1) = u.at(1) + projectionComponentNodeOne.at(1) * u.at(13) + projectionComponentNodeOne.at(2) * u.at(18) + projectionComponentNodeOne.at(3) * u.at(17);
251 u.at(2) = u.at(2) + projectionComponentNodeOne.at(2) * u.at(14) + projectionComponentNodeOne.at(3) * u.at(16);
252 u.at(3) = u.at(3) + projectionComponentNodeOne.at(3) * u.at(15);
253
254 u.at(7) = u.at(7) + projectionComponentNodeTwo.at(1) * u.at(13) + projectionComponentNodeTwo.at(2) * u.at(18) + projectionComponentNodeTwo.at(3) * u.at(17);
255 u.at(8) = u.at(8) + projectionComponentNodeTwo.at(2) * u.at(14) + projectionComponentNodeTwo.at(3) * u.at(16);
256 u.at(9) = u.at(9) + projectionComponentNodeTwo.at(3) * u.at(15);
257
258
259 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
260 u.rotatedWith(rotationMatrix, 'n');
261 }
262 //define a temp displacement vector
263 FloatArray uTemp(12);
264
265 for ( int i = 1; i <= 12; i++ ) {
266 uTemp.at(i) = u.at(i);
267 }
268
269 answer.beProductOf(b, uTemp);
270}
271
272bool
273LatticeLink3dBoundary :: computeGtoLRotationMatrix(FloatMatrix &answer)
274{
275 FloatMatrix lcs;
276 int i, j;
277
278 answer.resize(18, 18);
279 answer.zero();
280
281 this->giveLocalCoordinateSystem(lcs);
282 for ( i = 1; i <= 3; i++ ) {
283 for ( j = 1; j <= 3; j++ ) {
284 answer.at(i, j) = lcs.at(i, j);
285 answer.at(i + 3, j + 3) = lcs.at(i, j);
286 answer.at(i + 6, j + 6) = lcs.at(i, j);
287 answer.at(i + 9, j + 9) = lcs.at(i, j);
288 }
289 }
290 answer.at(13, 13) = 1.;
291 answer.at(14, 14) = 1.;
292 answer.at(15, 15) = 1.;
293 answer.at(16, 16) = 1.;
294 answer.at(17, 17) = 1.;
295 answer.at(18, 18) = 1.;
296
297 return 1;
298}
299
300int
301LatticeLink3dBoundary :: giveLocalCoordinateSystem(FloatMatrix &answer)
302{
303 if ( geometryFlag == 0 ) {
305 }
306
307 answer = this->localCoordinateSystem;
308
309 return 1;
310}
311
312
313void
314LatticeLink3dBoundary :: giveDofManDofIDMask(int inode, IntArray &answer) const
315{
316 if ( inode == 3 ) {
317 answer = { E_xx, E_yy, E_zz, G_yz, G_xz, G_xy };
318 } else {
319 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
320 }
321}
322
323void
324LatticeLink3dBoundary :: initializeFrom(InputRecord &ir, int priority)
325{
326 ParameterManager &ppm = this->giveDomain()->elementPPM;
327 LatticeLink3d :: initializeFrom(ir, priority);
328
330}
331
332void
333LatticeLink3dBoundary :: postInitialize()
334{
335 ParameterManager &ppm = this->giveDomain()->elementPPM;
336 LatticeLink3d :: postInitialize();
338}
339
340
341
342void
343LatticeLink3dBoundary :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
344{
345 Material *mat = this->giveMaterial();
346
347 FloatMatrix b, bt, A, R, GNT;
348 FloatArray bs, TotalStressVector, u, strain;
349
350 this->computeVectorOf(VM_Total, tStep, u);
351
352 // subtract initial displacements, if defined
353 if ( initialDisplacements ) {
355 }
356
357 answer.resize(18);
358 answer.zero();
359
360 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), b);
361
362 bt.beTranspositionOf(b);
363 if ( useUpdatedGpRecord == 1 ) {
364 TotalStressVector = ( ( StructuralMaterialStatus * ) mat->giveStatus(integrationRulesArray [ 0 ]->getIntegrationPoint(0) ) )
365 ->giveStressVector();
366 } else
367 if ( !this->isActivated(tStep) ) {
368 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(integrationRulesArray [ 0 ]->getIntegrationPoint(0)->giveMaterialMode() ) );
369 strain.zero();
370 }
371 this->computeStrainVector(strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
372
373 this->computeStressVector(TotalStressVector, strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
374
375 bs.beProductOf(bt, TotalStressVector);
376
377 for ( int m = 1; m <= 12; m++ ) {
378 answer.at(m) = bs.at(m);
379 }
380
381 if ( this->giveRotationMatrix(R) ) {
382 answer.rotatedWith(R, 't');
383 }
384
385 IntArray projectionComponentNodeOne(3);
386 projectionComponentNodeOne.zero();
387 if ( location.at(1) != 0 ) {
388 giveSwitches(projectionComponentNodeOne, location.at(1) );
389 }
390
391 IntArray projectionComponentNodeTwo(3);
392 projectionComponentNodeTwo.zero();
393
394 if ( location.at(2) != 0 ) {
395 giveSwitches(projectionComponentNodeTwo, location.at(2) );
396 }
397
398 //Normal stresses
399 answer.at(13) = projectionComponentNodeOne.at(1) * answer.at(1) + projectionComponentNodeTwo.at(1) * answer.at(7);
400 answer.at(14) = projectionComponentNodeOne.at(2) * answer.at(2) + projectionComponentNodeTwo.at(2) * answer.at(8);
401 answer.at(15) = projectionComponentNodeOne.at(3) * answer.at(3) + projectionComponentNodeTwo.at(3) * answer.at(9);
402
403 //Shear stresses
404 answer.at(16) = projectionComponentNodeOne.at(3) * answer.at(2) + projectionComponentNodeTwo.at(3) * answer.at(8);
405 answer.at(17) = projectionComponentNodeOne.at(3) * answer.at(1) + projectionComponentNodeTwo.at(3) * answer.at(7);
406 answer.at(18) = projectionComponentNodeOne.at(2) * answer.at(1) + projectionComponentNodeTwo.at(2) * answer.at(7);
407
408 //Rotate to local system
409 if ( this->giveRotationMatrix(R) ) {
410 answer.rotatedWith(R, 'n');
411 }
412
413 return;
414}
415
416void
417LatticeLink3dBoundary :: giveSwitches(IntArray &answer, int location) {
418 int counter = 1;
419 for ( int x = -1; x < 2; x++ ) {
420 for ( int y = -1; y < 2; y++ ) {
421 for ( int z = -1; z < 2; z++ ) {
422 if ( !( z == 0 && y == 0 && x == 0 ) ) {
423 if ( counter == location ) {
424 answer(0) = x;
425 answer(1) = y;
426 answer(2) = z;
427 }
428 counter++;
429 }
430 }
431 }
432 }
433 return;
434}
435
436void
437LatticeLink3dBoundary :: computeGeometryProperties()
438{
439 Node *nodeA, *nodeB;
440
441 FloatArray coordsA(3);
442 FloatArray coordsB(3);
443
444 nodeA = this->giveNode(1);
445 nodeB = this->giveNode(2);
446
447 FloatArray specimenDimension(3);
448 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
449 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
450 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
451
452 IntArray projectionComponentNodeOne(3);
453 projectionComponentNodeOne.zero();
454 if ( location.at(1) != 0 ) {
455 giveSwitches(projectionComponentNodeOne, location.at(1) );
456 }
457
458 IntArray projectionComponentNodeTwo(3);
459 projectionComponentNodeTwo.zero();
460 if ( location.at(2) != 0 ) {
461 giveSwitches(projectionComponentNodeTwo, location.at(2) );
462 }
463
464 for ( int i = 0; i < 3; i++ ) {
465 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1) + projectionComponentNodeOne.at(i + 1) * specimenDimension.at(i + 1);
466 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1) + projectionComponentNodeTwo.at(i + 1) * specimenDimension.at(i + 1);
467 }
468
469 //From here on the normal functions should be valid
470
471 FloatArray rigidGlobal(3);
472
473 //Calculate normal vector
474 for ( int i = 0; i < 3; i++ ) {
475 rigidGlobal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
476 }
477
478 //Construct an initial temporary local coordinate system
479 FloatArray normal(3), s(3), t(3);
480
481 //Calculate normal vector
482 normal = this->directionVector;
483 normal.normalize();
484
485 //Construct two perpendicular axis so that n is normal to the plane which they create
486 //Check, if one of the components of the normal-direction is zero
487 if ( normal.at(1) == 0 ) {
488 s.at(1) = 0.;
489 s.at(2) = normal.at(3);
490 s.at(3) = -normal.at(2);
491 } else if ( normal.at(2) == 0 ) {
492 s.at(1) = normal.at(3);
493 s.at(2) = 0.;
494 s.at(3) = -normal.at(1);
495 } else {
496 s.at(1) = normal.at(2);
497 s.at(2) = -normal.at(1);
498 s.at(3) = 0.;
499 }
500
501 s.normalize();
502
503 t.beVectorProductOf(normal, s);
504 t.normalize();
505
506 //Set up rotation matrix
507 FloatMatrix lcs(3, 3);
508
509 this->localCoordinateSystem.resize(3, 3);
510 this->localCoordinateSystem.zero();
511
512 for ( int i = 1; i <= 3; i++ ) {
513 this->localCoordinateSystem.at(1, i) = normal.at(i);
514 this->localCoordinateSystem.at(2, i) = s.at(i);
515 this->localCoordinateSystem.at(3, i) = t.at(i);
516 }
517
518 // Rotate rigidarm vector into local coordinate system
519
520 this->rigid.beProductOf(localCoordinateSystem, rigidGlobal);
521
522 this->globalCentroid.resize(3);
523 for ( int i = 1; i <= 3; i++ ) {
524 this->globalCentroid.at(i) = coordsA.at(i);
525 }
526
527 this->geometryFlag = 1;
528
529 return;
530}
531
532
533void
534LatticeLink3dBoundary :: saveContext(DataStream &stream, ContextMode mode)
535{
536 LatticeLink3d :: saveContext(stream, mode);
537
539
540 if ( ( mode & CM_Definition ) ) {
541 if ( ( iores = location.storeYourself(stream) ) != CIO_OK ) {
542 THROW_CIOERR(iores);
543 }
544 }
545}
546
547
548
549void
550LatticeLink3dBoundary :: restoreContext(DataStream &stream, ContextMode mode)
551{
553
554 LatticeLink3d :: restoreContext(stream, mode);
555
556 if ( mode & CM_Definition ) {
557 if ( ( iores = this->location.restoreYourself(stream) ) != CIO_OK ) {
558 THROW_CIOERR(iores);
559 }
560 }
561}
562
563
564
565#ifdef __OOFEG
566
567void
568LatticeLink3dBoundary :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
569{
570 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
571
572 if ( mode == OGC_rawGeometry ) {
573 this->drawRawGeometry(gc, tStep);
574 } else if ( mode == OGC_deformedGeometry ) {
575 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
576 } else if ( mode == OGC_elemSpecial ) {
577 this->drawSpecial(gc, tStep);
578 } else {
579 OOFEM_ERROR("drawYourself : unsupported mode");
580 }
581}
582
583
584
585
586void LatticeLink3dBoundary :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
587{
588 GraphicObj *go;
589
590 WCRec p [ 2 ]; /* points */
591 if ( !gc.testElementGraphicActivity(this) ) {
592 return;
593 }
594
595 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
596 EASValsSetColor(gc.getActiveCrackColor() );
597 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
598
599
600 FloatArray specimenDimension(3);
601 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
602 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
603 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
604
605 //Compute projection vector
606 IntArray projectionComponentNodeOne(3);
607 projectionComponentNodeOne.zero();
608 if ( location.at(1) != 0 ) {
609 giveSwitches(projectionComponentNodeOne, location.at(1) );
610 }
611
612 IntArray projectionComponentNodeTwo(3);
613 projectionComponentNodeTwo.zero();
614
615 if ( location.at(2) != 0 ) {
616 giveSwitches(projectionComponentNodeTwo, location.at(2) );
617 }
618
619 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
620 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
621 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
622 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
623 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
624 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
625
626 go = CreateLine3D(p);
627 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
628 EGAttachObject(go, ( EObjectP ) this);
629 EMAddGraphicsToModel(ESIModel(), go);
630}
631
632void LatticeLink3dBoundary :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
633{
634 GraphicObj *go;
635
636 if ( !gc.testElementGraphicActivity(this) ) {
637 return;
638 }
639
640 double defScale = gc.getDefScale();
641
642 WCRec p [ 2 ]; /* points */
643
644 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
645 EASValsSetColor(gc.getDeformedElementColor() );
646 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
647
648 FloatArray specimenDimension(3);
649 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
650 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
651 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
652
653 FloatArray dispOne(6);
654 dispOne.at(1) = this->giveNode(1)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
655 dispOne.at(2) = this->giveNode(1)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
656 dispOne.at(3) = this->giveNode(1)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
657 dispOne.at(4) = this->giveNode(1)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
658 dispOne.at(5) = this->giveNode(1)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
659 dispOne.at(6) = this->giveNode(1)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
660
661 FloatArray dispTwo(6);
662 dispTwo.at(1) = this->giveNode(2)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
663 dispTwo.at(2) = this->giveNode(2)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
664 dispTwo.at(3) = this->giveNode(2)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
665 dispTwo.at(4) = this->giveNode(2)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
666 dispTwo.at(5) = this->giveNode(2)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
667 dispTwo.at(6) = this->giveNode(2)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
668
669 FloatArray dispThree(6);
670 dispThree.at(1) = this->giveNode(3)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
671 dispThree.at(2) = this->giveNode(3)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
672 dispThree.at(3) = this->giveNode(3)->giveDofWithID(D_w)->giveUnknown(VM_Total, tStep);
673 dispThree.at(4) = this->giveNode(3)->giveDofWithID(R_u)->giveUnknown(VM_Total, tStep);
674 dispThree.at(5) = this->giveNode(3)->giveDofWithID(R_v)->giveUnknown(VM_Total, tStep);
675 dispThree.at(6) = this->giveNode(3)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
676
677 IntArray projectionComponentNodeOne(3);
678 projectionComponentNodeOne.zero();
679 if ( location.at(1) != 0 ) {
680 giveSwitches(projectionComponentNodeOne, location.at(1) );
681 }
682
683 IntArray projectionComponentNodeTwo(3);
684 projectionComponentNodeTwo.zero();
685 if ( location.at(2) != 0 ) {
686 giveSwitches(projectionComponentNodeTwo, location.at(2) );
687 }
688
689
690 //Modify dispOne and dispTwo
691 dispOne.at(1) = dispOne.at(1) + projectionComponentNodeOne.at(1) * dispThree.at(1) + projectionComponentNodeOne.at(2) * dispThree.at(4) + projectionComponentNodeOne.at(3) * dispThree.at(5);
692 dispOne.at(2) = dispOne.at(2) + projectionComponentNodeOne.at(2) * dispThree.at(2) + projectionComponentNodeOne.at(3) * dispThree.at(6);
693 dispOne.at(3) = dispOne.at(3) + projectionComponentNodeOne.at(3) * dispThree.at(3);
694
695
696 dispTwo.at(1) = dispTwo.at(1) + projectionComponentNodeTwo.at(1) * dispThree.at(1) + projectionComponentNodeTwo.at(2) * dispThree.at(4) + projectionComponentNodeTwo.at(3) * dispThree.at(5);
697 dispTwo.at(2) = dispTwo.at(2) + projectionComponentNodeTwo.at(2) * dispThree.at(2) + projectionComponentNodeTwo.at(3) * dispThree.at(6);
698 dispTwo.at(3) = dispTwo.at(3) + projectionComponentNodeTwo.at(3) * dispThree.at(3);
699
700 double x1, y1, z1, x2, y2, z2;
701 x1 = this->giveNode(1)->giveCoordinate(1) + projectionComponentNodeOne.at(1) * specimenDimension.at(1);
702 y1 = this->giveNode(1)->giveCoordinate(2) + projectionComponentNodeOne.at(2) * specimenDimension.at(2);
703 z1 = this->giveNode(1)->giveCoordinate(3) + projectionComponentNodeOne.at(3) * specimenDimension.at(3);
704
705 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponentNodeTwo.at(1) * specimenDimension.at(1);
706 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponentNodeTwo.at(2) * specimenDimension.at(2);
707 z2 = this->giveNode(2)->giveCoordinate(3) + projectionComponentNodeTwo.at(3) * specimenDimension.at(3);
708
709 p [ 0 ].x = ( FPNum ) x1 + defScale * dispOne.at(1);
710 p [ 0 ].y = ( FPNum ) y1 + defScale * dispOne.at(2);
711 p [ 0 ].z = ( FPNum ) z1 + defScale * dispOne.at(3);
712
713 p [ 1 ].x = ( FPNum ) x2 + defScale * dispTwo.at(1);
714 p [ 1 ].y = ( FPNum ) y2 + defScale * dispTwo.at(2);
715 p [ 1 ].z = ( FPNum ) z2 + defScale * dispTwo.at(3);
716
717 go = CreateLine3D(p);
718 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
719 EMAddGraphicsToModel(ESIModel(), go);
720}
721
722#endif
723} // 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 beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
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 drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
int giveLocalCoordinateSystem(FloatMatrix &answer) override
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN) override
void giveSwitches(IntArray &answer, int location)
static ParamKey IPK_LatticeLink3dBoundary_location
bool computeGtoLRotationMatrix(FloatMatrix &) override
void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
FloatArray directionVector
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
FloatArray globalCentroid
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
FloatMatrix localCoordinateSystem
LatticeLink3d(int n, Domain *)
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
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