OOFEM 3.0
Loading...
Searching...
No Matches
trplanestressrotallman3d.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
50TrPlanestressRotAllman3d :: TrPlanestressRotAllman3d(int n, Domain *aDomain) : TrPlanestressRotAllman(n, aDomain)
51{
52 GtoLRotationMatrix = NULL;
53}
54
55
56void
57TrPlanestressRotAllman3d :: computeLocalNodalCoordinates(std :: vector< FloatArray > &lxy)
58// Returns global coordinates given in global vector
59// transformed into local coordinate system of the
60// receiver
61{
62 // test if pereviously computed
63 if ( GtoLRotationMatrix == NULL ) {
65 }
66
67
68 lxy.resize(6);
69 for ( int i = 0; i < 3; i++ ) {
70 const auto &nc = this->giveNode(i + 1)->giveCoordinates();
71 lxy [ i ].beProductOf(* GtoLRotationMatrix, nc);
72 }
73 lxy [ 3 ].resize(3);
74 lxy [ 4 ].resize(3);
75 lxy [ 5 ].resize(3);
76 for ( int i = 1; i <= 3; i++ ) {
77 lxy [ 3 ].at(i) = 0.5 * ( lxy [ 0 ].at(i) + lxy [ 1 ].at(i) );
78 lxy [ 4 ].at(i) = 0.5 * ( lxy [ 1 ].at(i) + lxy [ 2 ].at(i) );
79 lxy [ 5 ].at(i) = 0.5 * ( lxy [ 2 ].at(i) + lxy [ 0 ].at(i) );
80 }
81}
82
83
84double
85TrPlanestressRotAllman3d :: computeVolumeAround(GaussPoint *gp)
86{
87 double detJ, weight;
88
89 std :: vector< FloatArray >lc;
91
92 weight = gp->giveWeight();
93 detJ = fabs( this->interp.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIVertexListGeometryWrapper(lc, this->interp.giveGeometryType()) ) );
94 return detJ * weight * this->giveStructuralCrossSection()->give(CS_Thickness, gp);
95}
96
97
98void
99TrPlanestressRotAllman3d :: giveDofManDofIDMask(int inode, IntArray &answer) const
100{
101 answer = {
102 D_u, D_v, D_w, R_u, R_v, R_w
103 };
104}
105
106
107const FloatMatrix *
108TrPlanestressRotAllman3d :: computeGtoLRotationMatrix()
109// Returns the rotation matrix of the receiver of the size [3,3]
110// coords(local) = T * coords(global)
111//
112// local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
113//
114// e1' : [N2-N1] Ni - means i - th node
115// help : [N3-N1]
116// e3' : e1' x help
117// e2' : e3' x e1'
118{
119 if ( GtoLRotationMatrix == NULL ) {
120 FloatArray e1(3), e2(3), e3(3), help(3);
121
122 // compute e1' = [N2-N1] and help = [N3-N1]
123 e1.beDifferenceOf(this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates());
124 help.beDifferenceOf(this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates());
125
126 // let us normalize e1'
127 e1.normalize();
128
129 // compute e3' : vector product of e1' x help
130 e3.beVectorProductOf(e1, help);
131 // let us normalize
132 e3.normalize();
133
134 // now from e3' x e1' compute e2'
135 e2.beVectorProductOf(e3, e1);
136
137 //
139
140 for ( int i = 1; i <= 3; i++ ) {
141 GtoLRotationMatrix->at(1, i) = e1.at(i);
142 GtoLRotationMatrix->at(2, i) = e2.at(i);
143 GtoLRotationMatrix->at(3, i) = e3.at(i);
144 }
145 }
146
147 return GtoLRotationMatrix;
148}
149
150
151bool
152TrPlanestressRotAllman3d :: computeGtoLRotationMatrix(FloatMatrix &answer)
153// Returns the rotation matrix of the receiver of the size [9,18]
154// r(local) = T * r(global)
155// for one node (r written transposed): {u,v,r3} = T * {u,v,w,r1,r2,r3}
156{
157 // test if pereviously computed
158 if ( GtoLRotationMatrix == NULL ) {
160 }
161
162 answer.resize(9, 18);
163 answer.zero();
164
165 for ( int i = 1; i <= 3; i++ ) {
166 answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix->at(1, i);
167 answer.at(2, i) = answer.at(2 + 3, i + 6) = answer.at(2 + 6, i + 12) = GtoLRotationMatrix->at(2, i);
168 answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix->at(3, i);
169 }
170
171 return 1;
172}
173
174
175void
176TrPlanestressRotAllman3d :: giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
177// returns characteristic tensor of the receiver at given gp and tStep
178// strain vector = (Eps_X, Eps_y, Gamma_xy, Kappa_z)
179{
180 FloatArray charVect;
182
183 answer.resize(3, 3);
184 answer.zero();
185
186 if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
187 //this->computeStressVector(charVect, gp, tStep);
188 charVect = ms->giveStressVector();
189
190 answer.at(1, 1) = charVect.at(1);
191 answer.at(2, 2) = charVect.at(2);
192 answer.at(1, 2) = charVect.at(3);
193 answer.at(2, 1) = charVect.at(3);
194 } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {} else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
195 //this->computeStrainVector(charVect, gp, tStep);
196 charVect = ms->giveStrainVector();
197
198 answer.at(1, 1) = charVect.at(1);
199 answer.at(2, 2) = charVect.at(2);
200 answer.at(1, 2) = charVect.at(3) / 2.;
201 answer.at(2, 1) = charVect.at(3) / 2.;
202 } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {} else {
203 OOFEM_ERROR("unsupported tensor mode");
204 }
205
206 if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
207 ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
210 }
211}
212
213
214int
215TrPlanestressRotAllman3d :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
216{
217 FloatMatrix globTensor;
218 CharTensor cht;
219
220 if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
221 double c = 1.0;
222 if ( type == IST_ShellForceTensor ) {
223 cht = GlobalForceTensor;
224 } else {
225 cht = GlobalStrainTensor;
226 c = 2.0; // tensor components reported
227 }
228
229 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
230
231 answer.resize(6);
232 answer.at(1) = globTensor.at(1, 1); //xx
233 answer.at(2) = globTensor.at(2, 2); //yy
234 answer.at(3) = globTensor.at(3, 3); //zz
235 answer.at(4) = c * globTensor.at(2, 3); //yz
236 answer.at(5) = c * globTensor.at(1, 3); //xz
237 answer.at(6) = c * globTensor.at(1, 2); //xy
238 // mutiply stresses by thickness to get forces
239
240 return 1;
241 } else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
242 answer.clear();
243 return 1;
244 } else {
245 return TrPlanestressRotAllman :: giveIPValue(answer, gp, type, tStep);
246 }
247}
248
249int
250TrPlanestressRotAllman3d :: computeLoadGToLRotationMtrx(FloatMatrix &answer)
251// Returns the rotation matrix of the receiver of the size [6,6]
252// f(local) = T * f(global)
253{
254 // test if previously computed
255 if ( GtoLRotationMatrix == NULL ) {
257 }
258
259 answer.resize(6, 6);
260 answer.zero();
261
262 for ( int i = 1; i <= 3; i++ ) {
263 answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix->at(1, i);
264 answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix->at(2, i);
265 answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix->at(3, i);
266 }
267
268 return 1;
269}
270
271
272void
273TrPlanestressRotAllman3d :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
274// Computes numerically the load vector of the receiver due to the body loads, at tStep.
275// load is assumed to be in global cs.
276// load vector is then transformed to coordinate system in each node.
277// (should be global coordinate system, but there may be defined
278// different coordinate system in each node)
279{
280 double dens, dV, load;
281 GaussPoint *gp = NULL;
282 FloatArray force;
283 FloatMatrix T;
284
285 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
286 OOFEM_ERROR("unknown load type");
287 }
288
289 GaussIntegrationRule irule(0, this);
290 irule.SetUpPointsOnTriangle( 1, this->giveMaterialMode() );
291
292 // note: force is assumed to be in global coordinate system.
293 forLoad->computeComponentArrayAt(force, tStep, mode);
294
295 if ( force.giveSize() ) {
296 gp = irule.getIntegrationPoint(0);
297
298 dens = this->giveStructuralCrossSection()->give('d', gp);
299 dV = this->computeVolumeAround(gp);
300
301 answer.resize(18);
302 answer.zero();
303
304 load = force.at(1) * dens * dV / 3.0;
305 answer.at(1) = load;
306 answer.at(7) = load;
307 answer.at(13) = load;
308
309 load = force.at(2) * dens * dV / 3.0;
310 answer.at(2) = load;
311 answer.at(8) = load;
312 answer.at(14) = load;
313
314 load = force.at(3) * dens * dV / 3.0;
315 answer.at(3) = load;
316 answer.at(9) = load;
317 answer.at(15) = load;
318
319 // transform result from global cs to local element cs.
320 if ( this->computeGtoLRotationMatrix(T) ) {
321 answer.rotatedWith(T, 'n');
322 }
323 } else {
324 answer.clear(); // nil resultant
325 }
326}
327
328void
329TrPlanestressRotAllman3d :: printOutputAt(FILE *file, TimeStep *tStep)
330// Performs end-of-step operations.
331{
332 FloatArray v;
333
334 fprintf( file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
335
336 for ( int i = 0; i < ( int ) integrationRulesArray.size(); i++ ) {
337 for ( GaussPoint *gp : *integrationRulesArray [ i ] ) {
338 fprintf( file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
339
340 this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
341 fprintf(file, " strains ");
342 for ( auto &val : v ) {
343 fprintf(file, " %.4e", val);
344 }
345
346 this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
347 fprintf(file, "\n curvatures ");
348 for ( auto &val : v ) {
349 fprintf(file, " %.4e", val);
350 }
351
352 // Forces - Moments
353 this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
354 fprintf(file, "\n stresses ");
355 for ( auto &val : v ) {
356 fprintf(file, " %.4e", val);
357 }
358
359 this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
360 fprintf(file, "\n moments ");
361 for ( auto &val : v ) {
362 fprintf(file, " %.4e", val);
363 }
364
365 fprintf(file, "\n");
366 }
367 }
368}
369} // 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 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
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.
static FEI2dTrLin interp
Definition trplanstrss.h:70
void computeLocalNodalCoordinates(std::vector< FloatArray > &lxy) override
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
double computeVolumeAround(GaussPoint *gp) override
MaterialMode giveMaterialMode() override
#define OOFEM_ERROR(...)
Definition error.h:79
@ 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