OOFEM 3.0
Loading...
Searching...
No Matches
latticelink3d.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 "contextioerr.h"
48#include "datastream.h"
49#include "classfactory.h"
52#include "parametermanager.h"
53#include "paramkey.h"
54
55#ifdef __OOFEG
56 #include "oofeggraphiccontext.h"
57#endif
58
59namespace oofem {
65
66LatticeLink3d :: LatticeLink3d(int n, Domain *aDomain) : LatticeStructuralElement(n, aDomain)
67{
69 geometryFlag = 0;
70}
71
72LatticeLink3d :: ~LatticeLink3d()
73{}
74
75double
76LatticeLink3d :: computeVolumeAround(GaussPoint *aGaussPoint)
77{
78 //Returns artifical volume (bond area times bond length) so that general parts of post processing work
79 return pow(this->bondLength, 2.) * this->bondDiameter * M_PI;
80}
81
82
83void
84LatticeLink3d :: computeBmatrixAt(GaussPoint *aGaussPoint, FloatMatrix &answer, int li, int ui)
85// Returns the strain matrix of the receiver.
86{
87 if ( geometryFlag == 0 ) {
89 }
90
91 //Assemble Bmatrix based on three rigid arm components
92 //rigid.at(1) (tangential), rigid.at(2) (lateral), rigid.at(3) (lateral)
93 answer.resize(6, 12);
94 answer.zero();
95
96 //Normal displacement jump in x-direction
97 //First node
98 answer.at(1, 1) = -1.;
99 answer.at(1, 5) = -this->rigid.at(3);
100 answer.at(1, 6) = this->rigid.at(2);
101 //Second node
102 answer.at(1, 7) = 1.;
103
104 //Shear displacement jump in y-plane
105 //first node
106 answer.at(2, 2) = -1.;
107 answer.at(2, 4) = this->rigid.at(3);
108 answer.at(2, 6) = -this->rigid.at(1);
109 //Second node
110 answer.at(2, 8) = 1.;
111
112 //Shear displacement jump in z-plane
113 //first node
114 answer.at(3, 3) = -1.;
115 answer.at(3, 4) = -this->rigid.at(2);
116 answer.at(3, 5) = this->rigid.at(1);
117 //Second node
118 answer.at(3, 9) = 1.;
119
120 //Rotation around x-axis
121 //First node
122 answer.at(4, 4) = -1.;
123 //Second node
124 answer.at(4, 10) = 1.;
125
126 //Rotation around y-axis
127 //First node
128 answer.at(5, 5) = -1.;
129 //Second node
130 answer.at(5, 11) = 1.;
131
132 //Rotation around z-axis
133 //First node
134 answer.at(6, 6) = -1.;
135 //Second node
136 answer.at(6, 12) = 1.;
137
138 return;
139}
140
141void
142LatticeLink3d :: giveGPCoordinates(FloatArray &coords)
143{
144 if ( geometryFlag == 0 ) {
146 }
147 coords.resize(3);
148 coords = this->globalCentroid;
149 return;
150}
151
152
153void
154LatticeLink3d :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
155 TimeStep *tStep)
156// Computes numerically the stiffness matrix of the receiver.
157{
158 FloatMatrix d, b, bt, db;
159 FloatArray u, strain;
160
161 // This function can be quite costly to do inside the loops when one has many slave dofs.
162 this->computeVectorOf(VM_Total, tStep, u);
163 // subtract initial displacements, if defined
164 if ( initialDisplacements ) {
166 }
167
168 // zero answer will resize accordingly when adding first contribution
169 answer.clear();
170
171
172 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
173 this->computeBmatrixAt(gp, b);
174 bt.beTranspositionOf(b);
175
176 if ( !this->isActivated(tStep) ) {
177 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
178 strain.zero();
179 }
180 strain.beProductOf(b, u);
181
182
183 answer.resize(12, 12);
184 answer.zero();
185
186 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
187
188
189 db.beProductOf(d, b);
190 answer.beProductOf(bt, db);
191
192 //Introduce integration of bond strength
193 double area = this->computeVolumeAround(gp) / this->giveLength();
194 answer.times(area);
195
196 return;
197}
198
199void LatticeLink3d :: computeGaussPoints()
200// Sets up the array of Gauss Points of the receiver.
201{
202 integrationRulesArray.resize(1);
203 integrationRulesArray [ 0 ].reset(new GaussIntegrationRule(1, this, 1, 3) );
204 integrationRulesArray [ 0 ]->SetUpPointsOnLine(1, _3dLattice);
205}
206
207
208
209
210bool
211LatticeLink3d :: computeGtoLRotationMatrix(FloatMatrix &answer)
212{
213 FloatMatrix lcs;
214 int i, j;
215
216 answer.resize(12, 12);
217 answer.zero();
218
219 this->giveLocalCoordinateSystem(lcs);
220 for ( i = 1; i <= 3; i++ ) {
221 for ( j = 1; j <= 3; j++ ) {
222 answer.at(i, j) = lcs.at(i, j);
223 answer.at(i + 3, j + 3) = lcs.at(i, j);
224 answer.at(i + 6, j + 6) = lcs.at(i, j);
225 answer.at(i + 9, j + 9) = lcs.at(i, j);
226 }
227 }
228
229 return 1;
230}
231
232
233int
234LatticeLink3d :: giveLocalCoordinateSystem(FloatMatrix &answer)
235{
236 if ( geometryFlag == 0 ) {
238 }
239
240 answer = this->localCoordinateSystem;
241
242 return 1;
243}
244
245double
246LatticeLink3d :: giveBondLength()
247{
248 return this->bondLength;
249}
250
251double
252LatticeLink3d :: giveBondEndLength()
253{
254 return this->bondEndLength;
255}
256
257
258double
259LatticeLink3d :: giveBondDiameter()
260{
261 return this->bondDiameter;
262}
263
264
265
266
267void
268LatticeLink3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
269{
270 answer = {
271 D_u, D_v, D_w, R_u, R_v, R_w
272 };
273}
274
275void
276LatticeLink3d :: initializeFrom(InputRecord &ir, int priority)
277{
278 ParameterManager &ppm = this->giveDomain()->elementPPM;
279 // first call parent
280 LatticeStructuralElement :: initializeFrom(ir, priority);
281
286
287}
288
289void
290LatticeLink3d :: postInitialize()
291{
292 ParameterManager &ppm = this->giveDomain()->elementPPM;
293 LatticeStructuralElement :: postInitialize();
298}
299
300int
301LatticeLink3d :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
302{
303 if ( geometryFlag == 0 ) {
305 }
306
307 answer.resize(3);
308 answer = this->globalCentroid;
309
310 return 1;
311}
312
313
314void
315LatticeLink3d :: computeGeometryProperties()
316{
317 //coordinates of the two nodes
318 Node *nodeA, *nodeB;
319 FloatArray coordsA(3), coordsB(3);
320
321 //Order of nodes. However, might not matter.
322 //Reinforcement node
323 nodeA = this->giveNode(1);
324 //Lattice node
325 nodeB = this->giveNode(2);
326
327 //Calculate components of distance from reinforcement node to lattice node.
328
329 for ( int i = 0; i < 3; i++ ) {
330 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1);
331 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1);
332 }
333
334 FloatArray rigidGlobal(3);
335
336 //Calculate normal vector
337 for ( int i = 0; i < 3; i++ ) {
338 rigidGlobal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
339 }
340
341 //Construct an initial temporary local coordinate system
342 FloatArray normal(3), s(3), t(3);
343
344 //Calculate normal vector
345 normal = this->directionVector;
346 normal.normalize();
347
348 //Construct two perpendicular axis so that n is normal to the plane which they create
349 //Check, if one of the components of the normal-direction is zero
350 if ( normal.at(1) == 0 ) {
351 s.at(1) = 0.;
352 s.at(2) = normal.at(3);
353 s.at(3) = -normal.at(2);
354 } else if ( normal.at(2) == 0 ) {
355 s.at(1) = normal.at(3);
356 s.at(2) = 0.;
357 s.at(3) = -normal.at(1);
358 } else {
359 s.at(1) = normal.at(2);
360 s.at(2) = -normal.at(1);
361 s.at(3) = 0.;
362 }
363
364 s.normalize();
365
366 t.beVectorProductOf(normal, s);
367 t.normalize();
368
369 //Set up rotation matrix
370 FloatMatrix lcs(3, 3);
371
372 this->localCoordinateSystem.resize(3, 3);
373 this->localCoordinateSystem.zero();
374
375 for ( int i = 1; i <= 3; i++ ) {
376 this->localCoordinateSystem.at(1, i) = normal.at(i);
377 this->localCoordinateSystem.at(2, i) = s.at(i);
378 this->localCoordinateSystem.at(3, i) = t.at(i);
379 }
380
381 // Rotate rigidarm vector into local coordinate system
382 this->rigid.beProductOf(localCoordinateSystem, rigidGlobal);
383
384 this->globalCentroid.resize(3);
385 for ( int i = 1; i <= 3; i++ ) {
386 this->globalCentroid.at(i) = nodeA->giveCoordinate(i);
387 ;
388 }
389
390 this->geometryFlag = 1;
391
392 return;
393}
394void
395LatticeLink3d :: saveContext(DataStream &stream, ContextMode mode)
396{
397 LatticeStructuralElement :: saveContext(stream, mode);
398}
399
400
401void
402LatticeLink3d :: restoreContext(DataStream &stream, ContextMode mode)
403{
404 LatticeStructuralElement :: restoreContext(stream, mode);
405}
406
407
408void
409LatticeLink3d :: giveInternalForcesVector(FloatArray &answer,
410 TimeStep *tStep, int useUpdatedGpRecord)
411{
412 FloatMatrix b, bt;
413 FloatArray u, stress(6), strain(6);
414
415 // This function can be quite costly to do inside the loops when one has many slave dofs.
416 this->computeVectorOf(VM_Total, tStep, u);
417 // subtract initial displacements, if defined
418 if ( initialDisplacements ) {
420 }
421
422 // zero answer will resize accordingly when adding first contribution
423 answer.clear();
424
425 for ( GaussPoint *gp: * this->giveDefaultIntegrationRulePtr() ) {
426 LatticeMaterialStatus *matStat = static_cast< LatticeMaterialStatus * >( gp->giveMaterialStatus() );
427 this->computeBmatrixAt(gp, b);
428 bt.beTranspositionOf(b);
429
430 if ( useUpdatedGpRecord == 1 ) {
431 stress = matStat->giveLatticeStress();
432 } else {
433 if ( !this->isActivated(tStep) ) {
434 strain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
435 strain.zero();
436 }
437 strain.beProductOf(b, u);
438 this->computeStressVector(stress, strain, gp, tStep);
439 }
440
441 // updates gp stress and strain record acording to current
442 // increment of displacement
443 if ( stress.giveSize() == 0 ) {
444 break;
445 }
446
447 // now every gauss point has real stress vector
448 // compute nodal representation of internal forces using f = B^T*Sigma dV
449 // double dV = this->computeVolumeAround(gp);
450 // double dV = 1.;
451 if ( stress.giveSize() == 6 ) {
452 // It may happen that e.g. plane strain is computed
453 // using the default 3D implementation. If so,
454 // the stress needs to be reduced.
455 // (Note that no reduction will take place if
456 // the simulation is actually 3D.)
457 FloatArray stressTemp;
458 StructuralMaterial :: giveReducedSymVectorForm(stressTemp, stress, gp->giveMaterialMode() );
459
460 answer.beProductOf(bt, stressTemp);
461 } else {
462 answer.beProductOf(bt, stress);
463 }
464
465 //Introduce integration of bond strength
466 double area = this->computeVolumeAround(gp) / this->giveLength();
467 answer.times(area);
468 }
469
470
471 // if inactive update state, but no contribution to global system
472 if ( !this->isActivated(tStep) ) {
473 answer.zero();
474 return;
475 }
476}
477
478
479double
480LatticeLink3d :: giveLength()
481{
482 //returns the bond length, not the length between the two nodes
483 return this->bondLength;
484}
485
486void
487LatticeLink3d :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
488{
489 answer = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->give3dStiffnessMatrix(rMode, gp, tStep);
490}
491
492void
493LatticeLink3d :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
494{
495 answer = static_cast< LatticeCrossSection * >( this->giveCrossSection() )->giveLatticeStress3d(strain, gp, tStep);
496}
497
498
499
500#ifdef __OOFEG
501
502void
503LatticeLink3d :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
504{
505 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
506
507 if ( mode == OGC_rawGeometry ) {
508 this->drawRawGeometry(gc, tStep);
509 } else if ( mode == OGC_deformedGeometry ) {
510 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
511 } else if ( mode == OGC_eigenVectorGeometry ) {
512 this->drawDeformedGeometry(gc, tStep, EigenVector);
513 } else if ( mode == OGC_scalarPlot ) {
514 this->drawScalar(gc, tStep);
515 } else if ( mode == OGC_elemSpecial ) {
516 this->drawSpecial(gc, tStep);
517 } else {
518 OOFEM_ERROR("unsupported mode");
519 }
520}
521
522
523
524void LatticeLink3d :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
525{
526 GraphicObj *go;
527
528 WCRec p [ 2 ]; /* points */
529 if ( !gc.testElementGraphicActivity(this) ) {
530 return;
531 }
532
533 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
534 EASValsSetColor(gc.getElementColor() );
535 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
536
537 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
538 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
539 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
540 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
541 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
542 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
543
544 go = CreateLine3D(p);
545 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
546 EGAttachObject(go, ( EObjectP ) this);
547 EMAddGraphicsToModel(ESIModel(), go);
548}
549
550void LatticeLink3d :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
551{
552 GraphicObj *go;
553
554 if ( !gc.testElementGraphicActivity(this) ) {
555 return;
556 }
557
558 double defScale = gc.getDefScale();
559
560 WCRec p [ 2 ]; /* points */
561
562 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
563 EASValsSetColor(gc.getDeformedElementColor() );
564 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
565
566 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
567 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
568 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
569
570 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
571 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
572 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
573
574 go = CreateLine3D(p);
575 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
576 EMAddGraphicsToModel(ESIModel(), go);
577}
578
579#endif
580} // 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 isActivated(TimeStep *tStep)
Definition element.C:838
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition element.h:1077
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 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 zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
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 times(double s)
Definition floatarray.C:834
void times(double f)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
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
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition gausspoint.h:190
int giveLocalCoordinateSystem(FloatMatrix &answer) override
static ParamKey IPK_LatticeLink3d_l_end
double giveLength() override
static ParamKey IPK_LatticeLink3d_length
void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
FloatArray directionVector
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
FloatArray globalCentroid
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
double computeVolumeAround(GaussPoint *aGaussPoint) override
static ParamKey IPK_LatticeLink3d_dirvector
static ParamKey IPK_LatticeLink3d_diameter
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
FloatMatrix localCoordinateSystem
virtual void computeGeometryProperties()
const FloatArrayF< 6 > & giveLatticeStress() const
Returns lattice stress.
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
#define OOFEM_ERROR(...)
Definition error.h:79
#define M_PI
Definition mathfem.h:52
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