OOFEM 3.0
Loading...
Searching...
No Matches
cct3d.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 "fei2dtrlin.h"
38#include "node.h"
39#include "load.h"
40#include "mathfem.h"
42#include "gausspoint.h"
43#include "classfactory.h"
44
45#include <cstdlib>
46
47namespace oofem {
49
50CCTPlate3d::CCTPlate3d(int n, Domain *aDomain) : CCTPlate(n, aDomain)
51{}
52
53
54void
56{
57 if ( global.giveSize() != 3 ) {
58 OOFEM_ERROR("cannot transform coordinates - size mismatch");
59 }
60
61 // first ensure that receiver's GtoLRotationMatrix[3,3] is defined
63
64 FloatArray offset;
65 offset.beDifferenceOf( global, this->giveNode(1)->giveCoordinates() );
66 answer.beProductOf(GtoLRotationMatrix, offset);
67}
68
69
70void
71CCTPlate3d::giveNodeCoordinates(double &x1, double &x2, double &x3,
72 double &y1, double &y2, double &y3,
73 double &z1, double &z2, double &z3)
74{
75 FloatArray nc1(3), nc2(3), nc3(3);
76
77 this->giveLocalCoordinates(nc1, this->giveNode(1)->giveCoordinates() );
78 this->giveLocalCoordinates(nc2, this->giveNode(2)->giveCoordinates() );
79 this->giveLocalCoordinates(nc3, this->giveNode(3)->giveCoordinates() );
80
81 x1 = nc1.at(1);
82 x2 = nc2.at(1);
83 x3 = nc3.at(1);
84
85 y1 = nc1.at(2);
86 y2 = nc2.at(2);
87 y3 = nc3.at(2);
88
89 z1 = nc1.at(3);
90 z2 = nc2.at(3);
91 z3 = nc3.at(3);
92}
93
94
95void
97{
98 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
99}
100
101
102bool
104//converts global coordinates to local planar area coordinates,
105//does not return a coordinate in the thickness direction, but
106//does check that the point is in the element thickness
107{
108 // rotate the input point Coordinate System into the element CS
109 FloatArray inputCoords_ElCS;
110 std::vector< FloatArray >lc(3);
111 FloatArray llc;
112 this->giveLocalCoordinates(inputCoords_ElCS, coords);
113 for ( int _i = 0; _i < 3; _i++ ) {
114 this->giveLocalCoordinates(lc [ _i ], this->giveNode(_i + 1)->giveCoordinates() );
115 }
116 FEI2dTrLin _interp(1, 2);
117 bool inplane = _interp.global2local( llc, inputCoords_ElCS, FEIVertexListGeometryWrapper(lc, this->giveGeometryType()) ) > 0;
118 answer.resize(2);
119 answer.at(1) = inputCoords_ElCS.at(1);
120 answer.at(2) = inputCoords_ElCS.at(2);
121 GaussPoint _gp(NULL, 1, answer, 2.0, _2dPlate);
122 // now check if the third local coordinate is within the thickness of element
123 bool outofplane = ( fabs(inputCoords_ElCS.at(3) ) <= this->giveCrossSection()->give(CS_Thickness, & _gp) / 2. );
124
125 return inplane && outofplane;
126}
127
128
129int
131{
132 double l1 = lcoords.at(1);
133 double l2 = lcoords.at(2);
134 double l3 = 1. - l2 - l1;
135
136 answer.resize(3);
137 for ( int _i = 1; _i <= 3; _i++ ) {
138 answer.at(_i) = l1 * this->giveNode(1)->giveCoordinate(_i) + l2 * this->giveNode(2)->giveCoordinate(_i) + l3 * this->giveNode(3)->giveCoordinate(_i);
139 }
140 return true;
141}
142
143
144const FloatMatrix *
146// Returns the rotation matrix of the receiver of the size [3,3]
147// coords(local) = T * coords(global)
148//
149// local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
150//
151// e1' : [N2-N1] Ni - means i - th node
152// help : [N3-N1]
153// e3' : e1' x help
154// e2' : e3' x e1'
155{
156 if ( !GtoLRotationMatrix.isNotEmpty() ) {
157 FloatArray e1, e2, e3, help;
158
159 // compute e1' = [N2-N1] and help = [N3-N1]
160 e1.beDifferenceOf( this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
161 help.beDifferenceOf( this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
162
163 // let us normalize e1'
164 e1.normalize();
165
166 // compute e3' : vector product of e1' x help
167 e3.beVectorProductOf(e1, help);
168 // let us normalize
169 e3.normalize();
170
171 // now from e3' x e1' compute e2'
172 e2.beVectorProductOf(e3, e1);
173
174 //
175 GtoLRotationMatrix.resize(3, 3);
176
177 for ( int i = 1; i <= 3; i++ ) {
178 GtoLRotationMatrix.at(1, i) = e1.at(i);
179 GtoLRotationMatrix.at(2, i) = e2.at(i);
180 GtoLRotationMatrix.at(3, i) = e3.at(i);
181 }
182 }
183
184 return & GtoLRotationMatrix;
185}
186
187
188bool
190// Returns the rotation matrix of the receiver of the size [9,18]
191// r(local) = T * r(global)
192// for one node (r written transposed): {w,r1,r2} = T * {u,v,w,r1,r2,r3}
193{
195
196 answer.resize(9, 18);
197 answer.zero();
198
199 for ( int i = 1; i <= 3; i++ ) {
200 answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix.at(3, i);
201 answer.at(2, i + 3) = answer.at(2 + 3, i + 3 + 6) = answer.at(2 + 6, i + 3 + 12) = GtoLRotationMatrix.at(1, i);
202 answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix.at(2, i);
203 }
204
205 return 1;
206}
207
208void
210// returns characteristic tensor of the receiver at given gp and tStep
211// strain vector = (Kappa_x, Kappa_y, Kappa_xy, Gamma_zx, Gamma_zy)
212{
213 FloatArray charVect;
215
216 answer.resize(3, 3);
217 answer.zero();
218
219 if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
220 //this->computeStressVector(charVect, gp, tStep);
221 charVect = ms->giveStressVector();
222
223 answer.at(1, 3) = charVect.at(4);
224 answer.at(3, 1) = charVect.at(4);
225 answer.at(2, 3) = charVect.at(5);
226 answer.at(3, 2) = charVect.at(5);
227 } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {
228 //this->computeStressVector(charVect, gp, tStep);
229 charVect = ms->giveStressVector();
230
231 answer.at(1, 1) = charVect.at(1);
232 answer.at(2, 2) = charVect.at(2);
233 answer.at(1, 2) = charVect.at(3);
234 answer.at(2, 1) = charVect.at(3);
235 } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
236 //this->computeStrainVector(charVect, gp, tStep);
237 charVect = ms->giveStrainVector();
238
239 answer.at(1, 3) = charVect.at(4) / 2.;
240 answer.at(3, 1) = charVect.at(4) / 2.;
241 answer.at(2, 3) = charVect.at(5) / 2.;
242 answer.at(3, 2) = charVect.at(5) / 2.;
243 } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {
244 //this->computeStrainVector(charVect, gp, tStep);
245 charVect = ms->giveStrainVector();
246
247 answer.at(1, 1) = charVect.at(1);
248 answer.at(2, 2) = charVect.at(2);
249 answer.at(1, 2) = charVect.at(3) / 2.;
250 answer.at(2, 1) = charVect.at(3) / 2.;
251 } else {
252 OOFEM_ERROR("unsupported tensor mode");
253 }
254
255 if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
256 ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
259 }
260}
261
262
263int
265{
266 FloatMatrix globTensor;
267 CharTensor cht;
268
269 answer.resize(6);
270
271 if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
272 if ( type == IST_CurvatureTensor ) {
274 } else {
275 cht = GlobalStrainTensor;
276 }
277
278 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
279
280 answer.at(1) = globTensor.at(1, 1); //xx
281 answer.at(2) = globTensor.at(2, 2); //yy
282 answer.at(3) = globTensor.at(3, 3); //zz
283 answer.at(4) = 2 * globTensor.at(2, 3); //yz
284 answer.at(5) = 2 * globTensor.at(1, 3); //xz
285 answer.at(6) = 2 * globTensor.at(1, 2); //xy
286
287 return 1;
288 } else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
289 if ( type == IST_ShellMomentTensor ) {
290 cht = GlobalMomentTensor;
291 } else {
292 cht = GlobalForceTensor;
293 }
294
295 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
296
297 answer.at(1) = globTensor.at(1, 1); //xx
298 answer.at(2) = globTensor.at(2, 2); //yy
299 answer.at(3) = globTensor.at(3, 3); //zz
300 answer.at(4) = globTensor.at(2, 3); //yz
301 answer.at(5) = globTensor.at(1, 3); //xz
302 answer.at(6) = globTensor.at(1, 2); //xy
303
304 return 1;
305 } else {
306 return StructuralElement::giveIPValue(answer, gp, type, tStep);
307 }
308}
309
310int
312// Returns the rotation matrix of the receiver of the size [6,6]
313// f(local) = T * f(global)
314{
316
317 answer.resize(6, 6);
318 answer.zero();
319
320 for ( int i = 1; i <= 3; i++ ) {
321 answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix.at(1, i);
322 answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix.at(2, i);
323 answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix.at(3, i);
324 }
325
326 return 1;
327}
328
329void
330CCTPlate3d::computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
331// Computes numerically the load vector of the receiver due to the body loads, at tStep.
332// load is assumed to be in global cs.
333// load vector is then transformed to coordinate system in each node.
334// (should be global coordinate system, but there may be defined
335// different coordinate system in each node)
336{
337 double dens, dV, load;
338 FloatArray force;
339 FloatMatrix T;
340
341 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
342 OOFEM_ERROR("unknown load type");
343 }
344
345 GaussIntegrationRule irule(1, this, 1, 5);
346 irule.SetUpPointsOnTriangle(1, _2dPlate);
347
348 // note: force is assumed to be in global coordinate system, 6 components
349 forLoad->computeComponentArrayAt(force, tStep, mode);
350 // get it in local c.s.
352 force.rotatedWith(T, 'n');
353
354 if ( force.giveSize() ) {
355 GaussPoint *gp = irule.getIntegrationPoint(0);
356
357 dens = this->giveStructuralCrossSection()->give('d', gp); // constant density assumed
358 dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp); // constant thickness assumed
359
360 answer.resize(9);
361 answer.zero();
362
363 load = force.at(3) * dens * dV / 3.0;
364 answer.at(1) = load;
365 answer.at(4) = load;
366 answer.at(7) = load;
367
368 load = force.at(4) * dens * dV / 3.0;
369 answer.at(2) = load;
370 answer.at(5) = load;
371 answer.at(8) = load;
372
373 load = force.at(5) * dens * dV / 3.0;
374 answer.at(3) = load;
375 answer.at(6) = load;
376 answer.at(9) = load;
377
378 // transform result from global cs to local element cs.
379 //if ( this->computeGtoLRotationMatrix(T) ) {
380 // answer.rotatedWith(T, 'n');
381 //}
382 } else {
383 answer.clear(); // nil resultant
384 }
385}
386
387void
389{
390 FloatMatrix ne;
392
393 answer.resize(6, 18);
394 answer.zero();
395 int ri[] = {
396 2, 3, 4
397 };
398 int ci[] = {
399 2, 3, 4, 8, 9, 10, 14, 15, 16
400 };
401
402 for ( int i = 0; i < 3; i++ ) {
403 for ( int j = 0; j < 9; j++ ) {
404 answer(ri [ i ], ci [ j ]) = ne(i, j);
405 }
406 }
407}
408
409void
411{
412 answer.resize(18);
413 answer.zero();
414 if ( iSurf == 1 ) {
415 answer.at(3) = 1; // node 1
416 answer.at(4) = 2;
417 answer.at(5) = 3;
418
419 answer.at(9) = 4; // node 2
420 answer.at(10) = 5;
421 answer.at(11) = 6;
422
423 answer.at(15) = 7; // node 3
424 answer.at(16) = 8;
425 answer.at(17) = 9;
426 } else {
427 OOFEM_ERROR("wrong surface number");
428 }
429}
430
431
432double
434{
435 return this->computeVolumeAround(gp);
436}
437
438
439int
441{
442 return 0;
443}
444
445void
447// Performs end-of-step operations.
448{
449 FloatArray v;
450
451 fprintf(file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
452
453 for ( int i = 0; i < ( int ) integrationRulesArray.size(); i++ ) {
454 for ( GaussPoint *gp: * integrationRulesArray [ i ] ) {
455 fprintf(file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
456
457 this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
458 fprintf(file, " strains ");
459 for ( auto &val : v ) {
460 fprintf(file, " %.4e", val);
461 }
462
463 this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
464 fprintf(file, "\n curvatures ");
465 for ( auto &val : v ) {
466 fprintf(file, " %.4e", val);
467 }
468
469 // Forces - Moments
470 this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
471 fprintf(file, "\n stresses ");
472 for ( auto &val : v ) {
473 fprintf(file, " %.4e", val);
474 }
475
476 this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
477 fprintf(file, "\n moments ");
478 for ( auto &val : v ) {
479 fprintf(file, " %.4e", val);
480 }
481
482 fprintf(file, "\n");
483 }
484 }
485}
486} // end namespace oofem
#define REGISTER_Element(class)
void giveLocalCoordinates(FloatArray &answer, const FloatArray &global)
Definition cct3d.C:55
FloatMatrix GtoLRotationMatrix
Definition cct3d.h:79
void giveSurfaceDofMapping(IntArray &answer, int iSurf) const override
Definition cct3d.C:410
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
Definition cct3d.C:209
int computeLoadGToLRotationMtrx(FloatMatrix &answer) override
Definition cct3d.C:311
void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Definition cct3d.C:388
virtual const FloatMatrix * computeGtoLRotationMatrix()
Definition cct3d.C:145
int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) override
Definition cct3d.C:130
void printOutputAt(FILE *file, TimeStep *tStep) override
Definition cct3d.C:446
double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf) override
Definition cct3d.C:433
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition cct3d.C:96
void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3) override
Definition cct3d.C:71
bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords) override
Definition cct3d.C:103
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition cct3d.C:264
int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp) override
Definition cct3d.C:440
CCTPlate3d(int n, Domain *d)
Definition cct3d.C:50
void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode) override
Definition cct3d.C:330
CCTPlate(int n, Domain *d)
Definition cct.C:61
Element_Geometry_Type giveGeometryType() const override
Definition cct.h:80
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition cct.C:213
double computeVolumeAround(GaussPoint *gp) override
Definition cct.C:347
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
double giveCoordinate(int i) const
Definition dofmanager.h:383
Node * giveNode(int i) const
Definition element.h:629
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
CrossSection * giveCrossSection()
Definition element.C:534
int giveLabel() const
Definition element.h:1125
int global2local(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
Definition fei2dtrlin.C:130
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 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
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.
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
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.
#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