OOFEM 3.0
Loading...
Searching...
No Matches
lattice2dboundary.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 library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (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 GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
38#include "domain.h"
39#include "node.h"
40#include "material.h"
41#include "gausspoint.h"
43#include "floatmatrix.h"
44#include "intarray.h"
45#include "floatarray.h"
46#include "mathfem.h"
47#include "datastream.h"
48#include "contextioerr.h"
49#include "classfactory.h"
50#include "dof.h"
52#include "parametermanager.h"
53#include "paramkey.h"
54
55#ifdef __OOFEG
56 #include "oofeggraphiccontext.h"
57#endif
58
59namespace oofem {
62
63Lattice2dBoundary :: Lattice2dBoundary(int n, Domain *aDomain) : Lattice2d(n, aDomain)
64 // Constructor.
65{
67
68 kappa = -1; // set kappa to undef value (should be always > 0.)
69
70 pitch = 10.; // a dummy value
71}
72
73Lattice2dBoundary :: ~Lattice2dBoundary()
74{
75 // destructor
76}
77
78void
79Lattice2dBoundary :: computeBmatrixAt(GaussPoint *aGaussPoint, FloatMatrix &answer, int li, int ui)
80// Returns the strain matrix of the receiver.
81{
82 double length = this->giveLength();
83 double x1, x2, y1, y2, xp, yp;
84
85 FloatArray projectionComponent(2);
86 giveSwitches(projectionComponent);
87
88
89 xp = this->gpCoords.at(1);
90 yp = this->gpCoords.at(2);
91
92
93 //specimen dimension
94 FloatArray specimenDimension(2);
95 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
96 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
97
98
99 //Coordinates of the nodes
100 x1 = this->giveNode(1)->giveCoordinate(1);
101 y1 = this->giveNode(1)->giveCoordinate(2);
102 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
103 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
104
105
106 double ecc;
107
108 double areaHelp = 0.5 * ( x1 * y2 + x2 * yp + xp * y1 - ( xp * y2 + yp * x1 + x2 * y1 ) );
109
110 ecc = 2 * areaHelp / length;
111
112 // Assemble Bmatrix (used to compute strains and rotations
113 answer.resize(3, 6);
114 answer.zero();
115 answer.at(1, 1) = -1.;
116 answer.at(1, 2) = 0.;
117 answer.at(1, 3) = ecc;
118 answer.at(1, 4) = 1.;
119 answer.at(1, 5) = 0.;
120 answer.at(1, 6) = -ecc;
121
122 answer.at(2, 1) = 0.;
123 answer.at(2, 2) = -1.;
124 answer.at(2, 3) = -length / 2.;
125 answer.at(2, 4) = 0.;
126 answer.at(2, 5) = 1.;
127 answer.at(2, 6) = -length / 2.;
128
129 answer.at(3, 1) = 0.;
130 answer.at(3, 2) = 0.;
131 answer.at(3, 3) = -this->width / sqrt(12.);
132
133 answer.at(3, 4) = 0.;
134 answer.at(3, 5) = 0.;
135 answer.at(3, 6) = this->width / sqrt(12.);
136
137
138 answer.times(1. / length);
139
140 return;
141}
142
143const IntArray
144Lattice2dBoundary :: giveLocation()
145{
146 //Map 2D locations to 3D system
147 int newLocation = 0;
148 if ( this->location == 1 ) {
149 newLocation = 22;
150 } else if ( this->location == 2 ) {
151 newLocation = 25;
152 } else if ( this->location == 3 ) {
153 newLocation = 16;
154 } else if ( this->location == 4 ) {
155 newLocation = 8;
156 } else if ( this->location == 5 ) {
157 newLocation = 5;
158 } else if ( this->location == 6 ) {
159 newLocation = 2;
160 } else if ( this->location == 7 ) {
161 newLocation = 11;
162 } else if ( this->location == 8 ) {
163 newLocation = 19;
164 }
165
166 IntArray tempLocation(2);
167 tempLocation.at(1) = 0;
168 tempLocation.at(2) = newLocation;
169 return tempLocation;
170}
171
172void
173Lattice2dBoundary :: recalculateCoordinates(int nodeNumber, FloatArray &coords) {
174 coords.resize(3);
175 coords.zero();
176 Node *node;
177
178 FloatArray specimenDimension(3);
179 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
180 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
181 specimenDimension.at(3) = this->giveNode(3)->giveCoordinate(3);
182
183 FloatArray projectionComponent(3);
184 projectionComponent.zero();
185
186 node = this->giveNode(2);
187 if ( location != 0 ) {
188 giveSwitches(projectionComponent);
189 }
190
191 for ( int i = 0; i < 3; i++ ) {
192 coords.at(i + 1) = node->giveCoordinate(i + 1) + projectionComponent.at(i + 1) * specimenDimension.at(i + 1);
193 }
194 return;
195}
196
197
198void
199Lattice2dBoundary :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
200// Computes numerically the stiffness matrix of the receiver.
201{
202 double dV;
203 FloatMatrix d, bi, bj, dbj, dij;
204
205 FloatMatrix answerTemp(6, 6);
206 answerTemp.zero();
207
208
209 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), bj);
210 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
211 dV = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
212 dbj.beProductOf(d, bj);
213 answerTemp.plusProductUnsym(bj, dbj, dV);
214
216 answer.zero();
217
218 for ( int m = 1; m <= 6; m++ ) {
219 for ( int k = 1; k <= 6; k++ ) {
220 answer.at(m, k) = answerTemp.at(m, k);
221 }
222 }
223
224
225 //Rotate to global system
226 FloatMatrix R;
227 if ( this->giveRotationMatrix(R) ) {
228 answer.rotatedWith(R);
229 }
230
231 //Convert normal stiffness matrix to the modified one
232
233 FloatArray projectionComponent(2);
234 giveSwitches(projectionComponent);
235
236 answer.at(1, 7) = projectionComponent.at(1) * answer.at(1, 4);
237 answer.at(1, 8) = projectionComponent.at(2) * answer.at(1, 5);
238 answer.at(1, 9) = projectionComponent.at(1) * answer.at(1, 5);
239
240 answer.at(2, 7) = projectionComponent.at(1) * answer.at(2, 4);
241 answer.at(2, 8) = projectionComponent.at(2) * answer.at(2, 5);
242 answer.at(2, 9) = projectionComponent.at(1) * answer.at(2, 5);
243
244 answer.at(3, 7) = projectionComponent.at(1) * answer.at(3, 4);
245 answer.at(3, 8) = projectionComponent.at(2) * answer.at(3, 5);
246 answer.at(3, 9) = projectionComponent.at(1) * answer.at(3, 5);
247
248 answer.at(4, 7) = projectionComponent.at(1) * answer.at(4, 4);
249 answer.at(4, 8) = projectionComponent.at(2) * answer.at(4, 5);
250 answer.at(4, 9) = projectionComponent.at(1) * answer.at(4, 5);
251
252 answer.at(5, 7) = projectionComponent.at(1) * answer.at(5, 4);
253 answer.at(5, 8) = projectionComponent.at(2) * answer.at(5, 5);
254 answer.at(5, 9) = projectionComponent.at(1) * answer.at(5, 5);
255
256 answer.at(6, 7) = projectionComponent.at(1) * answer.at(6, 4);
257 answer.at(6, 8) = projectionComponent.at(2) * answer.at(6, 5);
258 answer.at(6, 9) = projectionComponent.at(1) * answer.at(6, 5);
259
260 answer.at(7, 1) = projectionComponent.at(1) * answer.at(4, 1);
261 answer.at(7, 2) = projectionComponent.at(1) * answer.at(4, 2);
262 answer.at(7, 3) = projectionComponent.at(1) * answer.at(4, 3);
263 answer.at(7, 4) = projectionComponent.at(1) * answer.at(4, 4);
264 answer.at(7, 5) = projectionComponent.at(1) * answer.at(4, 5);
265 answer.at(7, 6) = projectionComponent.at(1) * answer.at(4, 6);
266 answer.at(7, 7) = projectionComponent.at(1) * projectionComponent.at(1) * answer.at(4, 4);
267 answer.at(7, 8) = projectionComponent.at(1) * projectionComponent.at(2) * answer.at(4, 5);
268 answer.at(7, 9) = projectionComponent.at(1) * projectionComponent.at(1) * answer.at(4, 5);
269
270 answer.at(8, 1) = projectionComponent.at(2) * answer.at(5, 1);
271 answer.at(8, 2) = projectionComponent.at(2) * answer.at(5, 2);
272 answer.at(8, 3) = projectionComponent.at(2) * answer.at(5, 3);
273 answer.at(8, 4) = projectionComponent.at(2) * answer.at(5, 4);
274 answer.at(8, 5) = projectionComponent.at(2) * answer.at(5, 5);
275 answer.at(8, 6) = projectionComponent.at(2) * answer.at(5, 6);
276 answer.at(8, 7) = projectionComponent.at(1) * projectionComponent.at(2) * answer.at(5, 4);
277 answer.at(8, 8) = projectionComponent.at(2) * projectionComponent.at(2) * answer.at(5, 5);
278 answer.at(8, 9) = projectionComponent.at(1) * projectionComponent.at(2) * answer.at(5, 5);
279
280 answer.at(9, 1) = projectionComponent.at(1) * answer.at(5, 1);
281 answer.at(9, 2) = projectionComponent.at(1) * answer.at(5, 2);
282 answer.at(9, 3) = projectionComponent.at(1) * answer.at(5, 3);
283 answer.at(9, 4) = projectionComponent.at(1) * answer.at(5, 4);
284 answer.at(9, 5) = projectionComponent.at(1) * answer.at(5, 5);
285 answer.at(9, 6) = projectionComponent.at(1) * answer.at(5, 6);
286 answer.at(9, 7) = projectionComponent.at(1) * projectionComponent.at(1) * answer.at(5, 4);
287 answer.at(9, 8) = projectionComponent.at(1) * projectionComponent.at(2) * answer.at(5, 5);
288 answer.at(9, 9) = projectionComponent.at(1) * projectionComponent.at(1) * answer.at(5, 5);
289 //Rotate back to local system
290 FloatMatrix Rtranspose;
291 Rtranspose.beTranspositionOf(R);
292 answer.rotatedWith(Rtranspose);
293}
294
295void
296Lattice2dBoundary :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN)
297// Computes the vector containing the strains at the Gauss point gp of
298// the receiver, at time step stepN. The nature of these strains depends
299// on the element's type.
300{
301 FloatMatrix b;
302 FloatArray u;
303
304 //Compute strain vector
305 //Get the 9 components of the displacement vector of this element
306 this->computeVectorOf(VM_Total, stepN, u);
307 this->computeBmatrixAt(gp, b);
308
309
310 //Compute projection vector
311 FloatArray projectionComponent(2);
312 giveSwitches(projectionComponent);
313
314 FloatMatrix rotationMatrix;
315 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
316 u.rotatedWith(rotationMatrix, 't');
317 }
318
319 //General expressions for the corrected displacements
320 u.at(4) = u.at(4) + projectionComponent.at(1) * u.at(7);
321 u.at(5) = u.at(5) + projectionComponent.at(2) * u.at(8) + projectionComponent.at(1) * u.at(9);
322
323
324 if ( this->computeGtoLRotationMatrix(rotationMatrix) ) {
325 u.rotatedWith(rotationMatrix, 'n');
326 }
327
328
329 //define a temp displacement vector
330 FloatArray uTemp(6);
331
332 for ( int i = 1; i <= 6; i++ ) {
333 uTemp.at(i) = u.at(i);
334 }
335
336 answer.beProductOf(b, uTemp);
337}
338
339
340bool
341Lattice2dBoundary :: computeGtoLRotationMatrix(FloatMatrix &answer)
342{
343 double sine, cosine;
344 answer.resize(9, 9);
345 answer.zero();
346
347 sine = sin(this->givePitch() );
348 cosine = cos(pitch);
349 answer.at(1, 1) = cosine;
350 answer.at(1, 2) = sine;
351 answer.at(2, 1) = -sine;
352 answer.at(2, 2) = cosine;
353 answer.at(3, 3) = 1.;
354 answer.at(4, 4) = cosine;
355 answer.at(4, 5) = sine;
356 answer.at(5, 4) = -sine;
357 answer.at(5, 5) = cosine;
358 answer.at(6, 6) = 1.;
359 answer.at(7, 7) = 1.;
360 answer.at(8, 8) = 1.;
361 answer.at(9, 9) = 1.;
362
363 return 1;
364}
365
366
367double
368Lattice2dBoundary :: computeVolumeAround(GaussPoint *aGaussPoint)
369{
370 double weight = aGaussPoint->giveWeight();
371 return weight * 0.5 * this->giveLength() * this->width * this->thickness;
372}
373
374
375void
376Lattice2dBoundary :: giveDofManDofIDMask(int inode, IntArray &answer) const {
377 if ( inode == 3 ) {
378 answer = { E_xx, E_yy, G_xy };
379 } else {
380 answer = { D_u, D_v, R_w };
381 }
382}
383
384
385
386double Lattice2dBoundary :: giveLength()
387// Returns the length of the receiver.
388{
389 //Compute strain vector
390 FloatArray projectionComponent(2);
391 giveSwitches(projectionComponent);
392
393 //specimen dimension
394 FloatArray specimenDimension(2);
395 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
396 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
397
398 double dx, dy;
399 Node *nodeA, *nodeB;
400
401 nodeA = this->giveNode(1);
402 nodeB = this->giveNode(2);
403
404 dx = nodeB->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1) - nodeA->giveCoordinate(1);
405 dy = nodeB->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2) - nodeA->giveCoordinate(2);
406
407 return sqrt(dx * dx + dy * dy);
408}
409
410
411double Lattice2dBoundary :: givePitch()
412// Returns the pitch of the receiver.
413{
414 double xA, xB, yA, yB;
415 Node *nodeA, *nodeB;
416
417 //specimen dimension
418 FloatArray specimenDimension(2);
419 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
420 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
421
422 FloatArray projectionComponent(2);
423 giveSwitches(projectionComponent);
424
425 if ( pitch == 10. ) { // 10. : dummy initialization value
426 nodeA = this->giveNode(1);
427 nodeB = this->giveNode(2);
428 xA = nodeA->giveCoordinate(1);
429 xB = nodeB->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
430 yA = nodeA->giveCoordinate(2);
431 yB = nodeB->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
432 pitch = atan2(yB - yA, xB - xA);
433 }
434 return pitch;
435}
436
437
438
439void
440Lattice2dBoundary :: initializeFrom(InputRecord &ir, int priority)
441{
442 ParameterManager &ppm = this->giveDomain()->elementPPM;
443 // first call parent
444 Lattice2d :: initializeFrom(ir, priority);
446}
447
448void
449Lattice2dBoundary :: postInitialize()
450{
451 ParameterManager &ppm = this->giveDomain()->elementPPM;
452 // first call parent
453 Lattice2d :: postInitialize();
455 this->computeGaussPoints();
456}
457
458
459
460void
461Lattice2dBoundary :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
462{
463 Material *mat = this->giveMaterial();
464
465 FloatMatrix b, bt, A, R, GNT;
466 FloatArray bs, TotalStressVector, u, strain;
467 double dV;
468
469
470 // This function can be quite costly to do inside the loops when one has many slave dofs.
471 this->computeVectorOf(VM_Total, tStep, u);
472 // subtract initial displacements, if defined
473 if ( initialDisplacements ) {
475 }
476
477 answer.resize(9);
478 answer.zero();
479
480 this->computeBmatrixAt(integrationRulesArray [ 0 ]->getIntegrationPoint(0), b);
481
482 bt.beTranspositionOf(b);
483 // TotalStressVector = gp->giveStressVector() ;
484 if ( useUpdatedGpRecord == 1 ) {
485 TotalStressVector = ( ( LatticeMaterialStatus * ) mat->giveStatus(integrationRulesArray [ 0 ]->getIntegrationPoint(0) ) )
486 ->giveLatticeStress();
487 } else
488 if ( !this->isActivated(tStep) ) {
489 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(integrationRulesArray [ 0 ]->getIntegrationPoint(0)->giveMaterialMode() ) );
490 strain.zero();
491 }
492 this->computeStrainVector(strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
493
494 // strain.beProductOf(b, u);
495
496 this->computeStressVector(TotalStressVector, strain, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
497 //
498 // compute nodal representation of internal forces using f = B^T*Sigma dV
499 //
500 dV = this->computeVolumeAround(integrationRulesArray [ 0 ]->getIntegrationPoint(0) );
501 bs.beProductOf(bt, TotalStressVector);
502 bs.times(dV);
503
504 for ( int m = 1; m <= 6; m++ ) {
505 answer.at(m) = bs.at(m);
506 }
507
508 //Rotate to global system
509
510 if ( this->giveRotationMatrix(R) ) {
511 answer.rotatedWith(R, 't');
512 }
513 //Compute strain vector
514
515 FloatArray projectionComponent(2);
516 giveSwitches(projectionComponent);
517
518 answer.at(7) = projectionComponent.at(1) * answer.at(4);
519 answer.at(8) = projectionComponent.at(2) * answer.at(5);
520 answer.at(9) = projectionComponent.at(1) * answer.at(5);
521
522 //Rotate to local system
523 if ( this->giveRotationMatrix(R) ) {
524 answer.rotatedWith(R, 'n');
525 }
526
527 return;
528}
529
530
531void
532Lattice2dBoundary :: giveSwitches(FloatArray &answer) {
533 if ( this->location == 1 ) {
534 answer(0) = 1;
535 answer(1) = 0;
536 } else if ( this->location == 2 ) {
537 answer(0) = 1;
538 answer(1) = 1;
539 } else if ( this->location == 3 ) {
540 answer(0) = 0;
541 answer(1) = 1;
542 } else if ( this->location == 4 ) {
543 answer(0) = -1;
544 answer(1) = 1;
545 } else if ( this->location == 5 ) {
546 answer(0) = -1;
547 answer(1) = 0;
548 } else if ( this->location == 6 ) {
549 answer(0) = -1;
550 answer(1) = -1;
551 } else if ( this->location == 7 ) {
552 answer(0) = 0;
553 answer(1) = -1;
554 } else if ( this->location == 8 ) {
555 answer(0) = 1;
556 answer(1) = -1;
557 }
558
559 return;
560}
561
562
563void Lattice2dBoundary :: saveContext(DataStream &stream, ContextMode mode)
564{
565 Lattice2d :: saveContext(stream, mode);
566
567
568 if ( ( mode & CM_Definition ) ) {
569 if ( !stream.write(location) ) {
571 }
572 }
573}
574
575
576void
577Lattice2dBoundary :: restoreContext(DataStream &stream, ContextMode mode)
578{
579 Lattice2d :: restoreContext(stream, mode);
580
581 if ( mode & CM_Definition ) {
582 if ( !stream.read(location) ) {
584 }
585 }
586}
587
588
589#ifdef __OOFEG
590
591void
592Lattice2dBoundary :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
593{
594 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
595
596 if ( mode == OGC_rawGeometry ) {
597 this->drawRawGeometry(gc, tStep);
598 this->drawRawCrossSections(gc, tStep);
599 } else if ( mode == OGC_deformedGeometry ) {
600 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
601 } else if ( mode == OGC_elemSpecial ) {
602 this->drawSpecial(gc, tStep);
603 } else {
604 OOFEM_ERROR("unsupported mode");
605 }
606}
607
608
609void Lattice2dBoundary :: drawRawCrossSections(oofegGraphicContext &gc, TimeStep *tStep)
610{
611 GraphicObj *go;
612
613 // if (!go) { // create new one
614 WCRec p [ 2 ]; /* poin */
615 if ( !gc.testElementGraphicActivity(this) ) {
616 return;
617 }
618
619 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
620 EASValsSetColor(gc.getCrossSectionColor() );
621 EASValsSetLayer(OOFEG_RAW_CROSSSECTION_LAYER);
622
623 //The cs definition is not needed anymore
624
627
628 FloatArray projectionComponent(2);
629 giveSwitches(projectionComponent);
630
631 FloatArray specimenDimension(2);
632
633 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
634 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
635
636 // double width = cs->give(WIDTH);
637
638 double x1, y1, x2, y2;
639 x1 = this->giveNode(1)->giveCoordinate(1);
640 y1 = this->giveNode(1)->giveCoordinate(2);
641
642 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
643 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
644
645 //Compute normal and shear direction
646 FloatArray normalDirection;
647 FloatArray shearDirection;
648 normalDirection.resize(2);
649 normalDirection.zero();
650 shearDirection.resize(2);
651 shearDirection.zero();
652 normalDirection.at(1) = x2 - x1;
653 normalDirection.at(2) = y2 - y1;
654 normalDirection.normalize();
655 if ( normalDirection.at(2) == 0. ) {
656 shearDirection.at(1) = 0.;
657 shearDirection.at(2) = 1.;
658 } else {
659 shearDirection.at(1) = 1.0;
660 shearDirection.at(2) =
661 -normalDirection.at(1) / normalDirection.at(2);
662 }
663 shearDirection.normalize();
664
665 p [ 0 ].x = ( FPNum ) this->gpCoords.at(1) - shearDirection.at(1) * this->width / 2.;
666 p [ 0 ].y = ( FPNum ) this->gpCoords.at(2) - shearDirection.at(2) * this->width / 2.;
667 p [ 0 ].z = 0.;
668
669 p [ 1 ].x = ( FPNum ) this->gpCoords.at(1) + shearDirection.at(1) * this->width / 2.;
670 p [ 1 ].y = ( FPNum ) this->gpCoords.at(2) + shearDirection.at(2) * this->width / 2.;
671 p [ 1 ].z = 0.;
672
673 go = CreateLine3D(p);
674 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
675 EGAttachObject(go, ( EObjectP ) this);
676 EMAddGraphicsToModel(ESIModel(), go);
677}
678
679
680void Lattice2dBoundary :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
681{
682 GraphicObj *go;
683
684 // if (!go) { // create new one
685 WCRec p [ 2 ]; /* poin */
686 if ( !gc.testElementGraphicActivity(this) ) {
687 return;
688 }
689
690 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
691 EASValsSetColor(gc.getElementColor() );
692 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
693
694 FloatArray projectionComponent(2);
695 giveSwitches(projectionComponent);
696
697
698 FloatArray specimenDimension(2);
699
700 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
701 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
702
703 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
704 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
705 p [ 0 ].z = 0.;
706
707 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
708 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
709 p [ 1 ].z = 0.;
710
711 go = CreateLine3D(p);
712 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
713 EGAttachObject(go, ( EObjectP ) this);
714 EMAddGraphicsToModel(ESIModel(), go);
715}
716
717
718
719void Lattice2dBoundary :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
720{
721 GraphicObj *go;
722
723 if ( !gc.testElementGraphicActivity(this) ) {
724 return;
725 }
726
727 double defScale = gc.getDefScale();
728 // if (!go) { // create new one
729 WCRec p [ 2 ]; /* poin */
730 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
731 EASValsSetColor(gc.getDeformedElementColor() );
732 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
733
734
735 FloatArray projectionComponent(2);
736 giveSwitches(projectionComponent);
737
738 //specimen dimension
739 FloatArray specimenDimension(2);
740 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
741 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
742
743 double x1, y1, x2, y2;
744 x1 = this->giveNode(1)->giveCoordinate(1);
745 y1 = this->giveNode(1)->giveCoordinate(2);
746 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
747 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
748
749 //Node One
750 FloatArray dispOne(3);
751 dispOne.at(1) = this->giveNode(1)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
752 dispOne.at(2) = this->giveNode(1)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
753 dispOne.at(3) = this->giveNode(1)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
754
755 //Node Two
756 FloatArray dispTwo(3);
757 dispTwo.at(1) = this->giveNode(2)->giveDofWithID(D_u)->giveUnknown(VM_Total, tStep);
758 dispTwo.at(2) = this->giveNode(2)->giveDofWithID(D_v)->giveUnknown(VM_Total, tStep);
759 dispTwo.at(3) = this->giveNode(2)->giveDofWithID(R_w)->giveUnknown(VM_Total, tStep);
760
761 //Node Three
762 FloatArray dispThree(3);
763 dispThree.at(1) = this->giveNode(3)->giveDofWithID(E_xx)->giveUnknown(VM_Total, tStep);
764 dispThree.at(2) = this->giveNode(3)->giveDofWithID(E_yy)->giveUnknown(VM_Total, tStep);
765 dispThree.at(3) = this->giveNode(3)->giveDofWithID(G_xy)->giveUnknown(VM_Total, tStep);
766
767 //Modify dispOne and dispTwo
768 dispTwo.at(1) = dispTwo.at(1) + projectionComponent.at(1) * dispThree.at(1);
769 dispTwo.at(2) = dispTwo.at(2) + projectionComponent.at(2) * dispThree.at(2) + projectionComponent.at(1) * dispThree.at(3);
770
771 p [ 0 ].x = ( FPNum ) x1 + defScale * dispOne.at(1);
772 p [ 0 ].y = ( FPNum ) y1 + defScale * dispOne.at(2);
773 p [ 0 ].z = 0.;
774
775 p [ 1 ].x = ( FPNum ) x2 + defScale * dispTwo.at(1);
776 p [ 1 ].y = ( FPNum ) y2 + defScale * dispTwo.at(2);
777 p [ 1 ].z = 0.;
778
779 go = CreateLine3D(p);
780 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
781 EMAddGraphicsToModel(ESIModel(), go);
782}
783
784
785
786
787
788
789void
790Lattice2dBoundary :: drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
791{
792 WCRec l [ 2 ];
793 GraphicObj *tr;
795 GaussPoint *gp;
796 FloatArray crackStatuses, cf;
797
798 if ( !gc.testElementGraphicActivity(this) ) {
799 return;
800 }
801 if ( gc.giveIntVarType() == IST_CrackState ) {
802 gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
803 mat->giveIPValue(crackStatuses, gp, IST_CrackStatuses, tStep); // printf("crackFlag= %e\n",crackStatuses(0));
804 if ( crackStatuses(0) == 1. || crackStatuses(0) == 2. ) {
807
808 //specimen dimension
809 FloatArray specimenDimension(2);
810 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
811 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
812
813
814 FloatArray projectionComponent(2);
815 giveSwitches(projectionComponent);
816
817
818 double x1, y1, x2, y2;
819 x1 = this->giveNode(1)->giveCoordinate(1);
820 y1 = this->giveNode(1)->giveCoordinate(2);
821 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
822 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
823
824 //Compute normal and shear direction
825 FloatArray normalDirection;
826 FloatArray shearDirection;
827 normalDirection.resize(2);
828 normalDirection.zero();
829 shearDirection.resize(2);
830 shearDirection.zero();
831 normalDirection.at(1) = x2 - x1;
832 normalDirection.at(2) = y2 - y1;
833 normalDirection.normalize();
834 if ( normalDirection.at(2) == 0. ) {
835 shearDirection.at(1) = 0.;
836 shearDirection.at(2) = 1.;
837 } else {
838 shearDirection.at(1) = 1.0;
839 shearDirection.at(2) =
840 -normalDirection.at(1) / normalDirection.at(2);
841 }
842 shearDirection.normalize();
843
844 l [ 0 ].x = ( FPNum ) this->gpCoords.at(1) - shearDirection.at(1) * this->width / 2.;
845 l [ 0 ].y = ( FPNum ) this->gpCoords.at(2) - shearDirection.at(2) * this->width / 2.;
846 l [ 0 ].z = 0.;
847 l [ 1 ].x = ( FPNum ) this->gpCoords.at(1) + shearDirection.at(1) * this->width / 2.;
848 ;
849 l [ 1 ].y = ( FPNum ) this->gpCoords.at(2) + shearDirection.at(2) * this->width / 2.;
850 l [ 1 ].z = 0.;
851
852 EASValsSetLayer(OOFEG_CRACK_PATTERN_LAYER);
853 EASValsSetLineWidth(OOFEG_CRACK_PATTERN_WIDTH);
854 if ( ( crackStatuses(0) == 1. ) ) {
855 EASValsSetColor(gc.getActiveCrackColor() );
856 } else if ( crackStatuses(0) == 2. ) {
857 EASValsSetColor(gc.getCrackPatternColor() );
858 }
859 tr = CreateLine3D(l);
860 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
861 EMAddGraphicsToModel(ESIModel(), tr);
862 }
863 }
864}
865
866void
867Lattice2dBoundary :: giveCrossSectionCoordinates(FloatArray &coords)
868{
871 FloatArray projectionComponent(2);
872 giveSwitches(projectionComponent);
873
874 FloatArray specimenDimension(2);
875
876 specimenDimension.at(1) = this->giveNode(3)->giveCoordinate(1);
877 specimenDimension.at(2) = this->giveNode(3)->giveCoordinate(2);
878
879 double x1, y1, x2, y2;
880 x1 = this->giveNode(1)->giveCoordinate(1);
881 y1 = this->giveNode(1)->giveCoordinate(2);
882
883 x2 = this->giveNode(2)->giveCoordinate(1) + projectionComponent.at(1) * specimenDimension.at(1);
884 y2 = this->giveNode(2)->giveCoordinate(2) + projectionComponent.at(2) * specimenDimension.at(2);
885
886 //Compute normal and shear direction
887 FloatArray normalDirection;
888 FloatArray shearDirection;
889 normalDirection.resize(2);
890 normalDirection.zero();
891 shearDirection.resize(2);
892 shearDirection.zero();
893 normalDirection.at(1) = x2 - x1;
894 normalDirection.at(2) = y2 - y1;
895 normalDirection.normalize();
896 if ( normalDirection.at(2) == 0. ) {
897 shearDirection.at(1) = 0.;
898 shearDirection.at(2) = 1.;
899 } else {
900 shearDirection.at(1) = 1.0;
901 shearDirection.at(2) =
902 -normalDirection.at(1) / normalDirection.at(2);
903 }
904 shearDirection.normalize();
905
906 coords.resize(6);
907 coords.at(1) = this->gpCoords.at(1) - shearDirection.at(1) * this->width / 2.;
908 coords.at(2) = this->gpCoords.at(2) - shearDirection.at(2) * this->width / 2.;
909 coords.at(3) = 0.;
910 coords.at(4) = this->gpCoords.at(1) + shearDirection.at(1) * this->width / 2.;
911 coords.at(5) = this->gpCoords.at(2) + shearDirection.at(2) * this->width / 2.;
912 coords.at(6) = 0.;
913
914 return;
915}
916
917#endif
918} // 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 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
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 plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
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
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
int & at(std::size_t i)
Definition intarray.h:104
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *stepN) override
void giveSwitches(FloatArray &answer)
static ParamKey IPK_Lattice2dBoundary_location
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep) override
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
bool computeGtoLRotationMatrix(FloatMatrix &) override
void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
void drawRawCrossSections(oofegGraphicContext &, TimeStep *tStep)
int computeNumberOfDofs() override
double computeVolumeAround(GaussPoint *) override
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition lattice2d.C:168
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition lattice2d.C:174
void computeGaussPoints() override
Definition lattice2d.C:196
FloatArray gpCoords
Definition lattice2d.h:61
void giveGpCoordinates(FloatArray &coords) override
Definition lattice2d.C:365
Lattice2d(int n, Domain *d)
Definition lattice2d.C:69
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.
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
#define THROW_CIOERR(e)
#define CM_Definition
Definition contextmode.h:47
#define OOFEM_ERROR(...)
Definition error.h:79
long ContextMode
Definition contextmode.h:43
@ CIO_IOERR
General IO error.
#define OOFEG_RAW_CROSSSECTION_LAYER
#define OOFEG_CRACK_PATTERN_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 OOFEG_CRACK_PATTERN_WIDTH
#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