OOFEM 3.0
Loading...
Searching...
No Matches
trplanrot.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 "node.h"
39#include "material.h"
40#include "gausspoint.h"
42#include "floatmatrix.h"
43#include "floatarray.h"
44#include "intarray.h"
45#include "domain.h"
46#include "verbose.h"
47#include "load.h"
48#include "mathfem.h"
49#include "classfactory.h"
50#include "paramkey.h"
51#include "parametermanager.h"
52
53namespace oofem {
56
57TrPlaneStrRot :: TrPlaneStrRot(int n, Domain *aDomain) :
58 TrPlaneStress2d(n, aDomain)
59{
63}
64
65
66void
67TrPlaneStrRot :: computeGaussPoints()
68// Sets up the array containing the four Gauss points of the receiver.
69{
70 if ( integrationRulesArray.size() == 0 ) {
71 integrationRulesArray.resize(1);
72 integrationRulesArray [ 0 ] = std::make_unique<GaussIntegrationRule>(1, this, 1, 3);
73 this->giveCrossSection()->setupIntegrationPoints(* integrationRulesArray [ 0 ], numberOfGaussPoints, this);
74 }
75}
76
77
78void
79TrPlaneStrRot :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
80// Returns part of strain-displacement matrix {B} of the receiver,
81// (epsilon_x,epsilon_y,gamma_xy) = B . r
82// type of this part is [3,9] r=(u1,w1,fi1,u2,w2,fi2,u3,w3,fi3)
83// evaluated at gp.
84{
85#if 1
86 // New version (13-09-2014 /JB)
87 // Computes the B-matrix, directly taking into account the reduced
88 // integration of the fourth strain component.
89
90 // get node coordinates
91 FloatArray x(3), y(3);
92 this->giveNodeCoordinates(x, y);
93
94 FloatArray b(3), c(3);
95
96 for ( int i = 1; i <= 3; i++ ) {
97 int j = i + 1 - i / 3 * 3;
98 int k = j + 1 - j / 3 * 3;
99 b.at(i) = y.at(j) - y.at(k);
100 c.at(i) = x.at(k) - x.at(j);
101 }
102
103 // Derivatives of the shape functions (for this special interpolation)
106
107 FloatArray center = Vec2(0.33333333, 0.33333333);
108 FloatArray nxRed = this->GiveDerivativeVX( center );
109 FloatArray nyRed = this->GiveDerivativeUY( center );
110
111 // These are the regular shape functions of a linear triangle evaluated at the center
112 FloatArray shapeFunct = Vec3( center.at(1), center.at(2), 1.0 - center.at(1) - center.at(2) ); // = {0, 0, 1}
113
114 double area = this->giveArea();
115 double detJ = 1.0 / ( 2.0 * area );
116 answer.resize(4, 9);
117 for ( int i = 1; i <= 3; i++ ) {
118 answer.at(1, 3 * i - 2) = b.at(i) * detJ;
119 answer.at(1, 3 * i - 0) = nx.at(i) * detJ;
120
121 answer.at(2, 3 * i - 1) = c.at(i) * detJ;
122 answer.at(2, 3 * i - 0) = ny.at(i) * detJ;
123
126 answer.at(3, 3 * i - 2) = c.at(i) * detJ;
127 answer.at(3, 3 * i - 1) = b.at(i) * detJ;
128 answer.at(3, 3 * i - 0) = ( nxRed.at(i) + nyRed.at(i) ) * detJ;
129
130 // Reduced integration of the fourth strain component
131 answer.at(4, 3 * i - 2) = -1. * c.at(i) * 1.0 / 4.0 / area;
132 answer.at(4, 3 * i - 1) = b.at(i) * 1.0 / 4.0 / area;
133 answer.at(4, 3 * i - 0) = ( -4. * area * shapeFunct.at(i) + nxRed.at(i) - nyRed.at(i) ) * 1.0 / 4.0 / area;
134
135 }
136
137#else
138 // OLD version - commented out 13-09-2014 // JB
139
140 // get node coordinates
141 FloatArray x(3), y(3);
142 this->giveNodeCoordinates(x, y);
143
144 //
145 FloatArray b(3), c(3);
146
147 for ( int i = 1; i <= 3; i++ ) {
148 int j = i + 1 - i / 3 * 3;
149 int k = j + 1 - j / 3 * 3;
150 b.at(i) = y.at(j) - y.at(k);
151 c.at(i) = x.at(k) - x.at(j);
152 }
153
154 //
155 double area;
156 area = this->giveArea();
157
158 //
159 int size;
160 if ( ui == ALL_STRAINS ) {
161 size = 4;
162 ui = 4;
163 } else {
164 size = ui - li + 1;
165 }
166
167 if ( ( size < 0 ) || ( size > 4 ) ) {
168 OOFEM_ERROR("ComputeBmatrixAt size mismatch");
169 }
170
171 //
172 int ind = 1;
173
174 answer.resize(size, 9);
175
176 if ( ( li <= 2 ) ) {
179
180 if ( ( li <= 1 ) && ( ui >= 1 ) ) {
181 for ( int i = 1; i <= 3; i++ ) {
182 answer.at(1, 3 * i - 2) = b.at(i) * 1. / ( 2. * area );
183 answer.at(1, 3 * i - 0) = nx.at(i) * 1. / ( 2. * area );
184 }
185
186 ind++;
187 }
188
189 if ( ( li <= 2 ) && ( ui >= 2 ) ) {
190 for ( int i = 1; i <= 3; i++ ) {
191 answer.at(2, 3 * i - 1) = c.at(i) * 1. / ( 2. * area );
192 answer.at(2, 3 * i - 0) = ny.at(i) * 1. / ( 2. * area );
193 }
194
195 ind++;
196 }
197 }
198
199 if ( ( li <= 3 ) && ( ui >= 3 ) ) {
200 GaussIntegrationRule ir(1, this, 1, 3);
201 ir.SetUpPointsOnTriangle(1, _PlaneStress);
202
205
206 for ( int i = 1; i <= 3; i++ ) {
207 answer.at(3, 3 * i - 2) = c.at(i) * 1. / ( 2. * area );
208 answer.at(3, 3 * i - 1) = b.at(i) * 1. / ( 2. * area );
209 answer.at(3, 3 * i - 0) = ( nx.at(i) + ny.at(i) ) * 1. / ( 2. * area );
210 }
211
212 ind++;
213 }
214
215 if ( ( li <= 4 ) && ( ui >= 4 ) ) {
216 if ( numberOfRotGaussPoints == 1 ) {
217 //reduced integration, evaluate at coordinate (0,0)
218 FloatArray center = {0.0, 0.0};
219 FloatArray shapeFunct(3);
220
221 FloatArray nx = this->GiveDerivativeVX( center );
222 FloatArray ny = this->GiveDerivativeUY( center );
223
224 shapeFunct.at(1) = center.at(1);
225 shapeFunct.at(2) = center.at(2);
226 shapeFunct.at(3) = 1.0 - shapeFunct.at(1) - shapeFunct.at(2);
227
228 for ( int i = 1; i <= 3; i++ ) {
229 answer.at(ind, 3 * i - 2) = -1. * c.at(i) * 1.0 / 4.0 / area;
230 answer.at(ind, 3 * i - 1) = b.at(i) * 1.0 / 4.0 / area;
231 answer.at(ind, 3 * i - 0) = ( -4. * area * shapeFunct.at(i) + nx.at(i) - ny.at(i) ) * 1.0 / 4.0 / area;
232 }
233 }
234 }
235
236
237#endif
238
239
240}
241
242
243void
244TrPlaneStrRot :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
245// Returns the displacement interpolation matrix {N} of the receiver
246// evaluated at gp.
247{
248 // get node coordinates
249 FloatArray x(3), y(3);
250 this->giveNodeCoordinates(x, y);
251
252 //
253 FloatArray b(3), c(3), d(3);
254
255 for ( int i = 1; i <= 3; i++ ) {
256 int j = i + 1 - i / 3 * 3;
257 int k = j + 1 - j / 3 * 3;
258 b.at(i) = y.at(j) - y.at(k);
259 c.at(i) = x.at(k) - x.at(j);
260 d.at(i) = sqrt( b.at(i) * b.at(i) + c.at(i) * c.at(i) );
261 }
262
263 //
264 FloatArray angles = this->GivePitch();
265
266 //
267 double l1, l2, l3, f11, f12, f13, f21, f22, f23;
268
269 l1 = iLocCoord.at(1);
270 l2 = iLocCoord.at(2);
271 l3 = 1. - l1 - l2;
272
273 f11 = d.at(2) / 2. *l1 *l3 *sin( angles.at(2) ) - d.at(3) / 2. *l2 *l1 *sin( angles.at(3) );
274 f12 = d.at(3) / 2. *l2 *l1 *sin( angles.at(3) ) - d.at(1) / 2. *l3 *l2 *sin( angles.at(1) );
275 f13 = d.at(1) / 2. *l3 *l2 *sin( angles.at(1) ) - d.at(2) / 2. *l1 *l3 *sin( angles.at(2) );
276
277 f21 = d.at(3) / 2. *l2 *l1 *cos( angles.at(3) ) - d.at(2) / 2. *l1 *l3 *cos( angles.at(2) );
278 f22 = d.at(1) / 2. *l3 *l2 *cos( angles.at(1) ) - d.at(3) / 2. *l2 *l1 *cos( angles.at(3) );
279 f23 = d.at(2) / 2. *l1 *l3 *cos( angles.at(2) ) - d.at(1) / 2. *l3 *l2 *cos( angles.at(1) );
280
281 //
282 answer.resize(3, 9);
283 answer.zero();
284
285 answer.at(1, 1) = l1;
286 answer.at(1, 3) = f11;
287 answer.at(1, 4) = l2;
288 answer.at(1, 6) = f12;
289 answer.at(1, 7) = l3;
290 answer.at(1, 9) = f13;
291
292 answer.at(2, 2) = l1;
293 answer.at(2, 3) = f21;
294 answer.at(2, 5) = l2;
295 answer.at(2, 6) = f22;
296 answer.at(2, 8) = l3;
297 answer.at(2, 9) = f23;
298
299 answer.at(3, 3) = l1;
300 answer.at(3, 6) = l2;
301 answer.at(3, 9) = l3;
302}
303
304
305double
306TrPlaneStrRot :: giveArea()
307// returns the area occupied by the receiver
308{
310 if ( area > 0 ) { // check if previously computed
311 return area;
312 }
313
314 // get node coordinates
315 FloatArray x(3), y(3);
316 this->giveNodeCoordinates(x, y);
317
318 //
319 double x1, x2, x3, y1, y2, y3;
320 x1 = x.at(1);
321 x2 = x.at(2);
322 x3 = x.at(3);
323
324 y1 = y.at(1);
325 y2 = y.at(2);
326 y3 = y.at(3);
327
328 return ( area = fabs( 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 ) ) );
329 // return ( area = 0.5 * ( x2 * y3 + x1 * y2 + y1 * x3 - x2 * y1 - x3 * y2 - x1 * y3 ) );
330}
331
332
333void
334TrPlaneStrRot :: giveNodeCoordinates(FloatArray &x, FloatArray &y)
335{
336 const auto &nc1 = this->giveNode(1)->giveCoordinates();
337 const auto &nc2 = this->giveNode(2)->giveCoordinates();
338 const auto &nc3 = this->giveNode(3)->giveCoordinates();
339
340 x.at(1) = nc1.at(1);
341 x.at(2) = nc2.at(1);
342 x.at(3) = nc3.at(1);
343
344 y.at(1) = nc1.at(2);
345 y.at(2) = nc2.at(2);
346 y.at(3) = nc3.at(2);
347
348 //if (z) {
349 // z[0] = nc1.at(3);
350 // z[1] = nc2.at(3);
351 // z[2] = nc3.at(3);
352 //}
353}
354
355
357TrPlaneStrRot :: GivePitch()
358// Returns angles between each side and global x-axis
359{
360 // get node coordinates
361 FloatArray x(3), y(3);
362 this->giveNodeCoordinates(x, y);
363
364 //
365 FloatArray angles(3);
366
367 for ( int i = 0; i < 3; i++ ) {
368 int j = i + 1 - i / 2 * 3;
369 int k = j + 1 - j / 2 * 3;
370 if ( x(k) == x(j) ) {
371 if ( y(k) > y(j) ) {
372 angles.at(i + 1) = M_PI / 2.;
373 } else {
374 angles.at(i + 1) = M_PI * 3. / 2.;
375 }
376 }
377
378 if ( x(k) > x(j) ) {
379 if ( y(k) >= y(j) ) {
380 angles.at(i + 1) = atan( ( y(k) - y(j) ) / ( x(k) - x(j) ) );
381 } else {
382 angles.at(i + 1) = 2. * M_PI - atan( ( y(j) - y(k) ) / ( x(k) - x(j) ) );
383 }
384 }
385
386 if ( x(k) < x(j) ) {
387 if ( y(k) >= y(j) ) {
388 angles.at(i + 1) = M_PI - atan( ( y(k) - y(j) ) / ( x(j) - x(k) ) );
389 } else {
390 angles.at(i + 1) = M_PI + atan( ( y(j) - y(k) ) / ( x(j) - x(k) ) );
391 }
392 }
393 }
394
395 return angles;
396}
397
398
400TrPlaneStrRot :: GiveDerivativeUX(const FloatArray &lCoords)
401{
402 // get node coordinates
403 FloatArray x(3), y(3);
404 this->giveNodeCoordinates(x, y);
405
406 //
407 FloatArray b(3), c(3), d(3);
408
409 for ( int i = 1; i <= 3; i++ ) {
410 int j = i + 1 - i / 3 * 3;
411 int k = j + 1 - j / 3 * 3;
412 b.at(i) = y.at(j) - y.at(k);
413 c.at(i) = x.at(k) - x.at(j);
414 d.at(i) = sqrt( b.at(i) * b.at(i) + c.at(i) * c.at(i) );
415 }
416
417 //
418 FloatArray angles = this->GivePitch();
419
420 //
421 FloatArray shapeFunct(3);
422 shapeFunct.at(1) = lCoords.at(1);
423 shapeFunct.at(2) = lCoords.at(2);
424 shapeFunct.at(3) = 1.0 - shapeFunct.at(1) - shapeFunct.at(2);
425
426 //
427 FloatArray nx(3);
428
429 for ( int i = 1; i <= 3; i++ ) {
430 int j = i + 1 - i / 3 * 3;
431 int k = j + 1 - j / 3 * 3;
432 nx.at(i) = ( d.at(j) / 2. * ( b.at(k) * shapeFunct.at(i) + shapeFunct.at(k) * b.at(i) ) * sin( angles.at(j) ) -
433 d.at(k) / 2. * ( b.at(i) * shapeFunct.at(j) + shapeFunct.at(i) * b.at(j) ) * sin( angles.at(k) ) );
434 }
435
436 return nx;
437}
438
439
441TrPlaneStrRot :: GiveDerivativeVX(const FloatArray &lCoords)
442{
443 // get node coordinates
444 FloatArray x(3), y(3);
445 this->giveNodeCoordinates(x, y);
446
447 //
448 FloatArray b(3), c(3), d(3);
449
450 for ( int i = 1; i <= 3; i++ ) {
451 int j = i + 1 - i / 3 * 3;
452 int k = j + 1 - j / 3 * 3;
453 b.at(i) = y.at(j) - y.at(k);
454 c.at(i) = x.at(k) - x.at(j);
455 d.at(i) = sqrt( b.at(i) * b.at(i) + c.at(i) * c.at(i) );
456 }
457
458 //
459 FloatArray angles = this->GivePitch();
460
461 //
462 FloatArray shapeFunct(3);
463 shapeFunct.at(1) = lCoords.at(1);
464 shapeFunct.at(2) = lCoords.at(2);
465 shapeFunct.at(3) = 1.0 - shapeFunct.at(1) - shapeFunct.at(2);
466
467 //
468 FloatArray nx(3);
469
470 for ( int i = 1; i <= 3; i++ ) {
471 int j = i + 1 - i / 3 * 3;
472 int k = j + 1 - j / 3 * 3;
473 nx.at(i) = ( d.at(j) / 2. * ( b.at(k) * shapeFunct.at(i) + shapeFunct.at(k) * b.at(i) ) * cos( angles.at(j) ) -
474 d.at(k) / 2. * ( b.at(i) * shapeFunct.at(j) + shapeFunct.at(i) * b.at(j) ) * cos( angles.at(k) ) ) * ( -1.0 );
475 }
476
477 return nx;
478}
479
480
482TrPlaneStrRot :: GiveDerivativeUY(const FloatArray &lCoords)
483{
484 // get node coordinates
485 FloatArray x(3), y(3);
486 this->giveNodeCoordinates(x, y);
487
488 //
489 FloatArray b(3), c(3), d(3);
490
491 for ( int i = 1; i <= 3; i++ ) {
492 int j = i + 1 - i / 3 * 3;
493 int k = j + 1 - j / 3 * 3;
494 b.at(i) = y.at(j) - y.at(k);
495 c.at(i) = x.at(k) - x.at(j);
496 d.at(i) = sqrt( b.at(i) * b.at(i) + c.at(i) * c.at(i) );
497 }
498
499 //
500 FloatArray angles = this->GivePitch();
501
502 //
503 FloatArray shapeFunct(3);
504 shapeFunct.at(1) = lCoords.at(1);
505 shapeFunct.at(2) = lCoords.at(2);
506 shapeFunct.at(3) = 1.0 - shapeFunct.at(1) - shapeFunct.at(2);
507
508 //
509 FloatArray ny(3);
510
511 for ( int i = 1; i <= 3; i++ ) {
512 int j = i + 1 - i / 3 * 3;
513 int k = j + 1 - j / 3 * 3;
514 ny.at(i) = ( d.at(j) / 2. * ( c.at(k) * shapeFunct.at(i) + shapeFunct.at(k) * c.at(i) ) * sin( angles.at(j) ) -
515 d.at(k) / 2. * ( c.at(i) * shapeFunct.at(j) + shapeFunct.at(i) * c.at(j) ) * sin( angles.at(k) ) );
516 }
517
518 return ny;
519}
520
521
523TrPlaneStrRot :: GiveDerivativeVY(const FloatArray &lCoords)
524{
525 // get node coordinates
526 FloatArray x(3), y(3);
527 this->giveNodeCoordinates(x, y);
528
529 //
530 FloatArray b(3), c(3), d(3);
531
532 for ( int i = 1; i <= 3; i++ ) {
533 int j = i + 1 - i / 3 * 3;
534 int k = j + 1 - j / 3 * 3;
535 b.at(i) = y.at(j) - y.at(k);
536 c.at(i) = x.at(k) - x.at(j);
537 d.at(i) = sqrt( b.at(i) * b.at(i) + c.at(i) * c.at(i) );
538 }
539
540 //
541 FloatArray angles = this->GivePitch();
542
543 //
544 FloatArray shapeFunct(3);
545 shapeFunct.at(1) = lCoords.at(1);
546 shapeFunct.at(2) = lCoords.at(2);
547 shapeFunct.at(3) = 1.0 - shapeFunct.at(1) - shapeFunct.at(2);
548
549 //
550 FloatArray ny(3);
551
552 for ( int i = 1; i <= 3; i++ ) {
553 int j = i + 1 - i / 3 * 3;
554 int k = j + 1 - j / 3 * 3;
555 ny.at(i) = ( d.at(j) / 2. * ( c.at(k) * shapeFunct.at(i) + shapeFunct.at(k) * c.at(i) ) * cos( angles.at(j) ) -
556 d.at(k) / 2. * ( c.at(i) * shapeFunct.at(j) + shapeFunct.at(i) * c.at(j) ) * cos( angles.at(k) ) ) * ( -1.0 );
557 }
558
559 return ny;
560}
561
562
563void
564TrPlaneStrRot :: initializeFrom(InputRecord &ir, int priority)
565{
566 ParameterManager &pm = this->giveDomain()->elementPPM;
567 TrPlaneStress2d :: initializeFrom(ir, priority);
569}
570
571void
572TrPlaneStrRot :: postInitialize()
573{
574 TrPlaneStress2d :: postInitialize();
575
576 if ( !( ( numberOfGaussPoints == 1 ) ||
577 ( numberOfGaussPoints == 4 ) ||
578 ( numberOfGaussPoints == 7 ) ) ) {
580 }
581
582 // According to the implementation of the B-matrix only one gp is supported,
583 // so it shouldn't be an optional parameter //JB
584 if ( numberOfRotGaussPoints != 1 ) {
585 OOFEM_ERROR("numberOfRotGaussPoints size mismatch - must be equal to one");
586 }
587}
588
589void
590TrPlaneStrRot :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
591{
592 auto cs = this->giveStructuralCrossSection();
593 answer = cs->giveMembraneRotStiffMtrx(rMode, gp, tStep);
594}
595
596
597void
598TrPlaneStrRot :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
599{
600 auto cs = this->giveStructuralCrossSection();
601 answer = cs->giveGeneralizedStress_MembraneRot(strain, gp, tStep);
602}
603
604
605void
606TrPlaneStrRot :: giveDofManDofIDMask(int inode, IntArray &answer) const
607{
608 answer = {D_u, D_v, R_w};
609}
610
611
612void
613TrPlaneStrRot :: computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
614// Computes numerically the load vector of the receiver due to the body loads, at tStep.
615// load is assumed to be in global cs.
616// load vector is then transformed to coordinate system in each node.
617// (should be global coordinate system, but there may be defined
618// different coordinate system in each node)
619{
620 double dens, dV, load;
621 FloatArray force;
622 FloatMatrix T;
623
624 if ( ( forLoad->giveBCGeoType() != BodyLoadBGT ) || ( forLoad->giveBCValType() != ForceLoadBVT ) ) {
625 OOFEM_ERROR("unknown load type");
626 }
627
628 // note: force is assumed to be in global coordinate system.
629 forLoad->computeComponentArrayAt(force, tStep, mode);
630
631 if ( force.giveSize() ) {
632 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
633
634 dens = this->giveStructuralCrossSection()->give('d', gp);
635 dV = this->computeVolumeAround(gp) * this->giveCrossSection()->give(CS_Thickness, gp);
636
637 answer.resize(9);
638 answer.zero();
639
640 load = force.at(1) * dens * dV / 3.0;
641 answer.at(1) = load;
642 answer.at(4) = load;
643 answer.at(7) = load;
644
645 load = force.at(2) * dens * dV / 3.0;
646 answer.at(2) = load;
647 answer.at(5) = load;
648 answer.at(8) = load;
649
650 // transform result from global cs to local element cs.
651 if ( this->computeGtoLRotationMatrix(T) ) {
652 answer.rotatedWith(T, 'n');
653 }
654 } else {
655 answer.clear(); // nil resultant
656 }
657}
658
659
660double
661TrPlaneStrRot :: giveCharacteristicLength(const FloatArray &normalToCrackPlane)
662//
663// returns receiver's characteristic length for crack band models
664// for a crack formed in the plane with normal normalToCrackPlane.
665//
666{
667 return this->giveCharacteristicLengthForPlaneElements(normalToCrackPlane);
668}
669
670
671int
672TrPlaneStrRot :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
673{
674 if ( type == IST_StressTensor ) {
675 const FloatArray &help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
676 answer = Vec6(help.at(1), help.at(2), 0., 0., 0., help.at(3));
677 return 1;
678 } else if ( type == IST_StrainTensor ) {
679 const FloatArray &help = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
680 answer = Vec6(help.at(1), help.at(2), 0., 0., 0., help.at(3));
681 return 1;
682 } else {
683 return TrPlaneStress2d :: giveIPValue(answer, gp, type, tStep);
684 }
685}
686
687} // end namespace oofem
#define REGISTER_Element(class)
Node * giveNode(int i) const
Definition element.h:629
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Definition element.C:1710
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
Domain * giveDomain() const
Definition femcmpnn.h:97
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 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
GaussPoint * getIntegrationPoint(int n)
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition load.C:84
double computeVolumeAround(GaussPoint *gp) override
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver's strain vector.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
FloatArray GiveDerivativeUX(const FloatArray &lCoords)
Definition trplanrot.C:400
FloatArray GivePitch()
Definition trplanrot.C:357
static ParamKey IPK_TrPlaneStrRot_niprot
Definition trplanrot.h:61
virtual void giveNodeCoordinates(FloatArray &x, FloatArray &y)
Definition trplanrot.C:334
FloatArray GiveDerivativeVY(const FloatArray &lCoords)
Definition trplanrot.C:523
double giveArea() override
Definition trplanrot.C:306
FloatArray GiveDerivativeUY(const FloatArray &lCoords)
Definition trplanrot.C:482
FloatArray GiveDerivativeVX(const FloatArray &lCoords)
Definition trplanrot.C:441
TrPlaneStress2d(int n, Domain *d)
Definition trplanstrss.C:58
#define OOFEM_ERROR(...)
Definition error.h:79
#define M_PI
Definition mathfem.h:52
static FloatArray Vec2(const double &a, const double &b)
Definition floatarray.h:606
@ CS_Thickness
Thickness.
@ BodyLoadBGT
Distributed body load.
Definition bcgeomtype.h:43
static FloatArray Vec6(const double &a, const double &b, const double &c, const double &d, const double &e, const double &f)
Definition floatarray.h:610
@ ForceLoadBVT
Definition bcvaltype.h:43
static FloatArray Vec3(const double &a, const double &b, const double &c)
Definition floatarray.h:607
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)
#define ALL_STRAINS

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