OOFEM 3.0
Loading...
Searching...
No Matches
rershell.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 "floatarray.h"
44#include "intarray.h"
45#include "domain.h"
46#include "load.h"
47#include "mathfem.h"
48#include "classfactory.h"
49
50#ifdef __OOFEG
51 #include "oofeggraphiccontext.h"
52 #include "connectivitytable.h"
53#endif
54
55#include <cstdlib>
56
57namespace oofem {
59
60RerShell :: RerShell(int n, Domain *aDomain) :
61 CCTPlate3d(n, aDomain)
62{
63 Rx = 1.e+40;
64 Ry = 1.e+40;
65 Rxy = 1.e+40;
66
68}
69
70
72RerShell :: giveInterface(InterfaceType interface)
73{
74 if ( interface == LayeredCrossSectionInterfaceType ) {
75 return static_cast< LayeredCrossSectionInterface * >(this);
76 } else if ( interface == ZZNodalRecoveryModelInterfaceType ) {
77 return static_cast< ZZNodalRecoveryModelInterface * >(this);
78 } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
79 return static_cast< NodalAveragingRecoveryModelInterface * >(this);
80 }
81
82 return NULL;
83}
84
85
86void
87RerShell :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
88// Returns the [8x18] strain-displacement matrix {B} of the receiver,
89// evaluated at gp.
90{
91 FloatArray b(3), c(3), nodeCoords;
92 double x1, x2, x3, y1, y2, y3, area;
93
94 this->giveLocalCoordinates( nodeCoords, this->giveNode(1)->giveCoordinates() );
95 x1 = nodeCoords.at(1);
96 y1 = nodeCoords.at(2);
97
98 this->giveLocalCoordinates( nodeCoords, this->giveNode(2)->giveCoordinates() );
99 x2 = nodeCoords.at(1);
100 y2 = nodeCoords.at(2);
101
102 this->giveLocalCoordinates( nodeCoords, this->giveNode(3)->giveCoordinates() );
103 x3 = nodeCoords.at(1);
104 y3 = nodeCoords.at(2);
105
106
107 b.at(1) = y2 - y3;
108 b.at(2) = y3 - y1;
109 b.at(3) = y1 - y2;
110
111 c.at(1) = x3 - x2;
112 c.at(2) = x1 - x3;
113 c.at(3) = x2 - x1;
114
115 area = this->giveArea();
116
117 answer.resize(8, 18);
118 answer.zero();
119
120 for ( int i = 1; i <= 3; i++ ) {
121 int j = i + 1 - i / 3 * 3;
122 int k = j + 1 - j / 3 * 3;
123
124 int ii = 6 * ( i - 1 );
125 int ki = 6 * ( k - 1 );
126
127 answer.at(1, ii + 1) = b.at(i); // Eps_X
128 answer.at(1, ki + 4) = -( b.at(i) - b.at(j) ) / Rx * area / 12.;
129 answer.at(1, ki + 5) = -( c.at(i) - c.at(j) ) / Ry * area / 12.;
130
131 answer.at(2, ii + 2) = c.at(i); // Eps_y
132 answer.at(2, ki + 4) = -( b.at(i) - b.at(j) ) / Ry * area / 12.;
133 answer.at(2, ki + 5) = -( c.at(i) - c.at(j) ) / Ry * area / 12.;
134
135 answer.at(3, ii + 1) = c.at(i); // Gamma_xy;
136 answer.at(3, ii + 2) = b.at(i);
137 answer.at(3, ki + 4) = -( b.at(i) - b.at(j) ) / Rxy * area / 6.;
138 answer.at(3, ki + 5) = -( c.at(i) - c.at(j) ) / Rxy * area / 6.;
139
140 answer.at(4, ii + 5) = b.at(i); // Kappa_x
141
142 answer.at(5, ii + 4) = -c.at(i); // Kappa_y
143
144 answer.at(6, ii + 4) = -b.at(i); // Kappa_xy
145 answer.at(6, ii + 5) = c.at(i);
146
147 answer.at(7, ii + 3) = b.at(i); // Gamma_zx
148 answer.at(7, ii + 5) += area * ( 2. / 3. );
149 answer.at(7, ki + 5) += ( -2.0 * area + b.at(k) * ( c.at(i) - c.at(j) ) ) / 6.0;
150 answer.at(7, ki + 4) = b.at(k) * ( b.at(i) - b.at(j) ) / 6.0;
151
152 answer.at(8, ii + 3) = c.at(i); // Gamma_zy
153 answer.at(8, ii + 4) += area * ( -2.0 / 3.0 );
154 answer.at(8, ki + 4) += ( +2.0 * area + c.at(k) * ( b.at(i) - b.at(j) ) ) / 6.0;
155 answer.at(8, ki + 5) = c.at(k) * ( c.at(i) - c.at(j) ) / 6.0;
156 }
157
158 answer.times( 1. / ( 2. * area ) );
159}
160
161
162void
163RerShell :: computeGaussPoints()
164// Sets up the array containing the four Gauss points of the receiver.
165{
166 if ( integrationRulesArray.size() == 0 ) {
167 integrationRulesArray.resize( 1 );
168 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 8);
169 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
170 }
171}
172
173
174void
175RerShell :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
176// Returns the displacement interpolation matrix {N} of the receiver,
177// evaluated at gp.
178{
179 double x1, x2, x3, y1, y2, y3, l1, l2, l3, b1, b2, b3, c1, c2, c3;
180 FloatArray nodeCoords;
181
182 l1 = iLocCoord.at(1);
183 l2 = iLocCoord.at(2);
184 l3 = 1.0 - l1 - l2;
185
186 this->giveLocalCoordinates( nodeCoords, this->giveNode(1)->giveCoordinates() );
187 x1 = nodeCoords.at(1);
188 y1 = nodeCoords.at(2);
189
190 this->giveLocalCoordinates( nodeCoords, this->giveNode(2)->giveCoordinates() );
191 x2 = nodeCoords.at(1);
192 y2 = nodeCoords.at(2);
193
194 this->giveLocalCoordinates( nodeCoords, this->giveNode(3)->giveCoordinates() );
195 x3 = nodeCoords.at(1);
196 y3 = nodeCoords.at(2);
197
198 /*
199 * x1 = this -> giveNode(1) -> giveCoordinate(1);
200 * x2 = this -> giveNode(2) -> giveCoordinate(1);
201 * x3 = this -> giveNode(3) -> giveCoordinate(1);
202 *
203 * y1 = this -> giveNode(1) -> giveCoordinate(2);
204 * y2 = this -> giveNode(2) -> giveCoordinate(2);
205 * y3 = this -> giveNode(3) -> giveCoordinate(2);
206 */
207 b1 = y2 - y3;
208 b2 = y3 - y1;
209 b3 = y1 - y2;
210
211 c1 = x3 - x2;
212 c2 = x1 - x3;
213 c3 = x2 - x1;
214
215 answer.resize(5, 18);
216 answer.zero();
217
218
219 answer.at(1, 1) = l1;
220 answer.at(1, 7) = l2;
221 answer.at(1, 13) = l3;
222
223 answer.at(2, 2) = l1;
224 answer.at(2, 8) = l2;
225 answer.at(2, 14) = l3;
226
227 answer.at(3, 3) = l1;
228 answer.at(3, 4) = -( l1 * l2 * b3 - l1 * l3 * b2 ) * 0.5;
229 answer.at(3, 5) = -( l1 * l2 * c3 - l1 * l3 * c2 ) * 0.5;
230 answer.at(3, 9) = l2;
231 answer.at(3, 10) = -( l2 * l3 * b1 - l1 * l2 * b3 ) * 0.5;
232 answer.at(3, 11) = -( l2 * l3 * c1 - l1 * l2 * c3 ) * 0.5;
233 answer.at(3, 15) = l3;
234 answer.at(3, 16) = -( l3 * l1 * b2 - l2 * l3 * b1 ) * 0.5;
235 answer.at(3, 17) = -( l3 * l1 * c2 - l2 * l3 * c1 ) * 0.5;
236
237 answer.at(4, 4) = l1;
238 answer.at(4, 10) = l2;
239 answer.at(4, 16) = l3;
240
241 answer.at(5, 5) = l1;
242 answer.at(5, 11) = l2;
243 answer.at(5, 17) = l3;
244}
245
246
247double
248RerShell :: giveArea()
249// returns the area occupied by the receiver
250{
251 FloatArray nodeCoords;
252 if ( area > 0 ) {
253 return area; // check if previously computed
254 }
255
256 double x1, x2, x3, y1, y2, y3;
257
258 this->giveLocalCoordinates( nodeCoords, this->giveNode(1)->giveCoordinates() );
259 x1 = nodeCoords.at(1);
260 y1 = nodeCoords.at(2);
261
262 this->giveLocalCoordinates( nodeCoords, this->giveNode(2)->giveCoordinates() );
263 x2 = nodeCoords.at(1);
264 y2 = nodeCoords.at(2);
265
266 this->giveLocalCoordinates( nodeCoords, this->giveNode(3)->giveCoordinates() );
267 x3 = nodeCoords.at(1);
268 y3 = nodeCoords.at(2);
269
270 return ( area = 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 ) );
271}
272
273
274void
275RerShell :: initializeFrom(InputRecord &ir, int priority)
276{
277 StructuralElement :: initializeFrom(ir, priority);
279}
280
281
282void
283RerShell :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
284{
285 answer = this->giveStructuralCrossSection()->give3dShellStiffMtrx(rMode, gp, tStep);
286}
287
288void
289RerShell :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
290{
291 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Shell(strain, gp, tStep);
292}
293
294void
295RerShell :: computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
296// Returns the lumped mass matrix of the receiver.
297{
298 GaussPoint *gp;
299 double dV, mss1;
300
301 answer.resize(18, 18);
302 answer.zero();
303 gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
304
305 dV = this->computeVolumeAround(gp);
306 mss1 = dV * this->giveCrossSection()->give(CS_Thickness, gp) * this->giveStructuralCrossSection()->give('d', gp) / 3.;
307
308 answer.at(1, 1) = mss1;
309 answer.at(2, 2) = mss1;
310 answer.at(3, 3) = mss1;
311
312 answer.at(7, 7) = mss1;
313 answer.at(8, 8) = mss1;
314 answer.at(9, 9) = mss1;
315
316 answer.at(13, 13) = mss1;
317 answer.at(14, 14) = mss1;
318 answer.at(15, 15) = mss1;
319}
320
321
322void
323RerShell :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
324// Computes numerically the load vector of the receiver due to the body
325// loads, at tStep. - no support for momentum bodyload
326{
327 double dens, dV, load;
328 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
329 FloatArray f;
330 FloatMatrix T;
331
332
333 forLoad->computeComponentArrayAt(f, tStep, mode);
334
335 if ( f.giveSize() == 0 ) {
336 answer.clear();
337 return; // nil resultant
338 } else {
339 dens = this->giveStructuralCrossSection()->give('d', gp);
340 dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp);
341
342 answer.resize(18);
343
344 load = f.at(1) * dens * dV / 3.0;
345 answer.at(1) = load;
346 answer.at(7) = load;
347 answer.at(13) = load;
348
349 load = f.at(2) * dens * dV / 3.0;
350 answer.at(2) = load;
351 answer.at(8) = load;
352 answer.at(14) = load;
353
354 load = f.at(3) * dens * dV / 3.0;
355 answer.at(3) = load;
356 answer.at(9) = load;
357 answer.at(15) = load;
358
359 // transform result from global cs to local element cs.
360 if ( this->computeGtoLRotationMatrix(T) ) {
361 answer.rotatedWith(T, 'n');
362 }
363
364 return;
365 }
366}
367
368int
369RerShell :: giveLocalCoordinateSystem(FloatMatrix &answer)
370//
371// returns a unit vectors of local coordinate system at element
372// stored rowwise (mainly used by some materials with ortho and anisotrophy)
373//
374{
376 answer = GtoLRotationMatrix;
377 return 1;
378}
379
380
381//converts global coordinates to local planar area coordinates, does not return a coordinate in the thickness direction, but
382//does check that the point is in the element thickness
383#define POINT_TOL 1.e-3
384bool
385RerShell :: computeLocalCoordinates(FloatArray &answer, const FloatArray &coords)
386{
387 //set size of return value to 3 area coordinates
388 answer.resize(3);
389
390 //rotate the input point Coordinate System into the element CS
391 FloatArray inputCoords_ElCS;
392 this->giveLocalCoordinates( inputCoords_ElCS, coords );
393
394 //Nodes are defined in the global CS, so they also need to be rotated into the element CS, therefore get the node points and
395 //rotate them into the element CS
396 FloatArray nodeCoords;
397 double x1, x2, x3, y1, y2, y3, z1, z2, z3;
398
399 this->giveLocalCoordinates( nodeCoords, this->giveNode(1)->giveCoordinates() );
400 x1 = nodeCoords.at(1);
401 y1 = nodeCoords.at(2);
402 z1 = nodeCoords.at(3);
403
404 this->giveLocalCoordinates( nodeCoords, this->giveNode(2)->giveCoordinates() );
405 x2 = nodeCoords.at(1);
406 y2 = nodeCoords.at(2);
407 z2 = nodeCoords.at(3);
408
409 this->giveLocalCoordinates( nodeCoords, this->giveNode(3)->giveCoordinates() );
410 x3 = nodeCoords.at(1);
411 y3 = nodeCoords.at(2);
412 z3 = nodeCoords.at(3);
413
414 //Compute the area coordinates corresponding to this point
415 double area;
416 area = 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 );
417
418 answer.at(1) = ( ( x2 * y3 - x3 * y2 ) + ( y2 - y3 ) * inputCoords_ElCS.at(1) + ( x3 - x2 ) * inputCoords_ElCS.at(2) ) / 2. / area;
419 answer.at(2) = ( ( x3 * y1 - x1 * y3 ) + ( y3 - y1 ) * inputCoords_ElCS.at(1) + ( x1 - x3 ) * inputCoords_ElCS.at(2) ) / 2. / area;
420 answer.at(3) = ( ( x1 * y2 - x2 * y1 ) + ( y1 - y2 ) * inputCoords_ElCS.at(1) + ( x2 - x1 ) * inputCoords_ElCS.at(2) ) / 2. / area;
421
422 //get midplane location at this point
423 double midplZ;
424 midplZ = z1 * answer.at(1) + z2 *answer.at(2) + z3 *answer.at(3);
425
426 //check that the z is within the element
428 GaussPoint _gp(NULL, 1, answer, 1.0, _2dPlate);
429
430 double elthick;
431
432 elthick = cs->give(CS_Thickness, & _gp);
433
434 if ( elthick / 2.0 + midplZ - fabs( inputCoords_ElCS.at(3) ) < -POINT_TOL ) {
435 answer.zero();
436 return false;
437 }
438
439 //check that the point is in the element and set flag
440 for ( int i = 1; i <= 3; i++ ) {
441 if ( answer.at(i) < ( 0. - POINT_TOL ) ) {
442 return false;
443 }
444
445 if ( answer.at(i) > ( 1. + POINT_TOL ) ) {
446 return false;
447 }
448 }
449
450 return true;
451}
452
453
454bool
455RerShell :: computeGtoLRotationMatrix(FloatMatrix &answer)
456// Returns the rotation matrix of the receiver of the size [18,18]
457// r(local) = T * r(global)
458//
459// local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
460//
461// e1' : [N2-N1] Ni - means i - th node
462// help : [N3-N1]
463// e3' : e1' x help
464// e2' : e3' x e1'
465{
466 double val;
467
469
470
471 answer.resize(18, 18);
472 answer.zero();
473
474 for ( int i = 1; i <= 3; i++ ) {
475 for ( int j = 1; j <= 3; j++ ) {
476 val = GtoLRotationMatrix.at(i, j);
477 answer.at(i, j) = val;
478 answer.at(i + 3, j + 3) = val;
479 answer.at(i + 6, j + 6) = val;
480 answer.at(i + 9, j + 9) = val;
481 answer.at(i + 12, j + 12) = val;
482 answer.at(i + 15, j + 15) = val;
483 }
484 }
485
486 return 1;
487}
488
489
490void
491RerShell :: giveLocalCoordinates(FloatArray &answer, const FloatArray &global)
492{
493 if ( global.giveSize() != 3 ) {
494 OOFEM_ERROR("cannot transform coordinates- size mismatch");
495 }
496
497 // first ensure that receiver's GtoLRotationMatrix[3,3]
498 // is defined
499
501 answer.beProductOf(GtoLRotationMatrix, global);
502}
503
504
505void
506RerShell :: giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
507//
508// returns characteristic tensor of the receiver at given gp and tStep
509//
510{
511 answer.resize(3, 3);
512 answer.zero();
513
514 if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
515 FloatArray stress, strain;
516 this->computeStrainVector(strain, gp, tStep);
517 this->computeStressVector(stress, strain, gp, tStep);
518 answer.at(1, 1) = stress.at(1);
519 answer.at(2, 2) = stress.at(2);
520 answer.at(1, 2) = stress.at(3);
521 answer.at(2, 1) = stress.at(3);
522 answer.at(1, 3) = stress.at(7);
523 answer.at(3, 1) = stress.at(7);
524 answer.at(2, 3) = stress.at(8);
525 answer.at(3, 2) = stress.at(8);
526 } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {
527 FloatArray stress, strain;
528 this->computeStrainVector(strain, gp, tStep);
529 this->computeStressVector(stress, strain, gp, tStep);
530 answer.at(1, 1) = stress.at(4);
531 answer.at(2, 2) = stress.at(5);
532 answer.at(1, 2) = stress.at(6);
533 answer.at(2, 1) = stress.at(6);
534 } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
535 FloatArray strain;
536 this->computeStrainVector(strain, gp, tStep);
537
538 answer.at(1, 1) = strain.at(1);
539 answer.at(2, 2) = strain.at(2);
540 answer.at(1, 2) = strain.at(3) / 2.;
541 answer.at(2, 1) = strain.at(3) / 2.;
542 answer.at(1, 3) = strain.at(7) / 2.;
543 answer.at(3, 1) = strain.at(7) / 2.;
544 answer.at(2, 3) = strain.at(8) / 2.;
545 answer.at(3, 2) = strain.at(8) / 2.;
546 } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {
547 FloatArray curv;
548 this->computeStrainVector(curv, gp, tStep);
549 answer.at(1, 1) = curv.at(4);
550 answer.at(2, 2) = curv.at(5);
551 answer.at(1, 2) = curv.at(6) / 2.;
552 answer.at(2, 1) = curv.at(6) / 2.;
553 } else {
554 OOFEM_ERROR("unsupported tensor mode");
555 }
556
557 if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
558 ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
561 }
562}
563
564
565void
566RerShell :: computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain,
567 GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
568//
569// returns full 3d strain vector of given layer (whose z-coordinate from center-line is
570// stored in slaveGp) for given tStep
571//
572{
573 double layerZeta, layerZCoord, top, bottom;
574
575 top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
576 bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
577 layerZeta = slaveGp->giveNaturalCoordinate(3);
578 layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
579
580 answer.resize(5); // {Exx,Eyy,GMyz,GMzx,GMxy}
581
582 answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(4) * layerZCoord;
583 answer.at(2) = masterGpStrain.at(2) + masterGpStrain.at(5) * layerZCoord;
584 answer.at(5) = masterGpStrain.at(3) + masterGpStrain.at(6) * layerZCoord;
585 answer.at(3) = masterGpStrain.at(8);
586 answer.at(4) = masterGpStrain.at(7);
587}
588
589
590void
591RerShell :: printOutputAt(FILE *file, TimeStep *tStep)
592{
593 FloatArray v;
594
595 fprintf(file, "element %d (%8d):\n", this->giveLabel(), number);
596
597 for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
598
599 fprintf( file, " GP 1.%d :", gp->giveNumber() );
600 this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
601 fprintf(file, " strains ");
602 // eps_x, eps_y, eps_z, eps_yz, eps_xz, eps_xy (global)
603 fprintf( file,
604 " %.4e %.4e %.4e %.4e %.4e %.4e ",
605 v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
606
607 this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
608 fprintf(file, "\n curvatures ");
609 // k_x, k_y, k_z, k_yz, k_xz, k_xy (global)
610 fprintf( file,
611 " %.4e %.4e %.4e %.4e %.4e %.4e ",
612 v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
613
614 // Forces - Moments
615 this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
616 fprintf(file, "\n stresses ");
617 // n_x, n_y, n_z, v_yz, v_xz, v_xy (global)
618 fprintf( file,
619 " %.4e %.4e %.4e %.4e %.4e %.4e ",
620 v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
621
622 this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
623 fprintf(file, "\n moments ");
624 // m_x, m_y, m_z, m_yz, m_xz, m_xy (global)
625 fprintf( file,
626 " %.4e %.4e %.4e %.4e %.4e %.4e ",
627 v.at(1), v.at(2), v.at(3), v.at(4), v.at(5), v.at(6) );
628
629 fprintf(file, "\n");
630 }
631}
632
633
634void
635RerShell :: giveDofManDofIDMask(int inode, IntArray &answer) const
636{
637 answer = {D_u, D_v, D_w, R_u, R_v, R_w};
638}
639
640
641void
642RerShell :: NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node,
643 InternalStateType type, TimeStep *tStep)
644{
645 this->giveIPValue(answer, integrationRulesArray [ 0 ]->getIntegrationPoint(0), type, tStep);
646}
647
648
649int
650RerShell :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
651{
652 FloatMatrix globTensor;
653 CharTensor cht;
654
655 answer.resize(6);
656
657 if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
658 if ( type == IST_CurvatureTensor ) {
660 } else {
661 cht = GlobalStrainTensor;
662 }
663
664 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
665
666 answer.at(1) = globTensor.at(1, 1); //xx
667 answer.at(2) = globTensor.at(2, 2); //yy
668 answer.at(3) = globTensor.at(3, 3); //zz
669 answer.at(4) = 2*globTensor.at(2, 3); //yz
670 answer.at(5) = 2*globTensor.at(1, 3); //xz
671 answer.at(6) = 2*globTensor.at(1, 2); //xy
672
673 return 1;
674 } else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
675 if ( type == IST_ShellMomentTensor ) {
676 cht = GlobalMomentTensor;
677 } else {
678 cht = GlobalForceTensor;
679 }
680
681 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
682
683 answer.at(1) = globTensor.at(1, 1); //xx
684 answer.at(2) = globTensor.at(2, 2); //yy
685 answer.at(3) = globTensor.at(3, 3); //zz
686 answer.at(4) = globTensor.at(2, 3); //yz
687 answer.at(5) = globTensor.at(1, 3); //xz
688 answer.at(6) = globTensor.at(1, 2); //xy
689
690 return 1;
691 } else {
692 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
693 }
694}
695
696
697#ifdef __OOFEG
698/*
699 * void
700 * CCTPlate :: drawRawGeometry ()
701 * {
702 * WCRec p[3];
703 * GraphicObj *go;
704 *
705 * EASValsSetLineWidth(RAW_GEOMETRY_WIDTH);
706 * EASValsSetColor(elementColor);
707 * EASValsSetLayer(RAW_GEOMETRY_LAYER);
708 * p[0].x = (FPNum) this->giveNode(1)->giveCoordinate(1);
709 * p[0].y = (FPNum) this->giveNode(1)->giveCoordinate(2);
710 * p[0].z = (FPNum) this->giveNode(1)->giveCoordinate(3);
711 * p[1].x = (FPNum) this->giveNode(2)->giveCoordinate(1);
712 * p[1].y = (FPNum) this->giveNode(2)->giveCoordinate(2);
713 * p[1].z = (FPNum) this->giveNode(2)->giveCoordinate(3);
714 * p[2].x = (FPNum) this->giveNode(3)->giveCoordinate(1);
715 * p[2].y = (FPNum) this->giveNode(3)->giveCoordinate(2);
716 * p[2].z = (FPNum) this->giveNode(3)->giveCoordinate(3);
717 *
718 * go = CreateTriangle3D(p);
719 * EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
720 * EGAttachObject(go, (EObjectP) this);
721 * EMAddGraphicsToModel(ESIModel(), go);
722 * }
723 *
724 * void
725 * CCTPlate :: drawDeformedGeometry (UnknownType defMode)
726 * {
727 * WCRec p[3];
728 * GraphicObj *go;
729 *
730 * EASValsSetLineWidth(DEFORMED_GEOMETRY_WIDTH);
731 * EASValsSetColor(deformedElementColor);
732 * EASValsSetLayer(DEFORMED_GEOMETRY_LAYER);
733 * p[0].x = (FPNum) this->giveNode(1)->giveUpdatedCoordinate(1,tStep,defScale);
734 * p[0].y = (FPNum) this->giveNode(1)->giveUpdatedCoordinate(2,tStep,defScale);
735 * p[0].z = (FPNum) this->giveNode(1)->giveUpdatedCoordinate(3,tStep,defScale);
736 * p[1].x = (FPNum) this->giveNode(2)->giveUpdatedCoordinate(1,tStep,defScale);
737 * p[1].y = (FPNum) this->giveNode(2)->giveUpdatedCoordinate(2,tStep,defScale);
738 * p[1].z = (FPNum) this->giveNode(2)->giveUpdatedCoordinate(3,tStep,defScale);
739 * p[2].x = (FPNum) this->giveNode(3)->giveUpdatedCoordinate(1,tStep,defScale);
740 * p[2].y = (FPNum) this->giveNode(3)->giveUpdatedCoordinate(2,tStep,defScale);
741 * p[2].z = (FPNum) this->giveNode(3)->giveUpdatedCoordinate(3,tStep,defScale);
742 *
743 * go = CreateTriangle3D(p);
744 * EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
745 * EMAddGraphicsToModel(ESIModel(), go);
746 * }
747 *
748 */
749
750
751//void
752//RerShell :: drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
753//{}
754
755
756
757/*
758 * int
759 * RerShell :: giveInternalStateAtNode (FloatArray& answer, InternalValueRequest& r, int node, TimeStep* tStep)
760 *
761 * //
762 * // see eleemnt.h for description
763 * //
764 * {
765 * opravit;
766 *
767 * DrawMode mode = gc.getDrawMode();
768 *
769 * if (!gc.testElementGraphicActivity(this)) return 0;
770 *
771 * const FloatArray* nodval;
772 * int result = this->giveDomain()->giveSmoother()->giveNodalVector(nodval, this->giveNode(inode)->giveNumber(),
773 * this->giveDomain()->giveSmoother()->giveElementRegion(this));
774 * if (result) {
775 * if (mode == sxForce ) {
776 ***val = nodval->at(1);
777 * return 1;
778 * } else if (mode == syForce) {
779 ***val = nodval->at(2);
780 * return 1;
781 * } else if (mode == sxyForce) {
782 ***val = nodval->at(3);
783 * return 1;
784 * } else if (mode == mxForce ) {
785 ***val = nodval->at(4);
786 * return 1;
787 * } else if (mode == myForce) {
788 ***val = nodval->at(5);
789 * return 1;
790 * } else if (mode == mxyForce) {
791 ***val = nodval->at(6);
792 * return 1;
793 * } else if (mode == szxForce ) {
794 ***val = nodval->at(7);
795 * return 1;
796 * } else if (mode == syzForce) {
797 ***val = nodval->at(8);
798 * return 1;
799 * } else return 0;
800 * }
801 * return 0;
802 * }
803 */
804#endif
805} // end namespace oofem
#define REGISTER_Element(class)
FloatMatrix GtoLRotationMatrix
Definition cct3d.h:79
CCTPlate3d(int n, Domain *d)
Definition cct3d.C:50
double area
Definition cct.h:67
double computeVolumeAround(GaussPoint *gp) override
Definition cct.C:347
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
Node * giveNode(int i) const
Definition element.h:629
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
CrossSection * giveCrossSection()
Definition element.C:534
int giveLabel() const
Definition element.h:1125
int number
Component number.
Definition femcmpnn.h:77
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 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 times(double f)
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
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition gausspoint.h:136
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition load.C:84
double giveArea()
Definition rershell.C:248
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
Definition rershell.C:506
void giveLocalCoordinates(FloatArray &answer, const FloatArray &global)
Definition rershell.C:491
const FloatMatrix * computeGtoLRotationMatrix() override
Definition rershell.h:85
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition rershell.C:650
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition rershell.C:289
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
#define OOFEM_ERROR(...)
Definition error.h:79
#define POINT_TOL
@ CS_BottomZCoord
Bottom z coordinate.
@ CS_Thickness
Thickness.
@ CS_TopZCoord
Top z coordinate.
@ ZZNodalRecoveryModelInterfaceType
@ LayeredCrossSectionInterfaceType
@ NodalAveragingRecoveryModelInterfaceType

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