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