OOFEM 3.0
Loading...
Searching...
No Matches
dkt.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 "fei2dtrlin.h"
39#include "node.h"
40#include "material.h"
41#include "crosssection.h"
42#include "gausspoint.h"
44#include "floatmatrix.h"
45#include "floatarray.h"
46#include "intarray.h"
47#include "load.h"
48#include "mathfem.h"
49#include "classfactory.h"
50
51#ifdef __OOFEG
52 #include "oofeggraphiccontext.h"
53#endif
54
55namespace oofem {
57
59
69
70
73
74
77{
78 return & interp_lin;
79}
80
81
82void
84// Sets up the array containing the four Gauss points of the receiver.
85{
86 if ( integrationRulesArray.size() == 0 ) {
87 integrationRulesArray.resize(1);
88 integrationRulesArray [ 0 ] = std::make_unique< GaussIntegrationRule >(1, this, 1, 5);
90 }
91}
92
93
94void
95DKTPlate::computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
96// Computes numerically the load vector of the receiver due to the body loads, at tStep.
97// load is assumed to be in global cs.
98// load vector is then transformed to coordinate system in each node.
99// (should be global coordinate system, but there may be defined
100// different coordinate system in each node)
101{
102 FloatArray force;
103 FloatMatrix T;
104
105 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
106 OOFEM_ERROR("unknown load type");
107 }
108
109 GaussIntegrationRule irule(1, this, 1, 5);
110 irule.SetUpPointsOnTriangle(1, _2dPlate);
111
112 // note: force is assumed to be in global coordinate system.
113 forLoad->computeComponentArrayAt(force, tStep, mode);
114
115 if ( force.giveSize() ) {
116 GaussPoint *gp = irule.getIntegrationPoint(0);
117 double dens = this->giveStructuralCrossSection()->give('d', gp);
118 double dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp);
119
120 answer.resize(9);
121 answer.zero();
122
123 double load = force.at(1) * dens * dV / 3.0;
124 answer.at(1) = load;
125 answer.at(4) = load;
126 answer.at(7) = load;
127
128 // transform result from global cs to local element cs.
129 if ( this->computeGtoLRotationMatrix(T) ) {
130 answer.rotatedWith(T, 'n');
131 }
132 } else {
133 answer.clear(); // nil resultant
134 }
135}
136
137
138void
140// Returns the [5x9] strain-displacement matrix {B} of the receiver,
141// evaluated at gp.
142// strain components xx, yy, zz, xz, yz
143{
144 // get node coordinates
145 double x1, x2, x3, y1, y2, y3, z1, z2, z3;
146 this->giveNodeCoordinates(x1, x2, x3, y1, y2, y3, z1, z2, z3);
147
148
149 double ksi = gp->giveNaturalCoordinate(1);
150 double eta = gp->giveNaturalCoordinate(2);
151
152 double N1dk = 4.0 * ksi - 1.0;
153 double N2dk = 0.0;
154 double N3dk = -3.0 + 4.0 * ksi + 4.0 * eta;
155 double N4dk = 4.0 * eta;
156 double N5dk = -4.0 * eta;
157 double N6dk = 4.0 * ( 1.0 - 2.0 * ksi - eta );
158
159 double N1de = 0.0;
160 double N2de = 4.0 * eta - 1.0;
161 double N3de = -3.0 + 4.0 * eta + 4.0 * ksi;
162 double N4de = 4.0 * ksi;
163 double N5de = 4.0 * ( 1.0 - 2.0 * eta - ksi );
164 double N6de = -4.0 * ksi;
165
166 double A = N1dk * x1 + N2dk * x2 + N3dk * x3 + N4dk * ( x1 + x2 ) / 2 + N5dk * ( x2 + x3 ) / 2 + N6dk * ( x1 + x3 ) / 2;
167 double B = N1dk * y1 + N2dk * y2 + N3dk * y3 + N4dk * ( y1 + y2 ) / 2 + N5dk * ( y2 + y3 ) / 2 + N6dk * ( y1 + y3 ) / 2;
168 double C = N1de * x1 + N2de * x2 + N3de * x3 + N4de * ( x1 + x2 ) / 2 + N5de * ( x2 + x3 ) / 2 + N6de * ( x1 + x3 ) / 2;
169 double D = N1de * y1 + N2de * y2 + N3de * y3 + N4de * ( y1 + y2 ) / 2 + N5de * ( y2 + y3 ) / 2 + N6de * ( y1 + y3 ) / 2;
170
171 double dxdk = D / ( A * D - B * C );
172 double dydk = -C / ( A * D - B * C );
173 double dxde = -B / ( A * D - B * C );
174 double dyde = A / ( A * D - B * C );
175
176 double dN102 = N1dk * dxdk + N1de * dxde;
177 double dN104 = N2dk * dxdk + N2de * dxde;
178 double dN106 = N3dk * dxdk + N3de * dxde;
179 double dN108 = N4dk * dxdk + N4de * dxde;
180 double dN110 = N5dk * dxdk + N5de * dxde;
181 double dN112 = N6dk * dxdk + N6de * dxde;
182
183 double dN201 = -N1dk * dydk - N1de * dyde;
184 double dN203 = -N2dk * dydk - N2de * dyde;
185 double dN205 = -N3dk * dydk - N3de * dyde;
186 double dN207 = -N4dk * dydk - N4de * dyde;
187 double dN209 = -N5dk * dydk - N5de * dyde;
188 double dN211 = -N6dk * dydk - N6de * dyde;
189
190 // normals on element sides
191 double dx4 = x2 - x1;
192 double dy4 = y2 - y1;
193 double l4 = sqrt(dx4 * dx4 + dy4 * dy4);
194
195 double dx5 = x3 - x2;
196 double dy5 = y3 - y2;
197 double l5 = sqrt(dx5 * dx5 + dy5 * dy5);
198
199 double dx6 = x1 - x3;
200 double dy6 = y1 - y3;
201 double l6 = sqrt(dx6 * dx6 + dy6 * dy6);
202
203 double c4 = dy4 / l4;
204 double s4 = -dx4 / l4;
205
206 double c5 = dy5 / l5;
207 double s5 = -dx5 / l5;
208
209 double c6 = dy6 / l6;
210 double s6 = -dx6 / l6;
211
212 // transformation matrix between element DOFs (w, fi_x, fi_y) and DOFs of element with qudratic rotations
213 double T11 = -1.5 / l4 * c4;
214 double T12 = -0.25 * c4 * c4 + 0.5 * s4 * s4;
215 double T13 = -0.25 * c4 * s4 - 0.5 * c4 * s4;
216 double T14 = 1.5 / l4 * c4;
217 double T15 = -0.25 * c4 * c4 + 0.5 * s4 * s4;
218 double T16 = -0.25 * c4 * s4 - 0.5 * c4 * s4;
219
220 double T21 = -1.5 / l4 * s4;
221 double T22 = -0.25 * c4 * s4 - 0.5 * c4 * s4;
222 double T23 = -0.25 * s4 * s4 + 0.5 * c4 * c4;
223 double T24 = 1.5 / l4 * s4;
224 double T25 = -0.25 * c4 * s4 - 0.5 * c4 * s4;
225 double T26 = -0.25 * s4 * s4 + 0.5 * c4 * c4;
226
227 double T34 = -1.5 / l5 * c5;
228 double T35 = -0.25 * c5 * c5 + 0.5 * s5 * s5;
229 double T36 = -0.25 * c5 * s5 - 0.5 * c5 * s5;
230 double T37 = 1.5 / l5 * c5;
231 double T38 = -0.25 * c5 * c5 + 0.5 * s5 * s5;
232 double T39 = -0.25 * c5 * s5 - 0.5 * c5 * s5;
233
234 double T44 = -1.5 / l5 * s5;
235 double T45 = -0.25 * c5 * s5 - 0.5 * c5 * s5;
236 double T46 = -0.25 * s5 * s5 + 0.5 * c5 * c5;
237 double T47 = 1.5 / l5 * s5;
238 double T48 = -0.25 * c5 * s5 - 0.5 * c5 * s5;
239 double T49 = -0.25 * s5 * s5 + 0.5 * c5 * c5;
240
241 double T51 = 1.5 / l6 * c6;
242 double T52 = -0.25 * c6 * c6 + 0.5 * s6 * s6;
243 double T53 = -0.25 * c6 * s6 - 0.5 * c6 * s6;
244 double T57 = -1.5 / l6 * c6;
245 double T58 = -0.25 * c6 * c6 + 0.5 * s6 * s6;
246 double T59 = -0.25 * c6 * s6 - 0.5 * c6 * s6;
247
248 double T61 = 1.5 / l6 * s6;
249 double T62 = -0.25 * c6 * s6 - 0.5 * c6 * s6;
250 double T63 = -0.25 * s6 * s6 + 0.5 * c6 * c6;
251 double T67 = -1.5 / l6 * s6;
252 double T68 = -0.25 * c6 * s6 - 0.5 * c6 * s6;
253 double T69 = -0.25 * s6 * s6 + 0.5 * c6 * c6;
254
255 answer.resize(5, 9);
256 answer.zero();
257
258 answer.at(1, 1) = T21 * dN108 + T61 * dN112;
259 answer.at(1, 2) = T22 * dN108 + T62 * dN112;
260 answer.at(1, 3) = dN102 + T23 * dN108 + T63 * dN112;
261 answer.at(1, 4) = T24 * dN108 + T44 * dN110;
262 answer.at(1, 5) = T25 * dN108 + T45 * dN110;
263 answer.at(1, 6) = dN104 + T26 * dN108 + T46 * dN110;
264 answer.at(1, 7) = T47 * dN110 + T67 * dN112;
265 answer.at(1, 8) = T48 * dN110 + T68 * dN112;
266 answer.at(1, 9) = dN106 + T49 * dN110 + T69 * dN112;
267
268 answer.at(2, 1) = T11 * dN207 + T51 * dN211;
269 answer.at(2, 2) = dN201 + T12 * dN207 + T52 * dN211;
270 answer.at(2, 3) = T13 * dN207 + T53 * dN211;
271 answer.at(2, 4) = T14 * dN207 + T34 * dN209;
272 answer.at(2, 5) = dN203 + T15 * dN207 + T35 * dN209;
273 answer.at(2, 6) = T16 * dN207 + T36 * dN209;
274 answer.at(2, 7) = T37 * dN209 + T57 * dN211;
275 answer.at(2, 8) = dN205 + T38 * dN209 + T58 * dN211;
276 answer.at(2, 9) = T39 * dN209 + T59 * dN211;
277
278 answer.at(3, 1) = -T11 * dN108 - T51 * dN112 - T21 * dN207 - T61 * dN211;
279 answer.at(3, 2) = -dN102 - T12 * dN108 - T52 * dN112 - T22 * dN207 - T62 * dN211;
280 answer.at(3, 3) = -dN201 - T13 * dN108 - T53 * dN112 - T23 * dN207 - T63 * dN211;
281 answer.at(3, 4) = -T14 * dN108 - T34 * dN110 - T24 * dN207 - T44 * dN209;
282 answer.at(3, 5) = -dN104 - T15 * dN108 - T35 * dN110 - T25 * dN207 - T45 * dN209;
283 answer.at(3, 6) = -dN203 - T16 * dN108 - T36 * dN110 - T26 * dN207 - T46 * dN209;
284 answer.at(3, 7) = -T37 * dN110 - T57 * dN112 - T47 * dN209 - T67 * dN211;
285 answer.at(3, 8) = -dN106 - T38 * dN110 - T58 * dN112 - T48 * dN209 - T68 * dN211;
286 answer.at(3, 9) = -dN205 - T39 * dN110 - T59 * dN112 - T49 * dN209 - T69 * dN211;
287
288 // Note: no shear strains, no shear forces => the 4th and 5th rows are zero
289}
290
291
292void
294// Returns the [3x9] displacement interpolation matrix {N} of the receiver,
295// evaluated at gp.
296// Note: this interpolation is not available, as the deflection is cubic along the edges,
297// but not define in the interior of the element
298// Note: the interpolation of rotations is quadratic
299// NOTE: linear interpolation returned instead
300{
302
303 answer.resize(3, 9);
304 answer.zero();
305 giveInterpolation()->evalN(N, iLocCoord, FEIElementGeometryWrapper(this) );
306
307 answer.beNMatrixOf(N, 3);
308}
309
310
311void
313{
314 answer = this->giveStructuralCrossSection()->giveGeneralizedStress_Plate(strain, gp, tStep);
315}
316
317
318void
319DKTPlate::computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
320{
321 answer = this->giveStructuralCrossSection()->give2dPlateStiffMtrx(rMode, gp, tStep);
322}
323
324
325void
326DKTPlate::giveNodeCoordinates(double &x1, double &x2, double &x3,
327 double &y1, double &y2, double &y3,
328 double &z1, double &z2, double &z3)
329{
330 const auto &nc1 = this->giveNode(1)->giveCoordinates();
331 const auto &nc2 = this->giveNode(2)->giveCoordinates();
332 const auto &nc3 = this->giveNode(3)->giveCoordinates();
333
334 x1 = nc1.at(1);
335 x2 = nc2.at(1);
336 x3 = nc3.at(1);
337
338 y1 = nc1.at(2);
339 y2 = nc2.at(2);
340 y3 = nc3.at(2);
341
342 z1 = nc1.at(3);
343 z2 = nc2.at(3);
344 z3 = nc3.at(3);
345}
346
347
348void
350{
352}
353
354
355void
357{
358 answer = { D_w, R_u, R_v };
359}
360
361
362void
364// returns normal vector to midPlane in GaussPoinr gp of receiver
365{
366 FloatArray u, v;
367 u.beDifferenceOf(this->giveNode(2)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
368 v.beDifferenceOf(this->giveNode(3)->giveCoordinates(), this->giveNode(1)->giveCoordinates() );
369
370 answer.beVectorProductOf(u, v);
371 answer.normalize();
372}
373
374
375double
377//
378// returns receiver's characteristic length for crack band models
379// for a crack formed in the plane with normal normalToCrackPlane.
380//
381{
382 return this->giveCharacteristicLengthForPlaneElements(normalToCrackPlane);
383}
384
385
386double
388// Returns the portion of the receiver which is attached to gp.
389{
390 double detJ, weight;
391 std::vector< FloatArray >lc = { FloatArray(3), FloatArray(3), FloatArray(3) };
392 this->giveNodeCoordinates( lc [ 0 ].at(1), lc [ 1 ].at(1), lc [ 2 ].at(1),
393 lc [ 0 ].at(2), lc [ 1 ].at(2), lc [ 2 ].at(2),
394 lc [ 0 ].at(3), lc [ 1 ].at(3), lc [ 2 ].at(3) );
395
396 weight = gp->giveWeight();
397 detJ = fabs(this->interp_lin.giveTransformationJacobian(gp->giveNaturalCoordinates(), FEIVertexListGeometryWrapper(lc, this->giveGeometryType()) ) );
398 return detJ * weight;
399}
400
401
402void
404// Returns the lumped mass matrix of the receiver.
405{
406 answer.resize(9, 9);
407 answer.zero();
408
409 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
410
411 double dV = this->computeVolumeAround(gp);
412 double density = this->giveStructuralCrossSection()->give('d', gp);
413 double mss1 = dV * this->giveCrossSection()->give(CS_Thickness, gp) * density / 3.;
414
415 answer.at(1, 1) = mss1;
416 answer.at(4, 4) = mss1;
417 answer.at(7, 7) = mss1;
418}
419
420
421Interface *
423{
424 if ( interface == LayeredCrossSectionInterfaceType ) {
425 return static_cast< LayeredCrossSectionInterface * >( this );
426 } else if ( interface == ZZNodalRecoveryModelInterfaceType ) {
427 return static_cast< ZZNodalRecoveryModelInterface * >( this );
428 } else if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
429 return static_cast< NodalAveragingRecoveryModelInterface * >( this );
430 } else if ( interface == SPRNodalRecoveryModelInterfaceType ) {
431 return static_cast< SPRNodalRecoveryModelInterface * >( this );
432 } else if ( interface == ZZErrorEstimatorInterfaceType ) {
433 return static_cast< ZZErrorEstimatorInterface * >( this );
434 }
435
436
437 return NULL;
438}
439
440
441#define POINT_TOL 1.e-3
442
443bool
445//converts global coordinates to local planar area coordinates,
446//does not return a coordinate in the thickness direction, but
447//does check that the point is in the element thickness
448{
449 // get node coordinates
450 double x1, x2, x3, y1, y2, y3, z1, z2, z3;
451 this->giveNodeCoordinates(x1, x2, x3, y1, y2, y3, z1, z2, z3);
452
453 // Fetch local coordinates.
454 bool ok = this->interp_lin.global2local(answer, coords, FEIElementGeometryWrapper(this) ) > 0;
455
456 //check that the point is in the element and set flag
457 for ( int i = 1; i <= 3; i++ ) {
458 if ( answer.at(i) < ( 0. - POINT_TOL ) ) {
459 return false;
460 }
461
462 if ( answer.at(i) > ( 1. + POINT_TOL ) ) {
463 return false;
464 }
465 }
466
467 //get midplane location at this point
468 double midplZ;
469 midplZ = z1 * answer.at(1) + z2 * answer.at(2) + z3 * answer.at(3);
470
471 //check that the z is within the element
473 double elthick = cs->give(CS_Thickness, answer, this);
474
475 if ( elthick / 2.0 + midplZ - fabs(coords.at(3) ) < -POINT_TOL ) {
476 answer.zero();
477 return false;
478 }
479
480
481 return ok;
482}
483
484
485int
487{
488 FloatArray help;
489 answer.resize(6);
490 if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
491 if ( type == IST_ShellForceTensor ) {
492 help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
493 } else {
494 help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
495 }
496 answer.at(1) = 0.0; // nx
497 answer.at(2) = 0.0; // ny
498 answer.at(3) = 0.0; // nz
499 answer.at(4) = help.at(5); // vyz in answer
500 answer.at(5) = help.at(4); // vxz in answer
501 answer.at(6) = 0.0; // vxy
502 return 1;
503 } else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
504 if ( type == IST_ShellMomentTensor ) {
505 help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
506 } else {
507 help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
508 }
509 answer.at(1) = help.at(1); // mx
510 answer.at(2) = help.at(2); // my
511 answer.at(3) = 0.0; // mz
512 answer.at(4) = 0.0; // myz
513 answer.at(5) = 0.0; // mxz
514 answer.at(6) = help.at(3); // mxy
515 return 1;
516 } else {
517 return StructuralElement::giveIPValue(answer, gp, type, tStep);
518 }
519}
520
521
522//
523// The element interface required by NodalAveragingRecoveryModel
524//
525void
527 InternalStateType type, TimeStep *tStep)
528{
529 GaussPoint *gp;
530 if ( ( type == IST_ShellForceTensor ) || ( type == IST_ShellMomentTensor ) ||
531 ( type == IST_ShellStrainTensor ) || ( type == IST_CurvatureTensor ) ) {
532 gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
533 this->giveIPValue(answer, gp, type, tStep);
534 } else {
535 answer.clear();
536 }
537}
538
539
540//
541// The element interface required by SPRNodalRecoveryModelInterface
542//
543void
545{
546 pap.resize(3);
547 pap.at(1) = this->giveNode(1)->giveNumber();
548 pap.at(2) = this->giveNode(2)->giveNumber();
549 pap.at(3) = this->giveNode(3)->giveNumber();
550}
551
552
553void
555{
556 answer.resize(1);
557 if ( ( pap == this->giveNode(1)->giveNumber() ) ||
558 ( pap == this->giveNode(2)->giveNumber() ) ||
559 ( pap == this->giveNode(3)->giveNumber() ) ) {
560 answer.at(1) = pap;
561 } else {
562 OOFEM_ERROR("node unknown");
563 }
564}
565
566
572
573
574void
575DKTPlate::computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords)
576{
577 FloatArray n;
578
579 this->interp_lin.edgeEvalN(n, boundaryID, lcoords, FEIElementGeometryWrapper(this) );
580
581 answer.resize(3, 6);
582 answer.at(1, 1) = n.at(1);
583 answer.at(1, 4) = n.at(2);
584 answer.at(2, 2) = answer.at(3, 3) = n.at(1);
585 answer.at(2, 5) = answer.at(3, 6) = n.at(2);
586}
587
588
589
590//
591// layered cross section support functions
592//
593void
595 GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
596// returns full 3d strain vector of given layer (whose z-coordinate from center-line is
597// stored in slaveGp) for given tStep
598{
599 double layerZeta, layerZCoord, top, bottom;
600
601 top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
602 bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
603 layerZeta = slaveGp->giveNaturalCoordinate(3);
604 layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
605
606 answer.resize(5); // {Exx,Eyy,GMyz,GMzx,GMxy}
607
608 answer.at(1) = masterGpStrain.at(1) * layerZCoord;
609 answer.at(2) = masterGpStrain.at(2) * layerZCoord;
610 answer.at(5) = masterGpStrain.at(3) * layerZCoord;
611 answer.at(3) = masterGpStrain.at(5);
612 answer.at(4) = masterGpStrain.at(4);
613}
614
615void
616DKTPlate::giveEdgeDofMapping(IntArray &answer, int iEdge) const
617{
618 /*
619 * provides dof mapping of local edge dofs (only nonzero are taken into account)
620 * to global element dofs
621 */
622
623 answer.resize(6);
624 if ( iEdge == 1 ) { // edge between nodes 1,2
625 answer.at(1) = 1;
626 answer.at(2) = 2;
627 answer.at(3) = 3;
628 answer.at(4) = 4;
629 answer.at(5) = 5;
630 answer.at(6) = 6;
631 } else if ( iEdge == 2 ) { // edge between nodes 2 3
632 answer.at(1) = 4;
633 answer.at(2) = 5;
634 answer.at(3) = 6;
635 answer.at(4) = 7;
636 answer.at(5) = 8;
637 answer.at(6) = 9;
638 } else if ( iEdge == 3 ) { // edge between nodes 3 1
639 answer.at(1) = 7;
640 answer.at(2) = 8;
641 answer.at(3) = 9;
642 answer.at(4) = 1;
643 answer.at(5) = 2;
644 answer.at(6) = 3;
645 } else {
646 OOFEM_ERROR("wrong edge number");
647 }
648}
649
650double
652{
653 double detJ = this->interp_lin.edgeGiveTransformationJacobian(iEdge, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
654 return detJ * gp->giveWeight();
655}
656
657
658int
660{
661 // returns transformation matrix from
662 // edge local coordinate system
663 // to element local c.s
664 // (same as global c.s in this case)
665 //
666 // i.e. f(element local) = T * f(edge local)
667 //
668 const auto &edgeNodes = this->interp_lin.computeLocalEdgeMapping(iEdge);
669
670 auto nodeA = this->giveNode(edgeNodes.at(1) );
671 auto nodeB = this->giveNode(edgeNodes.at(2) );
672
673 double dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
674 double dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
675 double length = sqrt(dx * dx + dy * dy);
676
677 answer.resize(3, 3);
678 answer.zero();
679 answer.at(1, 1) = 1.0;
680 answer.at(2, 2) = dx / length;
681 answer.at(2, 3) = -dy / length;
682 answer.at(3, 2) = dy / length;
683 answer.at(3, 3) = dx / length;
684
685 return 1;
686}
687
688
689
690void
692{
693 this->computeNmatrixAt(sgp->giveNaturalCoordinates(), answer);
694}
695
696void
698{
699 answer.resize(9);
700 answer.zero();
701 if ( iSurf == 1 ) {
702 for ( int i = 1; i <= 9; i++ ) {
703 answer.at(i) = i;
704 }
705 } else {
706 OOFEM_ERROR("wrong surface number");
707 }
708}
709
710
711double
713{
714 return this->computeVolumeAround(gp);
715}
716
717
718int
720{
721 return 0;
722}
723
724void
726{
727#ifdef DKT_EnableVertexMomentsCache
728 if ( stateCounter == tStep->giveSolutionStateCounter() ) {
729 answer = vertexMoments;
730 return;
731 }
732#endif
733
734 // the results should be cached somehow, as computing on the fly is highly inefficient
735 // due to multiple requests
736 answer.resize(5, 3);
737
738 FloatArray eps, m;
739 FloatArray coords[ 3 ]; // vertex local coordinates
740 coords [ 0 ] = Vec2(
741 1.0, 0.0
742 );
743 coords [ 1 ] = Vec2(
744 0.0, 1.0
745 );
746 coords [ 2 ] = Vec2(
747 0.0, 0.0
748 );
749
750 GaussIntegrationRule iRule = GaussIntegrationRule(1, this, 1, 1); // dummy rule used to evaluate B at vertices
751 iRule.SetUpPointsOnTriangle(1, _Unknown);
752 GaussPoint *vgp = iRule.getIntegrationPoint(0);
753
754 for ( int i = 1; i <= this->numberOfDofMans; i++ ) {
755 vgp->setNaturalCoordinates(coords [ i - 1 ]);
756 this->computeStrainVector(eps, vgp, tStep);
757 m = this->giveStructuralCrossSection()->giveGeneralizedStress_Plate(eps, vgp, tStep);
758 answer.setColumn(m, i);
759 }
760
761#ifdef DKT_EnableVertexMomentsCache
762 this->vertexMoments = answer;
763 this->stateCounter = tStep->giveSolutionStateCounter();
764#endif
765}
766
767void
769{
770 // as shear strains are enforced to be zero (art least on element edges) the shear forces are computed from equlibrium
771 FloatMatrix m, dndx;
772 answer.resize(5);
773
774 this->computeVertexBendingMoments(m, tStep);
775 this->interp_lin.evaldNdx(dndx, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
776 for ( int i = 1; i <= this->numberOfDofMans; i++ ) {
777 answer.at(4) += m.at(1, i) * dndx.at(i, 1) + m.at(3, i) * dndx.at(i, 2); //dMxdx + dMxydy
778 answer.at(5) += m.at(2, i) * dndx.at(i, 2) + m.at(3, i) * dndx.at(i, 1); //dMydy + dMxydx
779 }
780}
781
782
783//
784// io routines
785//
786#ifdef __OOFEG
787void
789{
790 WCRec p[ 3 ];
791 GraphicObj *go;
792
793 if ( !gc.testElementGraphicActivity(this) ) {
794 return;
795 }
796
797 if ( this->giveMaterial()->isActivated(tStep) ) {
798 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
799 EASValsSetColor(gc.getElementColor() );
800 EASValsSetEdgeColor(gc.getElementEdgeColor() );
801 EASValsSetEdgeFlag(true);
802 EASValsSetFillStyle(FILL_SOLID);
803 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
804 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
805 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
806 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
807 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
808 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
809 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
810 p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
811 p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
812 p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveCoordinate(3);
813
814 go = CreateTriangle3D(p);
815 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
816 EGAttachObject(go, ( EObjectP ) this);
817 EMAddGraphicsToModel(ESIModel(), go);
818 }
819}
820
821
822void
824{
825 WCRec p[ 3 ];
826 GraphicObj *go;
827 double defScale = gc.getDefScale();
828
829 if ( !gc.testElementGraphicActivity(this) ) {
830 return;
831 }
832
833 if ( this->giveMaterial()->isActivated(tStep) ) {
834 EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
835 EASValsSetColor(gc.getDeformedElementColor() );
836 EASValsSetEdgeColor(gc.getElementEdgeColor() );
837 EASValsSetEdgeFlag(true);
838 EASValsSetFillStyle(FILL_SOLID);
839 EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
840 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
841 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
842 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
843 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
844 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
845 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
846 p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale);
847 p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale);
848 p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(3, tStep, defScale);
849
850 go = CreateTriangle3D(p);
851 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
852 EMAddGraphicsToModel(ESIModel(), go);
853 }
854}
855
856
857void
859{
860 int i, indx, result = 0;
861 WCRec p[ 3 ];
862 GraphicObj *tr;
863 FloatArray v1, v2, v3;
864 double s[ 3 ], defScale;
865
866 if ( !gc.testElementGraphicActivity(this) ) {
867 return;
868 }
869
870 if ( !this->giveMaterial()->isActivated(tStep) ) {
871 return;
872 }
873
874 if ( gc.giveIntVarMode() == ISM_recovered ) {
875 result += this->giveInternalStateAtNode(v1, gc.giveIntVarType(), gc.giveIntVarMode(), 1, tStep);
876 result += this->giveInternalStateAtNode(v2, gc.giveIntVarType(), gc.giveIntVarMode(), 2, tStep);
877 result += this->giveInternalStateAtNode(v3, gc.giveIntVarType(), gc.giveIntVarMode(), 3, tStep);
878 } else if ( gc.giveIntVarMode() == ISM_local ) {
879 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
880 result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
881 v2 = v1;
882 v3 = v1;
883 result *= 3;
884 }
885
886 if ( result != 3 ) {
887 return;
888 }
889
890 indx = gc.giveIntVarIndx();
891
892 s [ 0 ] = v1.at(indx);
893 s [ 1 ] = v2.at(indx);
894 s [ 2 ] = v3.at(indx);
895
896 EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
897
898 if ( gc.getScalarAlgo() == SA_ISO_SURF ) {
899 for ( i = 0; i < 3; i++ ) {
900 if ( gc.getInternalVarsDefGeoFlag() ) {
901 // use deformed geometry
902 defScale = gc.getDefScale();
903 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(1, tStep, defScale);
904 p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(2, tStep, defScale);
905 p [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveUpdatedCoordinate(3, tStep, defScale);
906 } else {
907 p [ i ].x = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(1);
908 p [ i ].y = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(2);
909 p [ i ].z = ( FPNum ) this->giveNode(i + 1)->giveCoordinate(3);
910 }
911 }
912
913 //EASValsSetColor(gc.getYieldPlotColor(ratio));
914 gc.updateFringeTableMinMax(s, 3);
915 tr = CreateTriangleWD3D(p, s [ 0 ], s [ 1 ], s [ 2 ]);
916 EGWithMaskChangeAttributes(LAYER_MASK, tr);
917 EMAddGraphicsToModel(ESIModel(), tr);
918 }
919}
920
921#endif
922} // end namespace oofem
double length(const Vector &a)
Definition CSG.h:88
#define N(a, b)
#define REGISTER_Element(class)
virtual double give(CrossSectionProperty a, GaussPoint *gp) const
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
SPRPatchType SPRNodalRecoveryMI_givePatchType() override
Definition dkt.C:568
bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords) override
Definition dkt.C:444
double computeVolumeAround(GaussPoint *gp) override
Definition dkt.C:387
FEInterpolation * giveInterpolation() const override
Definition dkt.C:72
void drawScalar(oofegGraphicContext &gc, TimeStep *tStep) override
Definition dkt.C:858
StateCounterType stateCounter
Definition dkt.h:78
DKTPlate(int n, Domain *d)
Definition dkt.C:60
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
Definition dkt.C:319
void giveSurfaceDofMapping(IntArray &answer, int iSurf) const override
Definition dkt.C:697
void computeShearForces(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Definition dkt.C:768
void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType type) override
Definition dkt.C:823
void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) override
Definition dkt.C:95
Interface * giveInterface(InterfaceType it) override
Definition dkt.C:422
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
Definition dkt.C:293
double computeEdgeVolumeAround(GaussPoint *gp, int iEdge) override
Definition dkt.C:651
void computeEdgeNMatrix(FloatMatrix &answer, int boundaryID, const FloatArray &lcoords) override
computes edge interpolation matrix
Definition dkt.C:575
void giveDofManDofIDMask(int inode, IntArray &) const override
Definition dkt.C:356
void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep) override
Definition dkt.C:403
void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap) override
Definition dkt.C:544
FloatMatrix vertexMoments
Definition dkt.h:77
int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp) override
Definition dkt.C:719
double giveCharacteristicLength(const FloatArray &normalToCrackPlane) override
Definition dkt.C:376
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Definition dkt.C:139
void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) override
Definition dkt.C:526
void giveEdgeDofMapping(IntArray &answer, int iEdge) const override
Definition dkt.C:616
void initializeFrom(InputRecord &ir, int priority) override
Definition dkt.C:349
double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf) override
Definition dkt.C:712
double area
Definition dkt.h:75
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition dkt.C:486
virtual void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3)
Definition dkt.C:326
void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep) override
Definition dkt.C:594
static FEI2dTrLin interp_lin
Element geometry approximation.
Definition dkt.h:74
Element_Geometry_Type giveGeometryType() const override
Definition dkt.h:89
void computeMidPlaneNormal(FloatArray &answer, const GaussPoint *gp) override
Definition dkt.C:363
void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap) override
Definition dkt.C:554
int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp) override
Definition dkt.C:659
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
Definition dkt.C:312
void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep) override
Definition dkt.C:788
void computeGaussPoints() override
Definition dkt.C:83
void computeVertexBendingMoments(FloatMatrix &answer, TimeStep *tStep)
Definition dkt.C:725
void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Definition dkt.C:691
double giveCoordinate(int i) const
Definition dofmanager.h:383
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
Node * giveNode(int i) const
Definition element.h:629
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Definition element.C:1710
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
void initializeFrom(InputRecord &ir, int priority) override
Definition element.C:687
virtual Material * giveMaterial()
Definition element.C:523
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
double giveCharacteristicLengthForPlaneElements(const FloatArray &normalToCrackPlane)
Definition element.C:1188
CrossSection * giveCrossSection()
Definition element.C:534
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
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 beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beNMatrixOf(const FloatArray &n, int nsd)
void setColumn(const FloatArray &src, int c)
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
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition gausspoint.h:136
void setNaturalCoordinates(const FloatArray &c)
Definition gausspoint.h:139
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
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Definition node.C:234
virtual FloatMatrixF< 5, 5 > give2dPlateStiffMtrx(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0
virtual FloatArrayF< 5 > giveGeneralizedStress_Plate(const FloatArrayF< 5 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
StructuralElement(int n, Domain *d)
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, 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.
StateCounterType giveSolutionStateCounter()
Definition timestep.h:211
ZZErrorEstimatorInterface(Element *element)
Constructor.
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
#define OOFEM_ERROR(...)
Definition error.h:79
#define POINT_TOL
static FloatArray Vec2(const double &a, const double &b)
Definition floatarray.h:606
@ CS_BottomZCoord
Bottom z coordinate.
@ CS_Thickness
Thickness.
@ CS_TopZCoord
Top z coordinate.
@ BodyLoadBGT
Distributed body load.
Definition bcgeomtype.h:43
@ ForceLoadBVT
Definition bcvaltype.h:43
@ SPRNodalRecoveryModelInterfaceType
@ ZZNodalRecoveryModelInterfaceType
@ LayeredCrossSectionInterfaceType
@ ZZErrorEstimatorInterfaceType
@ NodalAveragingRecoveryModelInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_VARPLOT_PATTERN_LAYER
#define OOFEG_DEFORMED_GEOMETRY_LAYER
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_LAYER

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