OOFEM 3.0
Loading...
Searching...
No Matches
bondlink3dboundary.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 {
62
63BondLink3dBoundary :: BondLink3dBoundary(int n, Domain *aDomain) : BondLink3d(n, aDomain), location(2)
64{
66 geometryFlag = 0;
67}
68
69BondLink3dBoundary :: ~BondLink3dBoundary()
70{}
71
72void
73BondLink3dBoundary :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode,
74 TimeStep *tStep)
75// Computes numerically the stiffness matrix of the receiver.
76{
77 FloatMatrix d, b, bt, db;
78 FloatArray u, slip;
79
80 this->computeVectorOf(VM_Total, tStep, u);
81
82 // subtract initial displacements, if defined
85 }
86
87 answer.clear();
88
89 GaussPoint *gp = this->giveDefaultIntegrationRulePtr()->getIntegrationPoint(0);
90 this->computeBmatrixAt(gp, b);
92
93 if ( !this->isActivated(tStep) ) {
94 slip.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
95 slip.zero();
96 }
97 slip.beProductOf(b, u);
98
99 answer.resize(9, 9);
100 answer.zero();
101
102 this->computeConstitutiveMatrixAt(d, rMode, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
103
104 db.beProductOf(d, b);
105 answer.beProductOf(bt, db);
106
107 //Introduce integration of bond strength
108 double area = this->computeVolumeAround(gp) / this->giveLength();
109 answer.times(area);
110
111
112 return;
113}
114
115void
116BondLink3dBoundary :: computeTransformationMatrix(FloatMatrix &answer, TimeStep *tStep)
117{
118 FloatArray unitCellSize;
119 unitCellSize.resize(3);
120 unitCellSize.at(1) = this->giveNode(3)->giveCoordinate(1);
121 unitCellSize.at(2) = this->giveNode(3)->giveCoordinate(2);
122 unitCellSize.at(3) = this->giveNode(3)->giveCoordinate(3);
123
124 IntArray switches1, switches2;
125 this->giveSwitches(switches1, this->location.at(1) );
126 this->giveSwitches(switches2, this->location.at(2) );
127
128 FloatMatrix k1, k2;
129 k1.resize(6, 9);
130 k2.resize(6, 9);
131
132 for ( int i = 1; i <= 3; i++ ) {
133 k1.at(i, 3 * i - 2) = unitCellSize.at(1) * switches1.at(1);
134 k1.at(i, 3 * i - 1) = unitCellSize.at(2) * switches1.at(2);
135 k1.at(i, 3 * i) = unitCellSize.at(3) * switches1.at(3);
136 }
137
138 for ( int i = 1; i <= 3; i++ ) {
139 k2.at(i, 3 * i - 2) = unitCellSize.at(1) * switches2.at(1);
140 k2.at(i, 3 * i - 1) = unitCellSize.at(2) * switches2.at(2);
141 k2.at(i, 3 * i) = unitCellSize.at(3) * switches2.at(3);
142 }
143
144 answer.resize(12, 12);
145 answer.beUnitMatrix();
146 answer.resizeWithData(12, 21);
147
148 answer.assemble(k1, { 1, 2, 3, 4, 5, 6 }, { 13, 14, 15, 16, 17, 18, 19, 20, 21 });
149 answer.assemble(k2, { 7, 8, 9, 10, 11, 12 }, { 13, 14, 15, 16, 17, 18, 19, 20, 21 });
150}
151
152
153bool
154BondLink3dBoundary :: computeGtoLRotationMatrix(FloatMatrix &answer)
155{
156 FloatMatrix lcs;
157 int i, j;
158
159 answer.resize(9, 9);
160 answer.zero();
161
162 this->giveLocalCoordinateSystem(lcs);
163 for ( i = 1; i <= 3; i++ ) {
164 for ( j = 1; j <= 3; j++ ) {
165 answer.at(i, j) = lcs.at(i, j);
166 answer.at(i + 3, j + 3) = lcs.at(i, j);
167 answer.at(i + 6, j + 6) = lcs.at(i, j);
168 // answer.at(i + 9, j + 9) = lcs.at(i, j);
169 }
170 }
171
172 return 1;
173}
174
175int
176BondLink3dBoundary :: giveLocalCoordinateSystem(FloatMatrix &answer)
177{
178 if ( geometryFlag == 0 ) {
180 }
181
182 answer = this->localCoordinateSystem;
183
184 return 1;
185}
186
187
188
189void
190BondLink3dBoundary :: giveDofManDofIDMask(int inode, IntArray &answer) const
191{
192 if ( inode == 1 ) {
193 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
194 } else if ( inode == 2 ) {
195 answer = { D_u, D_v, D_w };
196 } else if ( inode == 3 ) {
197 answer = { E_xx, E_xy, E_xz, E_yx, E_yy, E_yz, E_zx, E_zy, E_zz };
198 }
199}
200
201void
202BondLink3dBoundary :: initializeFrom(InputRecord &ir, int priority)
203{
204 // first call parent
205 StructuralElement :: initializeFrom(ir, priority);
206 ParameterManager &ppm = giveDomain()->elementPPM;
208}
209
210
211
212void
213BondLink3dBoundary :: computeGeometryProperties()
214{
215 //coordinates of the two nodes
216 Node *nodeA, *nodeB;
217 FloatArray coordsA(3), coordsB(3);
218
219 //Order of nodes. Important, because continuum node does not have rotational DOFs.
220 //Beam node
221 nodeA = this->giveNode(1);
222 //Continuum node
223 nodeB = this->giveNode(2);
224
225
226 //Calculate components of distance from continuum node to lattice node.
227 for ( int i = 0; i < 3; i++ ) {
228 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1);
229 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1);
230 }
231
232 FloatArray rigidGlobal(3);
233
234 //Calculate normal vector
235 for ( int i = 0; i < 3; i++ ) {
236 rigidGlobal.at(i + 1) = coordsA.at(i + 1) - coordsB.at(i + 1);
237 }
238
239 //Construct an initial temporary local coordinate system
240 FloatArray normal(3), s(3), t(3);
241
242 //Calculate normal vector
243 normal = this->directionVector;
244 normal.normalize();
245
246 //Construct two perpendicular axis so that n is normal to the plane which they create
247 //Check, if one of the components of the normal-direction is zero
248 if ( normal.at(1) == 0 ) {
249 s.at(1) = 0.;
250 s.at(2) = normal.at(3);
251 s.at(3) = -normal.at(2);
252 } else if ( normal.at(2) == 0 ) {
253 s.at(1) = normal.at(3);
254 s.at(2) = 0.;
255 s.at(3) = -normal.at(1);
256 } else {
257 s.at(1) = normal.at(2);
258 s.at(2) = -normal.at(1);
259 s.at(3) = 0.;
260 }
261
262 s.normalize();
263
264 t.beVectorProductOf(normal, s);
265 t.normalize();
266
267 //Set up rotation matrix
268 FloatMatrix lcs(3, 3);
269
270 this->localCoordinateSystem.resize(3, 3);
271 this->localCoordinateSystem.zero();
272
273 for ( int i = 1; i <= 3; i++ ) {
274 this->localCoordinateSystem.at(1, i) = normal.at(i);
275 this->localCoordinateSystem.at(2, i) = s.at(i);
276 this->localCoordinateSystem.at(3, i) = t.at(i);
277 }
278
279 // Rotate rigidarm vector into local coordinate system
280
281 this->rigid.beProductOf(localCoordinateSystem, rigidGlobal);
282
283
284 this->globalCentroid.resize(3);
285 for ( int i = 1; i <= 3; i++ ) {
286 this->globalCentroid.at(i) = nodeB->giveCoordinate(i);
287 ;
288 }
289
290 this->geometryFlag = 1;
291
292 return;
293}
294
295void
296BondLink3dBoundary :: saveContext(DataStream &stream, ContextMode mode)
297{
298 StructuralElement :: saveContext(stream, mode);
299}
300
301
302void
303BondLink3dBoundary :: restoreContext(DataStream &stream, ContextMode mode)
304{
305 StructuralElement :: restoreContext(stream, mode);
306}
307
308
309void
310BondLink3dBoundary :: giveInternalForcesVector(FloatArray &answer,
311 TimeStep *tStep, int useUpdatedGpRecord)
312{
313 FloatMatrix b, bt;
314 FloatArray u, stress(6), slip(6);
315
316 // This function can be quite costly to do inside the loops when one has many slave dofs.
317 this->computeVectorOf(VM_Total, tStep, u);
318 // subtract initial displacements, if defined
319 if ( initialDisplacements ) {
321 }
322
323 // zero answer will resize accordingly when adding first contribution
324 answer.clear();
325
326 for ( GaussPoint *gp: * this->giveDefaultIntegrationRulePtr() ) {
327 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
328 this->computeBmatrixAt(gp, b);
329 bt.beTranspositionOf(b);
330
331 if ( useUpdatedGpRecord == 1 ) {
332 stress = matStat->giveStressVector();
333 } else {
334 if ( !this->isActivated(tStep) ) {
335 slip.resize(6);
336 slip.zero();
337 }
338 slip.beProductOf(b, u);
339 this->computeStressVector(stress, slip, gp, tStep);
340 }
341
342 answer.beProductOf(bt, stress);
343
344 //Introduce integration of bond strength
345 double area = this->computeVolumeAround(gp) / this->giveLength();
346 answer.times(area);
347 }
348
349 // if inactive update state, but no contribution to global system
350 if ( !this->isActivated(tStep) ) {
351 answer.zero();
352 return;
353 }
354}
355
356
357
358void
359BondLink3dBoundary :: giveSwitches(IntArray &answer, int location) {
360 int counter = 1;
361 for ( int x = -1; x < 2; x++ ) {
362 for ( int y = -1; y < 2; y++ ) {
363 for ( int z = -1; z < 2; z++ ) {
364 if ( !( z == 0 && y == 0 && x == 0 ) ) {
365 if ( counter == location ) {
366 answer(0) = x;
367 answer(1) = y;
368 answer(2) = z;
369 }
370 counter++;
371 }
372 }
373 }
374 }
375 return;
376}
377
378
379void
380BondLink3dBoundary :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
381{
382 this->giveStructuralCrossSection()->giveCharMaterialStiffnessMatrix(answer, rMode, gp, tStep);
383}
384
385void
386BondLink3dBoundary :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
387{
388 answer = this->giveStructuralCrossSection()->giveRealStresses(strain, gp, tStep);
389}
390
391
392
393#ifdef __OOFEG
394
395void
396BondLink3dBoundary :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
397{
398 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
399
400 if ( mode == OGC_rawGeometry ) {
401 this->drawRawGeometry(gc, tStep);
402 } else if ( mode == OGC_deformedGeometry ) {
403 this->drawDeformedGeometry(gc, tStep, DisplacementVector);
404 } else if ( mode == OGC_eigenVectorGeometry ) {
405 this->drawDeformedGeometry(gc, tStep, EigenVector);
406 } else if ( mode == OGC_scalarPlot ) {
407 this->drawScalar(gc, tStep);
408 } else if ( mode == OGC_elemSpecial ) {
409 this->drawSpecial(gc, tStep);
410 } else {
411 OOFEM_ERROR("unsupported mode");
412 }
413}
414
415
416
417void BondLink3dBoundary :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
418{
419 GraphicObj *go;
420
421 WCRec p [ 2 ]; /* points */
422 if ( !gc.testElementGraphicActivity(this) ) {
423 return;
424 }
425
426 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
427 EASValsSetColor(gc.getElementColor() );
428 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
429
430 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
431 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
432 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
433 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
434 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
435 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
436
437 go = CreateLine3D(p);
438 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
439 EGAttachObject(go, ( EObjectP ) this);
440 EMAddGraphicsToModel(ESIModel(), go);
441}
442
443void BondLink3dBoundary :: drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type)
444{
445 GraphicObj *go;
446
447 if ( !gc.testElementGraphicActivity(this) ) {
448 return;
449 }
450
451 double defScale = gc.getDefScale();
452
453 WCRec p [ 2 ]; /* points */
454
455 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
456 EASValsSetColor(gc.getDeformedElementColor() );
457 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
458
459 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
460 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
461 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
462
463 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
464 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
465 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
466
467 go = CreateLine3D(p);
468 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
469 EMAddGraphicsToModel(ESIModel(), go);
470}
471
472#endif
473} // end namespace oofem
#define REGISTER_Element(class)
void giveSwitches(IntArray &answer, int location)
static ParamKey IPK_BondLink3dBoundary_location
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
void computeGeometryProperties() override
int giveLocalCoordinateSystem(FloatMatrix &answer) override
void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
FloatArray directionVector
Definition bondlink3d.h:69
FloatArray rigid
Definition bondlink3d.h:72
double giveLength()
Definition bondlink3d.C:476
double computeVolumeAround(GaussPoint *aGaussPoint) override
Definition bondlink3d.C:76
BondLink3d(int n, Domain *)
Definition bondlink3d.C:69
FloatArray globalCentroid
Definition bondlink3d.h:73
FloatMatrix localCoordinateSystem
Definition bondlink3d.h:67
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
Definition bondlink3d.C:85
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
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 resizeWithData(Index, Index)
Definition floatmatrix.C:91
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.
void assemble(const FloatMatrix &src, const IntArray &loc)
double at(std::size_t i, std::size_t j) const
void beUnitMatrix()
Sets receiver to unity matrix.
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition gausspoint.h:190
int & at(std::size_t i)
Definition intarray.h:104
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
#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_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