OOFEM 3.0
Loading...
Searching...
No Matches
libeam3dboundary.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
37#include "node.h"
38#include "material.h"
39#include "crosssection.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "intarray.h"
44#include "floatarray.h"
45#include "mathfem.h"
46#include "classfactory.h"
47#include "parametermanager.h"
48#include "paramkey.h"
49
50namespace oofem {
54
55LIBeam3dBoundary :: LIBeam3dBoundary(int n, Domain *aDomain) : LIBeam3d(n, aDomain)
56 // Constructor.
57{
59 referenceNode = 0;
60 length = 0.;
61}
62
63
64void
65LIBeam3dBoundary :: initializeFrom(InputRecord &ir, int priority)
66{
67 LIBeam3d :: initializeFrom(ir, priority);
68 ParameterManager &ppm = giveDomain()->elementPPM;
71
72}
73
74void
75LIBeam3dBoundary :: postInitialize()
76{
77 LIBeam3d :: postInitialize();
78 ParameterManager &ppm = giveDomain()->elementPPM;
81 if ( referenceNode == 0 ) {
82 OOFEM_ERROR("wrong reference node specified");
83 }
84}
85
86void
87LIBeam3dBoundary :: giveDofManDofIDMask(int inode, IntArray &answer) const
88{
89 if ( inode == 3 ) {
90 answer = { E_xx, E_xy, E_xz, E_yx, E_yy, E_yz, E_zx, E_zy, E_zz };
91 } else {
92 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
93 }
94}
95
96
97void
98LIBeam3dBoundary :: giveSwitches(IntArray &answer, int location) {
99 answer.resize(3);
100 int counter = 1;
101 for ( int x = -1; x < 2; x++ ) {
102 for ( int y = -1; y < 2; y++ ) {
103 for ( int z = -1; z < 2; z++ ) {
104 if ( !( z == 0 && y == 0 && x == 0 ) ) {
105 if ( counter == location ) {
106 answer(0) = x;
107 answer(1) = y;
108 answer(2) = z;
109 }
110 counter++;
111 }
112 }
113 }
114 }
115 return;
116}
117
118
119void
120LIBeam3dBoundary :: recalculateCoordinates(int nodeNumber, FloatArray &coords)
121{
122 FloatArray unitCellSize;
123 unitCellSize.resize(3);
124 unitCellSize.at(1) = this->giveNode(3)->giveCoordinate(1);
125 unitCellSize.at(2) = this->giveNode(3)->giveCoordinate(2);
126 unitCellSize.at(3) = this->giveNode(3)->giveCoordinate(3);
127
128 IntArray switches;
129 this->giveSwitches(switches, this->location.at(nodeNumber) );
130
131 coords.resize(3);
132 coords.at(1) = this->giveNode(nodeNumber)->giveCoordinate(1) + switches.at(1) * unitCellSize.at(1);
133 coords.at(2) = this->giveNode(nodeNumber)->giveCoordinate(2) + switches.at(2) * unitCellSize.at(2);
134 coords.at(3) = this->giveNode(nodeNumber)->giveCoordinate(3) + switches.at(3) * unitCellSize.at(3);
135
136 return;
137}
138
139
140double
141LIBeam3dBoundary :: computeLength()
142// Returns the length of the receiver.
143{
144 double dx, dy, dz;
145 FloatArray nodeA, nodeB;
146
147 if ( length == 0. ) {
148 recalculateCoordinates(1, nodeA);
149 recalculateCoordinates(2, nodeB);
150 dx = nodeB.at(1) - nodeA.at(1);
151 dy = nodeB.at(2) - nodeA.at(2);
152 dz = nodeB.at(3) - nodeA.at(3);
153 length = sqrt(dx * dx + dy * dy + dz * dz);
154 }
155
156 return length;
157}
158
159
160int
161LIBeam3dBoundary :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
162{
163 double ksi, n1, n2;
164 FloatArray coordsNodeA, coordsNodeB;
165
166 ksi = lcoords.at(1);
167 n1 = ( 1. - ksi ) * 0.5;
168 n2 = ( 1. + ksi ) * 0.5;
169 recalculateCoordinates(1, coordsNodeA);
170 recalculateCoordinates(2, coordsNodeB);
171
172 answer.resize(3);
173 answer.at(1) = n1 * coordsNodeA.at(1) + n2 * coordsNodeB.at(1);
174 answer.at(2) = n1 * coordsNodeA.at(2) + n2 * coordsNodeB.at(2);
175 answer.at(3) = n1 * coordsNodeA.at(3) + n2 * coordsNodeB.at(3);
176
177 return 1;
178}
179
180
181int
182LIBeam3dBoundary :: giveLocalCoordinateSystem(FloatMatrix &answer)
183//
184// returns a unit vectors of local coordinate system at element
185// stored rowwise (mainly used by some materials with ortho and anisotrophy)
186//
187{
188 FloatArray lx(3), ly(3), lz(3), help(3);
189 double length = this->computeLength();
190 Node *refNode;
191 FloatArray coordsNodeA, coordsNodeB;
192
193 answer.resize(3, 3);
194 answer.zero();
195 refNode = this->giveReferenceNode(referenceNode);
196 recalculateCoordinates(1, coordsNodeA);
197 recalculateCoordinates(2, coordsNodeB);
198
199 for ( int i = 1; i <= 3; i++ ) {
200 lx.at(i) = ( coordsNodeB.at(i) - coordsNodeA.at(i) ) / length;
201 help.at(i) = ( refNode->giveCoordinate(i) - coordsNodeA.at(i) );
202 }
203
204 lz.beVectorProductOf(lx, help);
205 lz.normalize();
206 ly.beVectorProductOf(lz, lx);
207 ly.normalize();
208
209 for ( int i = 1; i <= 3; i++ ) {
210 answer.at(1, i) = lx.at(i);
211 answer.at(2, i) = ly.at(i);
212 answer.at(3, i) = lz.at(i);
213 }
214
215 return 1;
216}
217
218
219bool
220LIBeam3dBoundary :: computeGtoLRotationMatrix(FloatMatrix &answer)
221{
222 FloatMatrix lcs;
223
224 answer.resize(this->computeNumberOfDofs(), this->computeNumberOfDofs() );
225 answer.zero();
226
227 this->giveLocalCoordinateSystem(lcs);
228 for ( int i = 1; i <= 3; i++ ) {
229 for ( int j = 1; j <= 3; j++ ) {
230 answer.at(i, j) = lcs.at(i, j);
231 answer.at(i + 3, j + 3) = lcs.at(i, j);
232 answer.at(i + 6, j + 6) = lcs.at(i, j);
233 answer.at(i + 9, j + 9) = lcs.at(i, j);
234 }
235 }
236 //do not rotate extra DOFs
237 for ( int i = 13; i <= this->computeNumberOfDofs(); i++ ) {
238 for ( int j = 13; j <= this->computeNumberOfDofs(); j++ ) {
239 if ( i == j ) {
240 answer.at(i, j) = 1.;
241 }
242 }
243 }
244 return true;
245}
246
247
248void
249LIBeam3dBoundary :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
250// Returns the stiffness matrix of the receiver, expressed in the global
251// axes.
252{
253 //get the raw stiffness matrix in local cs
254 FloatMatrix Korig, GtoL, R, Rt, T, Tt, TtK;
255 StructuralElement :: computeStiffnessMatrix(Korig, rMode, tStep);
256
257 //rotate to global cs
258 this->giveRotationMatrix(GtoL);
259 R.beSubMatrixOf(GtoL, 1, 12, 1, 12);
260 Korig.rotatedWith(R);
261
262 //transform
263 this->computeTransformationMatrix(T, tStep);
264 Tt.beTranspositionOf(T);
265 TtK.beProductOf(Tt, Korig);
266 answer.beProductOf(TtK, T);
267
268 //rotate back to local cs
269 Rt.beTranspositionOf(GtoL);
270 answer.rotatedWith(Rt);
271}
272
273
274void
275LIBeam3dBoundary :: computeTransformationMatrix(FloatMatrix &answer, TimeStep *tStep)
276{
277 FloatArray unitCellSize;
278 unitCellSize.resize(3);
279 unitCellSize.at(1) = this->giveNode(3)->giveCoordinate(1);
280 unitCellSize.at(2) = this->giveNode(3)->giveCoordinate(2);
281 unitCellSize.at(3) = this->giveNode(3)->giveCoordinate(3);
282
283 IntArray switches1, switches2;
284 this->giveSwitches(switches1, this->location.at(1) );
285 this->giveSwitches(switches2, this->location.at(2) );
286
287 FloatMatrix k1, k2;
288 k1.resize(6, 9);
289 k2.resize(6, 9);
290
291 for ( int i = 1; i <= 3; i++ ) {
292 k1.at(i, 3 * i - 2) = unitCellSize.at(1) * switches1.at(1);
293 k1.at(i, 3 * i - 1) = unitCellSize.at(2) * switches1.at(2);
294 k1.at(i, 3 * i) = unitCellSize.at(3) * switches1.at(3);
295 }
296
297 for ( int i = 1; i <= 3; i++ ) {
298 k2.at(i, 3 * i - 2) = unitCellSize.at(1) * switches2.at(1);
299 k2.at(i, 3 * i - 1) = unitCellSize.at(2) * switches2.at(2);
300 k2.at(i, 3 * i) = unitCellSize.at(3) * switches2.at(3);
301 }
302
303 answer.resize(12, 12);
304 answer.beUnitMatrix();
305 answer.resizeWithData(12, 21);
306
307 answer.assemble(k1, { 1, 2, 3, 4, 5, 6 }, { 13, 14, 15, 16, 17, 18, 19, 20, 21 });
308 answer.assemble(k2, { 7, 8, 9, 10, 11, 12 }, { 13, 14, 15, 16, 17, 18, 19, 20, 21 });
309}
310
311
312void
313LIBeam3dBoundary :: computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
314{
315 FloatArray u, uG, dispVecG, dispVecL;
316 FloatMatrix B, T, GtoL, GtoLT, R;
317 this->computeVectorOf(VM_Total, tStep, u);
318 // subtract initial displacements, if defined
319 if ( initialDisplacements ) {
321 }
322
323 //rotate to global
324 this->giveRotationMatrix(GtoL);
325 GtoLT.beTranspositionOf(GtoL);
326 R.beSubMatrixOf(GtoL, 1, 12, 1, 12);
327 uG.beProductOf(GtoLT, u);
328 //transform
329 this->computeTransformationMatrix(T, tStep);
330 dispVecG.beProductOf(T, uG);
331 //rotate to local
332 dispVecL.beProductOf(R, dispVecG);
333
334 this->computeBmatrixAt(gp, B);
335 answer.beProductOf(B, dispVecL);
336}
337
338
339void
340LIBeam3dBoundary :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
341{
342 FloatMatrix B, T, Tt, GtoL, R, Rt;
343 FloatArray vStress, vStrain, fintsub, fintG, fintGT;
344
345 //get the raw internal force vector in local cs
346 answer.clear();
347 fintsub.resize(12);
348
349 for ( auto &gp : * this->giveDefaultIntegrationRulePtr() ) {
350 this->computeBmatrixAt(gp, B);
351
352 if ( useUpdatedGpRecord == 1 ) {
353 StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
354 vStress = matStat->giveStressVector();
355 } else {
356 if ( !this->isActivated(tStep) ) {
357 vStrain.resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
358 vStrain.zero();
359 }
360 this->computeStrainVector(vStrain, gp, tStep);
361 this->computeStressVector(vStress, vStrain, gp, tStep);
362 }
363
364 // Compute nodal internal forces at nodes as f = B^T*Stress dV
365 double dV = this->computeVolumeAround(gp);
366
367 if ( vStress.giveSize() == 6 ) {
368 FloatArray stressTemp;
369 StructuralMaterial :: giveReducedSymVectorForm(stressTemp, vStress, gp->giveMaterialMode() );
370 fintsub.plusProduct(B, stressTemp, dV);
371 } else {
372 fintsub.plusProduct(B, vStress, dV);
373 }
374 }
375
376 //rotate to global cs
377 this->giveRotationMatrix(GtoL);
378 R.beSubMatrixOf(GtoL, 1, 12, 1, 12);
379 Rt.beTranspositionOf(R);
380 fintG.beProductOf(Rt, fintsub);
381 //transform
382 this->computeTransformationMatrix(T, tStep);
383 Tt.beTranspositionOf(T);
384 fintGT.beProductOf(Tt, fintG);
385 //rotate back to local cs
386 answer.beProductOf(GtoL, fintGT);
387}
388
389
390int
391LIBeam3dBoundary :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
392{
393 if ( type == IST_DisplacementVector ) {
394 FloatArray u;
396 this->computeVectorOf(VM_Total, tStep, u);
397 u.resizeWithValues(12);
398 this->computeNmatrixAt(gp->giveSubPatchCoordinates(), N); //no special treatment needed for this interpolation
399 answer.beProductOf(N, u);
400 return 1;
401 }
402 return LIBeam3d :: giveIPValue(answer, gp, type, tStep);
403}
404} // end namespace oofem
#define N(a, b)
#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
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Definition element.h:886
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 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 resizeWithValues(Index s, std::size_t allocChunk=0)
Definition floatarray.C:103
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void subtract(const FloatArray &src)
Definition floatarray.C:320
void rotatedWith(const FloatMatrix &r, char mode='n')
void resizeWithData(Index, Index)
Definition floatmatrix.C:91
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beSubMatrixOf(const FloatMatrix &src, Index topRow, Index bottomRow, Index topCol, Index bottomCol)
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.
const FloatArray & giveSubPatchCoordinates() const
Returns local sub-patch coordinates of the receiver.
Definition gausspoint.h:142
void resize(int n)
Definition intarray.C:73
int & at(std::size_t i)
Definition intarray.h:104
void recalculateCoordinates(int nodeNumber, FloatArray &coords) override
static ParamKey IPK_LIBeam3dBoundary_refnode
double computeLength() override
virtual void computeTransformationMatrix(FloatMatrix &answer, TimeStep *tStep)
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override
void giveSwitches(IntArray &answer, int location)
static ParamKey IPK_LIBeam3dBoundary_location
int giveLocalCoordinateSystem(FloatMatrix &answer) override
int computeNumberOfDofs() override
double computeVolumeAround(GaussPoint *gp) override
Definition libeam3d.C:210
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition libeam3d.C:129
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition libeam3d.C:68
Node * giveReferenceNode(int refNode)
Definition libeam3d.C:386
LIBeam3d(int n, Domain *d)
Definition libeam3d.C:58
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition libeam3d.C:151
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
#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