OOFEM 3.0
Loading...
Searching...
No Matches
trplanrot3d.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 "material.h"
39#include "node.h"
40#include "load.h"
41#include "mathfem.h"
42#include "classfactory.h"
44#include "gausspoint.h"
45#include "fei2dtrlin.h"
46
47namespace oofem {
49
50TrPlaneStrRot3d :: TrPlaneStrRot3d(int n, Domain *aDomain) : TrPlaneStrRot(n, aDomain)
51{
52}
53
54
55void
56TrPlaneStrRot3d :: giveLocalCoordinates(FloatArray &answer, const FloatArray &global)
57// Returns global coordinates given in global vector
58// transformed into local coordinate system of the
59// receiver
60{
61 // test the parameter
62 if ( global.giveSize() != 3 ) {
63 OOFEM_ERROR("cannot transform coordinates - size mismatch");
64 }
65
66 // first ensure that receiver's GtoLRotationMatrix[3,3] is defined
67 if ( !GtoLRotationMatrix.isNotEmpty() ) {
69 }
70
71 FloatArray offset = global;
72 offset.subtract( this->giveNode(1)->giveCoordinates() );
73 answer.beProductOf(GtoLRotationMatrix, offset);
74}
75
76
77double
78TrPlaneStrRot3d :: computeVolumeAround(GaussPoint *gp)
79{
80 double detJ, weight;
81
82 FloatArray x, y;
83 this->giveNodeCoordinates(x, y);
84 std :: vector< FloatArray > lc = {Vec2(x[0], y[0]), Vec2(x[1], y[1]), Vec2(x[2], y[2])};
85
86 weight = gp->giveWeight();
87 detJ = fabs( this->interp.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIVertexListGeometryWrapper(lc, this->interp.giveGeometryType()) ) );
88 return detJ * weight * this->giveStructuralCrossSection()->give(CS_Thickness, gp);
89}
90
91
92void
93TrPlaneStrRot3d :: giveNodeCoordinates(FloatArray &x, FloatArray &y)
94{
95 FloatArray nc1(3), nc2(3), nc3(3);
96
97 this->giveLocalCoordinates( nc1, this->giveNode(1)->giveCoordinates() );
98 this->giveLocalCoordinates( nc2, this->giveNode(2)->giveCoordinates() );
99 this->giveLocalCoordinates( nc3, this->giveNode(3)->giveCoordinates() );
100
101 x.resize(3);
102 x.at(1) = nc1.at(1);
103 x.at(2) = nc2.at(1);
104 x.at(3) = nc3.at(1);
105
106 y.resize(3);
107 y.at(1) = nc1.at(2);
108 y.at(2) = nc2.at(2);
109 y.at(3) = nc3.at(2);
110
111 //if (z) {
112 // z[0] = nc1->at(3);
113 // z[1] = nc2->at(3);
114 // z[2] = nc3->at(3);
115 //}
116}
117
118
119void
120TrPlaneStrRot3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
121{
122 answer = {D_u, D_v, D_w, R_u, R_v, R_w};
123}
124
125
126const FloatMatrix *
127TrPlaneStrRot3d :: computeGtoLRotationMatrix()
128// Returns the rotation matrix of the receiver of the size [3,3]
129// coords(local) = T * coords(global)
130//
131// local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
132//
133// e1' : [N2-N1] Ni - means i - th node
134// help : [N3-N1]
135// e3' : e1' x help
136// e2' : e3' x e1'
137{
138 if ( !GtoLRotationMatrix.isNotEmpty() ) {
139 FloatArray e1, e2, e3, help;
140
141 // compute e1' = [N2-N1] and help = [N3-N1]
142 e1.beDifferenceOf( this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
143 help.beDifferenceOf( this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
144
145 // let us normalize e1'
146 e1.normalize();
147
148 // compute e3' : vector product of e1' x help
149 e3.beVectorProductOf(e1, help);
150 // let us normalize
151 e3.normalize();
152
153 // now from e3' x e1' compute e2'
154 e2.beVectorProductOf(e3, e1);
155
156 //
157 GtoLRotationMatrix.resize(3, 3);
158
159 for ( int i = 1; i <= 3; i++ ) {
160 GtoLRotationMatrix.at(1, i) = e1.at(i);
161 GtoLRotationMatrix.at(2, i) = e2.at(i);
162 GtoLRotationMatrix.at(3, i) = e3.at(i);
163 }
164 }
165
166 return &GtoLRotationMatrix;
167}
168
169
170bool
171TrPlaneStrRot3d :: computeGtoLRotationMatrix(FloatMatrix &answer)
172// Returns the rotation matrix of the receiver of the size [9,18]
173// r(local) = T * r(global)
174// for one node (r written transposed): {u,v,r3} = T * {u,v,w,r1,r2,r3}
175{
176 // test if pereviously computed
177 if ( !GtoLRotationMatrix.isNotEmpty() ) {
179 }
180
181 answer.resize(9, 18);
182 answer.zero();
183
184 for ( int i = 1; i <= 3; i++ ) {
185 answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix.at(1, i);
186 answer.at(2, i) = answer.at(2 + 3, i + 6) = answer.at(2 + 6, i + 12) = GtoLRotationMatrix.at(2, i);
187 answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix.at(3, i);
188 }
189
190 return 1;
191}
192
193
194void
195TrPlaneStrRot3d :: giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
196// returns characteristic tensor of the receiver at given gp and tStep
197// strain vector = (Eps_X, Eps_y, Gamma_xy, Kappa_z)
198{
199 FloatArray charVect;
201
202 answer.resize(3, 3);
203 answer.zero();
204
205 if ( type == LocalForceTensor || type == GlobalForceTensor ) {
206 //this->computeStressVector(charVect, gp, tStep);
207 charVect = ms->giveStressVector();
208
209 double h = this->giveStructuralCrossSection()->give(CS_Thickness, gp);
210 answer.at(1, 1) = charVect.at(1) * h;
211 answer.at(2, 2) = charVect.at(2) * h;
212 answer.at(1, 2) = charVect.at(3) * h;
213 answer.at(2, 1) = charVect.at(3) * h;
214 } else if ( type == LocalMomentTensor || type == GlobalMomentTensor ) {
215 //this->computeStressVector(charVect, gp, tStep);
216 charVect = ms->giveStressVector();
217
218 answer.at(3, 3) = charVect.at(4);
219 } else if ( type == LocalStrainTensor || type == GlobalStrainTensor ) {
220 //this->computeStrainVector(charVect, gp, tStep);
221 charVect = ms->giveStrainVector();
222
223 answer.at(1, 1) = charVect.at(1);
224 answer.at(2, 2) = charVect.at(2);
225 answer.at(1, 2) = charVect.at(3) / 2.;
226 answer.at(2, 1) = charVect.at(3) / 2.;
227 } else if ( type == LocalCurvatureTensor || type == GlobalCurvatureTensor ) {
228 //this->computeStrainVector(charVect, gp, tStep);
229 charVect = ms->giveStrainVector();
230
231 answer.at(3, 3) = charVect.at(4);
232 } else {
233 OOFEM_ERROR("unsupported tensor mode");
234 }
235
236 if ( type == GlobalForceTensor || type == GlobalMomentTensor ||
237 type == GlobalStrainTensor || type == GlobalCurvatureTensor ) {
240 }
241}
242
243
244int
245TrPlaneStrRot3d :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
246{
247 FloatMatrix globTensor;
248 CharTensor cht;
249
250 answer.resize(6);
251
252 if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
253 if ( type == IST_CurvatureTensor ) {
255 } else {
256 cht = GlobalStrainTensor;
257 }
258
259 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
260
261 answer.at(1) = globTensor.at(1, 1); //xx
262 answer.at(2) = globTensor.at(2, 2); //yy
263 answer.at(3) = globTensor.at(3, 3); //zz
264 answer.at(4) = 2 * globTensor.at(2, 3); //yz
265 answer.at(5) = 2 * globTensor.at(1, 3); //xz
266 answer.at(6) = 2 * globTensor.at(1, 2); //xy
267
268 return 1;
269 } else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
270 if ( type == IST_ShellMomentTensor ) {
271 cht = GlobalMomentTensor;
272 } else {
273 cht = GlobalForceTensor;
274 }
275
276 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
277
278 answer.at(1) = globTensor.at(1, 1); //xx
279 answer.at(2) = globTensor.at(2, 2); //yy
280 answer.at(3) = globTensor.at(3, 3); //zz
281 answer.at(4) = globTensor.at(2, 3); //yz
282 answer.at(5) = globTensor.at(1, 3); //xz
283 answer.at(6) = globTensor.at(1, 2); //xy
284
285 return 1;
286 } else {
287 return TrPlaneStrRot :: giveIPValue(answer, gp, type, tStep);
288 }
289}
290
291int
292TrPlaneStrRot3d :: computeLoadGToLRotationMtrx(FloatMatrix &answer)
293// Returns the rotation matrix of the receiver of the size [6,6]
294// f(local) = T * f(global)
295{
296 // test if previously computed
297 if ( !GtoLRotationMatrix.isNotEmpty() ) {
299 }
300
301 answer.resize(6, 6);
302 answer.zero();
303
304 for ( int i = 1; i <= 3; i++ ) {
305 answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix.at(1, i);
306 answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix.at(2, i);
307 answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix.at(3, i);
308 }
309
310 return 1;
311}
312
313
314void
315TrPlaneStrRot3d :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
316// Computes numerically the load vector of the receiver due to the body loads, at tStep.
317// load is assumed to be in global cs.
318// load vector is then transformed to coordinate system in each node.
319// (should be global coordinate system, but there may be defined
320// different coordinate system in each node)
321{
322 double dens, dV, load;
323 FloatArray force;
324 FloatMatrix T;
325
326 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
327 OOFEM_ERROR("unknown load type");
328 }
329
330 // note: force is assumed to be in global coordinate system; 6 components
331 forLoad->computeComponentArrayAt(force, tStep, mode);
332 // get it in local c.s.
334 force.rotatedWith(T, 'n');
335
336 if ( force.giveSize() ) {
337 GaussIntegrationRule iRule (1, this, 1, 1);
338 iRule.SetUpPointsOnTriangle(1, _Unknown);
339 GaussPoint *gp = iRule.getIntegrationPoint(0);
340
341 dens = this->giveStructuralCrossSection()->give('d', gp);
342 dV = this->computeVolumeAround(gp);// * this->giveCrossSection()->give(CS_Thickness, gp);
343
344 answer.resize(9);
345 answer.zero();
346
347 load = force.at(1) * dens * dV / 3.0;
348 answer.at(1) = load;
349 answer.at(4) = load;
350 answer.at(7) = load;
351
352 load = force.at(2) * dens * dV / 3.0;
353 answer.at(2) = load;
354 answer.at(5) = load;
355 answer.at(8) = load;
356
357 load = force.at(6) * dens * dV / 3.0;
358 answer.at(3) = load;
359 answer.at(6) = load;
360 answer.at(9) = load;
361
362 // transform result from global cs to local element cs.
363 //if ( this->computeGtoLRotationMatrix(T) ) {
364 // answer.rotatedWith(T, 'n');
365 //}
366 } else {
367 answer.clear(); // nil resultant
368 }
369}
370
371void
372TrPlaneStrRot3d :: computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *sgp)
373{
374 FloatMatrix ne;
376
377 answer.resize(6, 18);
378 answer.zero();
379 int ri[] = {
380 0, 1, 5
381 };
382 int ci[] = {
383 0, 1, 5, 6, 7, 11, 12, 13, 17
384 };
385
386 for ( int i = 0; i < 3; i++ ) {
387 for ( int j = 0; j < 9; j++ ) {
388 answer(ri [ i ], ci [ j ]) = ne(i, j);
389 }
390 }
391}
392
393void
394TrPlaneStrRot3d :: giveSurfaceDofMapping(IntArray &answer, int iSurf) const
395{
396 answer.resize(18);
397 answer.zero();
398 if ( iSurf == 1 ) {
399 answer.at(1) = 1; // node 1
400 answer.at(2) = 2;
401 answer.at(6) = 3;
402
403 answer.at(7) = 4; // node 2
404 answer.at(8) = 5;
405 answer.at(12) = 6;
406
407 answer.at(13) = 7; // node 3
408 answer.at(14) = 8;
409 answer.at(18) = 9;
410 } else {
411 OOFEM_ERROR("wrong surface number");
412 }
413}
414
415
416double
417TrPlaneStrRot3d :: computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
418{
419 return this->computeVolumeAround(gp);
420}
421
422
423int
424TrPlaneStrRot3d :: computeLoadLSToLRotationMatrix(FloatMatrix &answer, int isurf, GaussPoint *gp)
425{
426 return 0;
427}
428
429
430void
431TrPlaneStrRot3d :: printOutputAt(FILE *file, TimeStep *tStep)
432// Performs end-of-step operations.
433{
434 FloatArray v;
435
436 fprintf( file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
437
438 for ( int i = 0; i < (int)integrationRulesArray.size(); i++ ) {
439 for ( GaussPoint *gp: *integrationRulesArray [ i ] ) {
440
441 fprintf( file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
442
443 this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
444 fprintf(file, " strains ");
445 for ( auto &val : v ) fprintf(file, " %.4e", val);
446
447 this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
448 fprintf(file, "\n curvatures ");
449 for ( auto &val : v ) fprintf(file, " %.4e", val);
450
451 // Forces - Moments
452 this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
453 fprintf(file, "\n stresses ");
454 for ( auto &val : v ) fprintf(file, " %.4e", val);
455
456 this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
457 fprintf(file, "\n moments ");
458 for ( auto &val : v ) fprintf(file, " %.4e", val);
459
460 fprintf(file, "\n");
461 }
462 }
463}
464} // end namespace oofem
#define REGISTER_Element(class)
Node * giveNode(int i) const
Definition element.h:629
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int giveLabel() const
Definition element.h:1125
int giveNumber() const
Definition femcmpnn.h:104
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Definition floatarray.C:403
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 beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void subtract(const FloatArray &src)
Definition floatarray.C:320
void rotatedWith(const FloatMatrix &r, char mode='n')
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
int SetUpPointsOnTriangle(int nPoints, MaterialMode mode) override
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
Definition gausspoint.h:138
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
Definition gausspoint.h:204
double giveWeight()
Returns integration weight of receiver.
Definition gausspoint.h:180
virtual bcGeomType giveBCGeoType() const
void resize(int n)
Definition intarray.C:73
void zero()
Sets all component to zero.
Definition intarray.C:52
int & at(std::size_t i)
Definition intarray.h:104
GaussPoint * getIntegrationPoint(int n)
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition load.C:84
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver's strain vector.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
void giveNodeCoordinates(FloatArray &x, FloatArray &y) override
Definition trplanrot3d.C:93
const FloatMatrix * computeGtoLRotationMatrix()
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
void giveLocalCoordinates(FloatArray &answer, const FloatArray &global)
Definition trplanrot3d.C:56
double computeVolumeAround(GaussPoint *gp) override
Definition trplanrot3d.C:78
int computeLoadGToLRotationMtrx(FloatMatrix &answer) override
FloatMatrix GtoLRotationMatrix
Definition trplanrot3d.h:75
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition trplanrot.C:244
TrPlaneStrRot(int, Domain *)
Definition trplanrot.C:57
static FEI2dTrLin interp
Definition trplanstrss.h:70
#define OOFEM_ERROR(...)
Definition error.h:79
static FloatArray Vec2(const double &a, const double &b)
Definition floatarray.h:606
@ CS_Thickness
Thickness.
@ BodyLoadBGT
Distributed body load.
Definition bcgeomtype.h:43
@ ForceLoadBVT
Definition bcvaltype.h:43

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