OOFEM 3.0
Loading...
Searching...
No Matches
dkt3d.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
50DKTPlate3d::DKTPlate3d(int n, Domain *aDomain) : DKTPlate(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
71DKTPlate3d::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
102const FloatMatrix *
104// Returns the rotation matrix of the receiver of the size [3,3]
105// coords(local) = T * coords(global)
106//
107// local coordinate (described by vector triplet e1',e2',e3') is defined as follows:
108//
109// e1' : [N2-N1] Ni - means i - th node
110// help : [N3-N1]
111// e3' : e1' x help
112// e2' : e3' x e1'
113{
114 if ( !GtoLRotationMatrix.isNotEmpty() ) {
115 FloatArray e1, e2, e3, help;
116
117 // compute e1' = [N2-N1] and help = [N3-N1]
118 e1.beDifferenceOf( this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
119 help.beDifferenceOf( this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
120
121 // let us normalize e1'
122 e1.normalize();
123
124 // compute e3' : vector product of e1' x help
125 e3.beVectorProductOf(e1, help);
126 // let us normalize
127 e3.normalize();
128
129 // now from e3' x e1' compute e2'
130 e2.beVectorProductOf(e3, e1);
131
132 //
133 GtoLRotationMatrix.resize(3, 3);
134
135 for ( int i = 1; i <= 3; i++ ) {
136 GtoLRotationMatrix.at(1, i) = e1.at(i);
137 GtoLRotationMatrix.at(2, i) = e2.at(i);
138 GtoLRotationMatrix.at(3, i) = e3.at(i);
139 }
140 }
141
142 return & GtoLRotationMatrix;
143}
144
145
146bool
148// Returns the rotation matrix of the receiver of the size [9,18]
149// r(local) = T * r(global)
150// for one node (r written transposed): {w,r1,r2} = T * {u,v,w,r1,r2,r3}
151{
153
154 answer.resize(9, 18);
155 answer.zero();
156
157 for ( int i = 1; i <= 3; i++ ) {
158 answer.at(1, i) = answer.at(1 + 3, i + 6) = answer.at(1 + 6, i + 12) = GtoLRotationMatrix.at(3, i);
159 answer.at(2, i + 3) = answer.at(2 + 3, i + 3 + 6) = answer.at(2 + 6, i + 3 + 12) = GtoLRotationMatrix.at(1, i);
160 answer.at(3, i + 3) = answer.at(3 + 3, i + 3 + 6) = answer.at(3 + 6, i + 3 + 12) = GtoLRotationMatrix.at(2, i);
161 }
162
163 return 1;
164}
165
166void
168// returns characteristic tensor of the receiver at given gp and tStep
169// strain vector = (Kappa_x, Kappa_y, Kappa_xy, Gamma_zx, Gamma_zy)
170{
171 FloatArray charVect;
173
174 answer.resize(3, 3);
175 answer.zero();
176
177 if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) {
178 //this->computeStressVector(charVect, gp, tStep);
179 charVect = ms->giveStressVector();
180
181 answer.at(1, 3) = charVect.at(4);
182 answer.at(3, 1) = charVect.at(4);
183 answer.at(2, 3) = charVect.at(5);
184 answer.at(3, 2) = charVect.at(5);
185 } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) {
186 //this->computeStressVector(charVect, gp, tStep);
187 charVect = ms->giveStressVector();
188
189 answer.at(1, 1) = charVect.at(1);
190 answer.at(2, 2) = charVect.at(2);
191 answer.at(1, 2) = charVect.at(3);
192 answer.at(2, 1) = charVect.at(3);
193 } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) {
194 //this->computeStrainVector(charVect, gp, tStep);
195 charVect = ms->giveStrainVector();
196
197 answer.at(1, 3) = charVect.at(4) / 2.;
198 answer.at(3, 1) = charVect.at(4) / 2.;
199 answer.at(2, 3) = charVect.at(5) / 2.;
200 answer.at(3, 2) = charVect.at(5) / 2.;
201 } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) {
202 //this->computeStrainVector(charVect, gp, tStep);
203 charVect = ms->giveStrainVector();
204
205 answer.at(1, 1) = charVect.at(1);
206 answer.at(2, 2) = charVect.at(2);
207 answer.at(1, 2) = charVect.at(3) / 2.;
208 answer.at(2, 1) = charVect.at(3) / 2.;
209 } else {
210 OOFEM_ERROR("unsupported tensor mode");
211 }
212
213 if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) ||
214 ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) {
217 }
218}
219
220
221int
223{
224 FloatMatrix globTensor;
225 CharTensor cht;
226
227 answer.resize(6);
228
229 if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
230 if ( type == IST_CurvatureTensor ) {
232 } else {
233 cht = GlobalStrainTensor;
234 }
235
236 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
237
238 answer.at(1) = globTensor.at(1, 1); //xx
239 answer.at(2) = globTensor.at(2, 2); //yy
240 answer.at(3) = globTensor.at(3, 3); //zz
241 answer.at(4) = 2 * globTensor.at(2, 3); //yz
242 answer.at(5) = 2 * globTensor.at(1, 3); //xz
243 answer.at(6) = 2 * globTensor.at(1, 2); //xy
244
245 return 1;
246 } else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
247 if ( type == IST_ShellMomentTensor ) {
248 cht = GlobalMomentTensor;
249 } else {
250 cht = GlobalForceTensor;
251 }
252
253 this->giveCharacteristicTensor(globTensor, cht, gp, tStep);
254
255 answer.at(1) = globTensor.at(1, 1); //xx
256 answer.at(2) = globTensor.at(2, 2); //yy
257 answer.at(3) = globTensor.at(3, 3); //zz
258 answer.at(4) = globTensor.at(2, 3); //yz
259 answer.at(5) = globTensor.at(1, 3); //xz
260 answer.at(6) = globTensor.at(1, 2); //xy
261
262 return 1;
263 } else {
264 return StructuralElement::giveIPValue(answer, gp, type, tStep);
265 }
266}
267
268int
270// Returns the rotation matrix of the receiver of the size [6,6]
271// f(local) = T * f(global)
272{
274
275 answer.resize(6, 6);
276 answer.zero();
277
278 for ( int i = 1; i <= 3; i++ ) {
279 answer.at(1, i) = answer.at(4, i + 3) = GtoLRotationMatrix.at(1, i);
280 answer.at(2, i) = answer.at(5, i + 3) = GtoLRotationMatrix.at(2, i);
281 answer.at(3, i) = answer.at(6, i + 3) = GtoLRotationMatrix.at(3, i);
282 }
283
284 return 1;
285}
286
287void
289{
290 FloatMatrix ne;
292
293 answer.resize(6, 18);
294 answer.zero();
295 int ri[] = {
296 2, 3, 4
297 };
298 int ci[] = {
299 2, 3, 4, 8, 9, 10, 14, 15, 16
300 };
301
302 for ( int i = 0; i < 3; i++ ) {
303 for ( int j = 0; j < 9; j++ ) {
304 answer(ri [ i ], ci [ j ]) = ne(i, j);
305 }
306 }
307}
308
309void
311{
312 answer.resize(18);
313 answer.zero();
314 if ( iSurf == 1 ) {
315 answer.at(3) = 1; // node 1
316 answer.at(4) = 2;
317 answer.at(5) = 3;
318
319 answer.at(9) = 4; // node 2
320 answer.at(10) = 5;
321 answer.at(11) = 6;
322
323 answer.at(15) = 7; // node 3
324 answer.at(16) = 8;
325 answer.at(17) = 9;
326 } else {
327 OOFEM_ERROR("wrong surface number");
328 }
329}
330
331
332double
334{
335 return this->computeVolumeAround(gp);
336}
337
338
339int
341{
342 return 0;
343}
344
345void
347// Performs end-of-step operations.
348{
349 FloatArray v;
350
351 fprintf(file, "element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
352
353 for ( int i = 0; i < ( int ) integrationRulesArray.size(); i++ ) {
354 for ( GaussPoint *gp: * integrationRulesArray [ i ] ) {
355 fprintf(file, " GP %2d.%-2d :", i + 1, gp->giveNumber() );
356
357 this->giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
358 fprintf(file, " strains ");
359 // eps_x, eps_y, eps_z, eps_yz, eps_xz, eps_xy (global)
360 for ( auto &val : v ) {
361 fprintf(file, " %.4e", val);
362 }
363
364 this->giveIPValue(v, gp, IST_CurvatureTensor, tStep);
365 fprintf(file, "\n curvatures ");
366 for ( auto &val : v ) {
367 fprintf(file, " %.4e", val);
368 }
369
370 // Forces - Moments
371 this->giveIPValue(v, gp, IST_ShellForceTensor, tStep);
372 fprintf(file, "\n stresses ");
373 for ( auto &val : v ) {
374 fprintf(file, " %.4e", val);
375 }
376
377 this->giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
378 fprintf(file, "\n moments ");
379 for ( auto &val : v ) {
380 fprintf(file, " %.4e", val);
381 }
382
383 fprintf(file, "\n");
384 }
385 }
386}
387
388
389void
391{
392 /*
393 * provides dof mapping of local edge dofs (only nonzero are taken into account)
394 * to global element dofs
395 */
396
397 answer.resize(12);
398 answer.zero();
399 if ( iEdge == 1 ) { // edge between nodes 1,2
400 answer.at(3) = 3;
401 answer.at(4) = 4;
402 answer.at(5) = 5;
403 answer.at(9) = 9;
404 answer.at(10) = 10;
405 answer.at(11) = 11;
406 } else if ( iEdge == 2 ) { // edge between nodes 2 3
407 answer.at(3) = 9;
408 answer.at(4) = 10;
409 answer.at(5) = 11;
410 answer.at(9) = 15;
411 answer.at(10) = 16;
412 answer.at(11) = 17;
413 } else if ( iEdge == 3 ) { // edge between nodes 3 1
414 answer.at(3) = 15;
415 answer.at(4) = 16;
416 answer.at(5) = 17;
417 answer.at(9) = 4;
418 answer.at(10) = 5;
419 answer.at(11) = 6;
420 } else {
421 OOFEM_ERROR("wrong edge number");
422 }
423}
424
425double
427{
428 double detJ = this->interp_lin.edgeGiveTransformationJacobian(iEdge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
429 return detJ * gp->giveWeight();
430}
431
432
433int
435{
436 // returns transformation matrix from
437 // edge local coordinate system
438 // to element local c.s
439 // (same as global c.s in this case)
440 //
441 // i.e. f(element local) = T * f(edge local)
442 //
443
444 const auto &edgeNodes = this->interp_lin.computeLocalEdgeMapping(iEdge);
445
446 auto nodeA = this->giveNode(edgeNodes.at(1) );
447 auto nodeB = this->giveNode(edgeNodes.at(2) );
448
449 FloatArray cb(3), ca(3);
450 this->giveLocalCoordinates( ca, nodeA->giveCoordinates() );
451 this->giveLocalCoordinates( cb, nodeB->giveCoordinates() );
452
453 double dx = cb.at(1) - ca.at(1);
454 double dy = cb.at(2) - ca.at(2);
455 double length = sqrt(dx * dx + dy * dy);
456
457 answer.resize(3, 3);
458 answer.zero();
459 answer.at(1, 1) = 1.0;
460 answer.at(2, 2) = dx / length;
461 answer.at(2, 3) = -dy / length;
462 answer.at(3, 2) = dy / length;
463 answer.at(3, 3) = dx / length;
464
465 return 1;
466}
467
468bool
470//converts global coordinates to local planar area coordinates,
471//does not return a coordinate in the thickness direction, but
472//does check that the point is in the element thickness
473{
474 // rotate the input point Coordinate System into the element CS
475 FloatArray inputCoords_ElCS;
476 std::vector< FloatArray >lc(3);
477 FloatArray llc;
478 this->giveLocalCoordinates(inputCoords_ElCS, coords);
479 for ( int _i = 0; _i < 3; _i++ ) {
480 this->giveLocalCoordinates(lc [ _i ], this->giveNode(_i + 1)->giveCoordinates() );
481 }
482 FEI2dTrLin _interp(1, 2);
483 bool inplane = _interp.global2local( llc, inputCoords_ElCS, FEIVertexListGeometryWrapper(lc, this->giveGeometryType()) ) > 0;
484 answer.resize(2);
485 answer.at(1) = inputCoords_ElCS.at(1);
486 answer.at(2) = inputCoords_ElCS.at(2);
487 GaussPoint _gp(NULL, 1, answer, 2.0, _2dPlate);
488 // now check if the third local coordinate is within the thickness of element
489 bool outofplane = ( fabs(inputCoords_ElCS.at(3) ) <= this->giveCrossSection()->give(CS_Thickness, & _gp) / 2. );
490
491 return inplane && outofplane;
492}
493
494
495int
497{
498 double l1 = lcoords.at(1);
499 double l2 = lcoords.at(2);
500 double l3 = 1. - l2 - l1;
501
502 answer.resize(3);
503 for ( int _i = 1; _i <= 3; _i++ ) {
504 answer.at(_i) = l1 * this->giveNode(1)->giveCoordinate(_i) + l2 * this->giveNode(2)->giveCoordinate(_i) + l3 * this->giveNode(3)->giveCoordinate(_i);
505 }
506 return true;
507}
508
509void
510DKTPlate3d::computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
511// Computes numerically the load vector of the receiver due to the body loads, at tStep.
512// load is assumed to be in global cs.
513// load vector is then transformed to coordinate system in each node.
514// (should be global coordinate system, but there may be defined
515// different coordinate system in each node)
516{
517 double dens, dV, load;
518 FloatArray force;
519 FloatMatrix T;
520
521 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
522 OOFEM_ERROR("unknown load type");
523 }
524
525 GaussIntegrationRule irule(1, this, 1, 5);
526 irule.SetUpPointsOnTriangle(1, _2dPlate);
527
528 // note: force is assumed to be in global coordinate system.
529 forLoad->computeComponentArrayAt(force, tStep, mode);
530
531 if ( force.giveSize() ) {
532 GaussPoint *gp = irule.getIntegrationPoint(0);
533
534 dens = this->giveStructuralCrossSection()->give('d', gp); // constant density assumed
535 dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp); // constant thickness assumed
536
537 answer.resize(18);
538 answer.zero();
539
540 load = force.at(1) * dens * dV / 3.0;
541 answer.at(1) = load;
542 answer.at(7) = load;
543 answer.at(13) = load;
544
545 load = force.at(2) * dens * dV / 3.0;
546 answer.at(2) = load;
547 answer.at(8) = load;
548 answer.at(14) = load;
549
550 load = force.at(3) * dens * dV / 3.0;
551 answer.at(3) = load;
552 answer.at(9) = load;
553 answer.at(15) = load;
554
555 // transform result from global cs to local element cs.
556 if ( this->computeGtoLRotationMatrix(T) ) {
557 answer.rotatedWith(T, 'n');
558 }
559 } else {
560 answer.clear(); // nil resultant
561 }
562}
563} // end namespace oofem
double length(const Vector &a)
Definition CSG.h:88
#define REGISTER_Element(class)
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) override
Definition dkt3d.C:434
bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords) override
Definition dkt3d.C:469
FloatMatrix GtoLRotationMatrix
Definition dkt3d.h:77
void giveSurfaceDofMapping(IntArray &answer, int iSurf) const override
Definition dkt3d.C:310
int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords) override
Definition dkt3d.C:496
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
Definition dkt3d.C:167
int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp) override
Definition dkt3d.C:340
void giveLocalCoordinates(FloatArray &answer, const FloatArray &global)
Definition dkt3d.C:55
void printOutputAt(FILE *file, TimeStep *tStep) override
Definition dkt3d.C:346
const FloatMatrix * computeGtoLRotationMatrix()
Definition dkt3d.C:103
void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode) override
Definition dkt3d.C:510
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition dkt3d.C:96
double computeEdgeVolumeAround(GaussPoint *gp, int iEdge) override
Definition dkt3d.C:426
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition dkt3d.C:222
void giveEdgeDofMapping(IntArray &answer, int iEdge) const override
Definition dkt3d.C:390
double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf) override
Definition dkt3d.C:333
void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Definition dkt3d.C:288
int computeLoadGToLRotationMtrx(FloatMatrix &answer) override
Definition dkt3d.C:269
DKTPlate3d(int n, Domain *d)
Definition dkt3d.C:50
void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3) override
Definition dkt3d.C:71
double computeVolumeAround(GaussPoint *gp) override
Definition dkt.C:387
DKTPlate(int n, Domain *d)
Definition dkt.C:60
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition dkt.C:293
static FEI2dTrLin interp_lin
Element geometry approximation.
Definition dkt.h:74
Element_Geometry_Type giveGeometryType() const override
Definition dkt.h:89
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
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.
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