OOFEM 3.0
Loading...
Searching...
No Matches
stressvector.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
35#include "stressvector.h"
36#include "strainvector.h"
37#include "mathfem.h"
38#include "floatarray.h"
39#include "floatmatrix.h"
40#include "error.h"
41#include "mathfem.h"
42
43namespace oofem {
44StressVector :: StressVector(MaterialMode m) : StressStrainBaseVector(m)
45{ }
46
47StressVector :: StressVector(const FloatArray &src, MaterialMode m) : StressStrainBaseVector(src, m)
48{ }
49
50void
51StressVector :: computeDeviatoricVolumetricSplit(StressVector &dev, double &vol) const
52{
53 MaterialMode myMode = this->giveStressStrainMode();
54
55 if ( myMode == _1dMat ) {
56 // 1D model
57 OOFEM_ERROR("No Split for 1D!");
58 // dev.resize(1); dev.at(1) = 0.0;
59 // vol = this->at (1);
60 } else if ( myMode == _PlaneStress ) {
61 // plane stress problem
62 OOFEM_ERROR("No Split for plane stress!");
63 // dev = *this;
64 // vol = (this->at(1)+this->at(2))/2.0;
65 // dev.at(1) -= vol;
66 // dev.at(2) -= vol;
67 } else {
68 // 3d, plane strain or axisymmetric problem
69 dev = * this;
70 vol = ( this->at(1) + this->at(2) + this->at(3) ) / 3.0;
71 dev.at(1) -= vol;
72 dev.at(2) -= vol;
73 dev.at(3) -= vol;
74 }
75}
76
77void
78StressVector :: computeDeviatoricVolumetricSum(StressVector &answer, double vol) const
79{
80 MaterialMode myMode = this->giveStressStrainMode();
81
82 if ( myMode == _1dMat ) {
83 // 1D model
84 OOFEM_ERROR("No sum for 1D!");
85 // dev.resize(1); dev.at(1) = 0.0;
86 // vol = this->at (1);
87 } else if ( myMode == _PlaneStress ) {
88 // plane stress problem
89 OOFEM_ERROR("No sum for plane stress!");
90 // dev = *this;
91 // vol = (this->at(1)+this->at(2))/2.0;
92 // dev.at(1) -= vol;
93 // dev.at(2) -= vol;
94 } else {
95 // 3d, plane strain or axisymmetric problem
96 answer = * this;
97 for ( int i = 0; i < 3; i++ ) {
98 answer(i) += vol;
99 }
100 }
101}
102
103void
104StressVector :: computePrincipalValues(FloatArray &answer) const
105{
106 //
107 // This function computes sorted principal values of a stress vector.
108 // Engineering notation is used.
109 //
110 MaterialMode myMode = this->giveStressStrainMode();
111 int size = this->giveSize();
112
113 if ( myMode == _1dMat ) {
114 answer = * this;
115 } else if ( myMode == _PlaneStress ) {
116 // 2D problem
117 double ast, dst, D;
118 answer.resize(2);
119
120 ast = this->at(1) + this->at(2);
121 dst = this->at(1) - this->at(2);
122 D = std::sqrt( dst * dst + 4.0 * this->at(3) * this->at(3) );
123 answer.at(1) = 0.5 * ( ast + D );
124 answer.at(2) = 0.5 * ( ast - D );
125 } else {
126 // 3D problem
127 double I1 = 0.0, I2 = 0.0, I3 = 0.0, s1, s2, s3;
128 FloatArray s;
129 int nonzeroFlag = 0;
130
131 this->convertToFullForm(s);
132 for ( int i = 1; i <= size; i++ ) { // is this really needed?
133 if ( fabs( s.at(i) ) > 1.e-20 ) { // is there an underflow problem for small values?
134 nonzeroFlag = 1;
135 }
136 }
137
138 answer.resize(3);
139 answer.zero();
140 if ( nonzeroFlag == 0 ) {
141 return;
142 }
143
144 I1 = s.at(1) + s.at(2) + s.at(3);
145 I2 = s.at(1) * s.at(2) + s.at(2) * s.at(3) + s.at(3) * s.at(1) -
146 ( s.at(4) * s.at(4) + s.at(5) * s.at(5) + s.at(6) * s.at(6) );
147 I3 = s.at(1) * s.at(2) * s.at(3) + 2. * s.at(4) * s.at(5) * s.at(6) -
148 ( s.at(1) * s.at(4) * s.at(4) + s.at(2) * s.at(5) * s.at(5) +
149 s.at(3) * s.at(6) * s.at(6) );
150
151 /*
152 * Call cubic3r to ensure that all three real eigenvalues will be found, because we have a symmetric tensor.
153 * This allows to overcome various rounding errors when solving general cubic equation.
154 */
155 int num;
156 cubic3r( ( double ) -1., I1, -I2, I3, & s1, & s2, & s3, &num);
157
158 if (num > 0 ) {
159 answer.at(1) = s1;
160 }
161
162 if (num > 1 ) {
163 answer.at(2) = s2;
164 }
165
166 if (num > 2 ) {
167 answer.at(3) = s3;
168 }
169
170 // sort results
171 for ( int i = 1; i < 3; i++ ) {
172 for ( int j = 1; j < 3; j++ ) {
173 if ( answer.at(j + 1) > answer.at(j) ) {
174 std :: swap( answer.at(j + 1), answer.at(j) );
175 }
176 }
177 }
178 }
179}
180
181void
182StressVector :: computePrincipalValDir(FloatArray &answer, FloatMatrix &dir) const
183{
184 //
185 // This function cumputes Principal values & directions corresponding to receiver.
186 //
187 // Return Values:
188 //
189 // matrix dir -> principal directions of strains or stresses
190 // array sp -> principal strains or stresses
191 //
192
193 FloatMatrix ss;
194 FloatArray sp;
195 int nval;
196 int nonzeroFlag = 0;
197 int size = this->giveSize();
198 MaterialMode myMode = this->giveStressStrainMode();
199
200 if ( myMode == _1dMat ) {
201 answer = * this;
202 dir.resize(1, 1);
203 dir.at(1, 1) = 1.0;
204 return;
205 } else if ( myMode == _PlaneStress ) {
206 // 2D problem
207 ss.resize(2, 2);
208 answer.resize(2);
209
210 for ( int i = 1; i <= size; i++ ) {
211 if ( fabs( this->at(i) ) > 1.e-20 ) {
212 nonzeroFlag = 1;
213 }
214 }
215
216 if ( nonzeroFlag == 0 ) {
217 answer.zero();
218 ss.zero();
219 return;
220 }
221
222 ss.at(1, 1) = this->at(1);
223 ss.at(2, 2) = this->at(2);
224 ss.at(1, 2) = ss.at(2, 1) = this->at(3);
225 } else {
226 // 3D problem
227 ss.resize(3, 3);
228 FloatArray s;
229
230 answer.resize(3);
231
232 this->convertToFullForm(s);
233 for ( int i = 1; i <= size; i++ ) {
234 if ( fabs( s.at(i) ) > 1.e-20 ) {
235 nonzeroFlag = 1;
236 }
237 }
238
239 if ( nonzeroFlag == 0 ) {
240 answer.zero();
241 ss.zero();
242 return;
243 }
244
245 ss.at(1, 1) = s.at(1);
246 ss.at(2, 2) = s.at(2);
247 ss.at(3, 3) = s.at(3);
248 ss.at(1, 2) = ss.at(2, 1) = s.at(6);
249 ss.at(1, 3) = ss.at(3, 1) = s.at(5);
250 ss.at(2, 3) = ss.at(3, 2) = s.at(4);
251 }
252
253#if 0
254 ss.Jacobi(& answer, & dir, & i);
255#else
256 ss.jaco_(answer, dir, 10);
257#endif
258 // sort results
259 nval = 3;
260 if ( myMode == _PlaneStress ) {
261 nval = 2;
262 }
263
264 for ( int ii = 1; ii < nval; ii++ ) {
265 for ( int jj = 1; jj < nval; jj++ ) {
266 if ( answer.at(jj + 1) > answer.at(jj) ) {
267 // swap eigen values and eigen vectors
268 std :: swap( answer.at(jj + 1), answer.at(jj) );
269 for ( int kk = 1; kk <= nval; kk++ ) {
270 std :: swap( dir.at(kk, jj + 1), dir.at(kk, jj) );
271 }
272 }
273 }
274 }
275}
276
277void
278StressVector :: printYourself() const
279{
280 printf("StressVector (MaterialMode %d)\n", mode);
281 for ( double x : *this ) {
282 printf("%10.3e ", x);
283 }
284
285 printf("\n");
286}
287
288
289
290double
291StressVector :: computeFirstInvariant() const
292{
293 //
294 // This function computes the first invariant
295 // of the stress vector
296 //
297 #ifdef _USE_EIGEN
298 const StressVector& values=*this;
299 #endif
300 MaterialMode myMode = giveStressStrainMode();
301 if ( myMode == _1dMat ) {
302 // 1d problem
303 return values [ 0 ];
304 } else if ( myMode == _PlaneStress ) {
305 // 2d problem: plane stress
306 return values [ 0 ] + values [ 1 ];
307 } else {
308 // plane strain, axisymmetry or full 3d
309 return values [ 0 ] + values [ 1 ] + values [ 2 ];
310 }
311}
312
313double
314StressVector :: computeSecondInvariant() const
315{
316 //
317 // This function computes the second invariant of
318 // of the deviatoric stress vector
319 //
320 #ifdef _USE_EIGEN
321 const StressVector& values=*this;
322 #endif
323 MaterialMode myMode = giveStressStrainMode();
324 if ( myMode == _1dMat ) {
325 // 1d problem
326 return .5 * values [ 0 ] * values [ 0 ];
327 } else if ( myMode == _PlaneStress ) {
328 // 2d problem: plane stress
329 return .5 * ( values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] ) + values [ 2 ] * values [ 2 ];
330 } else if ( myMode == _PlaneStrain ) {
331 // plane strain or axisymmetry
332 return .5 * ( values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] + values [ 2 ] * values [ 2 ] ) +
333 values [ 3 ] * values [ 3 ];
334 } else {
335 // 3d problem
336 return .5 * ( values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] + values [ 2 ] * values [ 2 ] ) +
337 values [ 3 ] * values [ 3 ] + values [ 4 ] * values [ 4 ] + values [ 5 ] * values [ 5 ];
338 }
339}
340
341double
342StressVector :: computeThirdInvariant() const
343{
344 //
345 // This function computes the third invariant
346 // of the deviatoric stress vector.
347 //
348 #ifdef _USE_EIGEN
349 const StressVector& values=*this;
350 #endif
351 MaterialMode myMode = giveStressStrainMode();
352 if ( myMode == _1dMat ) {
353 // 1d problem
354 return ( 1. / 3. ) * values [ 0 ] * values [ 0 ] * values [ 0 ];
355 } else if ( myMode == _PlaneStress ) {
356 // 2d problem: plane stress
357 return ( 1. / 3. ) * ( values [ 0 ] * values [ 0 ] * values [ 0 ] + 3. * values [ 1 ] * values [ 2 ] * values [ 2 ]
358 + 3. * values [ 0 ] * values [ 2 ] * values [ 2 ] + values [ 1 ] * values [ 1 ] * values [ 1 ] );
359 } else if ( myMode == _PlaneStrain ) {
360 // plane strain or axisymmetry
361 return ( 1. / 3. ) * ( values [ 0 ] * values [ 0 ] * values [ 0 ] + 3. * values [ 0 ] * values [ 3 ] * values [ 3 ] +
362 3. * values [ 1 ] * values [ 3 ] * values [ 3 ] + values [ 1 ] * values [ 1 ] * values [ 1 ] +
363 values [ 2 ] * values [ 2 ] * values [ 2 ] );
364 } else {
365 // 3d problem
366 return ( 1. / 3. ) * ( values [ 0 ] * values [ 0 ] * values [ 0 ] + 3. * values [ 0 ] * values [ 5 ] * values [ 5 ] +
367 3. * values [ 0 ] * values [ 4 ] * values [ 4 ] + 6. * values [ 3 ] * values [ 5 ] * values [ 4 ] +
368 3. * values [ 1 ] * values [ 5 ] * values [ 5 ] + 3 * values [ 2 ] * values [ 4 ] * values [ 4 ] +
369 values [ 1 ] * values [ 1 ] * values [ 1 ] + 3. * values [ 1 ] * values [ 3 ] * values [ 3 ] +
370 3. * values [ 2 ] * values [ 3 ] * values [ 3 ] + values [ 2 ] * values [ 2 ] * values [ 2 ] );
371 }
372}
373
374
375void
376StressVector :: computeAllThreeHWCoordinates(double &xsi,
377 double &rho,
378 double &theta) const
379{
380 xsi = this->computeFirstCoordinate();
381
382 StressVector deviatoricStress( this->giveStressStrainMode() );
383 double volumetricStress;
384 this->computeDeviatoricVolumetricSplit(deviatoricStress, volumetricStress);
385 rho = deviatoricStress.computeSecondCoordinate();
386 theta = deviatoricStress.computeThirdCoordinate();
387}
388
389
390double
391StressVector :: computeFirstCoordinate() const
392{
393 //
394 // This function computes the first Haigh-Westergaard coordinate
395 // from the stress state
396 //
397 return computeFirstInvariant() / std::sqrt(3.);
398}
399
400double
401StressVector :: computeSecondCoordinate() const
402{
403 //
404 // This function computes the second Haigh-Westergaard coordinate
405 // from the deviatoric stress state
406 //
407 return std::sqrt( 2. * computeSecondInvariant() );
408}
409
410double
411StressVector :: computeThirdCoordinate() const
412{
413 //
414 // This function computes the third Haigh-Westergaard coordinate
415 // from the deviatoric stress state
416 //
417 double c1 = 0.0;
418 if ( computeSecondInvariant() == 0. ) {
419 c1 = 0.0;
420 } else {
421 c1 = ( 3. * std::sqrt(3.) / 2. ) * computeThirdInvariant() / ( std::pow( computeSecondInvariant(), ( 3. / 2. ) ) );
422 }
423
424 if ( c1 > 1.0 ) {
425 c1 = 1.0;
426 }
427
428 if ( c1 < -1.0 ) {
429 c1 = -1.0;
430 }
431
432 return 1. / 3. * acos(c1);
433}
434
435void
436StressVector :: applyElasticCompliance(StrainVector &strain, const double EModulus, const double nu) const
437{
438 //
439 // This function multiplies the receiver by the elastic compliance matrix
440 // and stores the result in strain
441 //
442 #ifdef _USE_EIGEN
443 const StressVector& values=*this;
444 #endif
445 MaterialMode myMode = giveStressStrainMode();
446 if ( myMode == _1dMat ) {
447 strain(0) = values [ 0 ] / EModulus;
448 } else if ( myMode == _PlaneStress ) {
449 strain(0) = ( values [ 0 ] - nu * values [ 1 ] ) / EModulus;
450 strain(1) = ( -nu * values [ 0 ] + values [ 1 ] ) / EModulus;
451 strain(2) = ( ( 2. + 2. * nu ) * values [ 2 ] ) / EModulus;
452 } else if ( myMode == _PlaneStrain ) {
453 strain(0) = ( values [ 0 ] - nu * values [ 1 ] - nu * values [ 2 ] ) / EModulus;
454 strain(1) = ( -nu * values [ 0 ] + values [ 1 ] - nu * values [ 2 ] ) / EModulus;
455 strain(2) = ( -nu * values [ 0 ] - nu * values [ 1 ] + values [ 2 ] ) / EModulus;
456 strain(3) = 2. * ( 1. + nu ) * values [ 3 ] / EModulus;
457 } else {
458 strain(0) = ( values [ 0 ] - nu * values [ 1 ] - nu * values [ 2 ] ) / EModulus;
459 strain(1) = ( -nu * values [ 0 ] + values [ 1 ] - nu * values [ 2 ] ) / EModulus;
460 strain(2) = ( -nu * values [ 0 ] - nu * values [ 1 ] + values [ 2 ] ) / EModulus;
461 strain(3) = ( 2. * ( 1. + nu ) * values [ 3 ] ) / EModulus;
462 strain(4) = ( 2. * ( 1. + nu ) * values [ 4 ] ) / EModulus;
463 strain(5) = ( 2. * ( 1. + nu ) * values [ 5 ] ) / EModulus;
464 }
465}
466
467void
468StressVector :: applyDeviatoricElasticCompliance(StrainVector &strain,
469 const double EModulus,
470 const double nu) const
471{
472 //
473 // This function applies the elastic compliance to the deviatoric strain vector
474 //
475 applyDeviatoricElasticCompliance( strain, EModulus / 2. / ( 1. + nu ) );
476}
477
478void
479StressVector :: applyDeviatoricElasticCompliance(StrainVector &strain,
480 const double GModulus) const
481{
482 //
483 // This function applies the elastic compliance to the deviatoric strain vector
484 //
485 #ifdef _USE_EIGEN
486 const StressVector& values=*this;
487 #endif
488 MaterialMode myMode = giveStressStrainMode();
489 if ( myMode == _1dMat ) {
490 OOFEM_ERROR("No Split for 1D");
491 } else if ( myMode == _PlaneStress ) {
492 OOFEM_ERROR("No Split for Plane Stress");
493 } else if ( myMode == _PlaneStrain ) {
494 strain(0) = 1. / ( 2. * GModulus ) * values [ 0 ];
495 strain(1) = 1. / ( 2. * GModulus ) * values [ 1 ];
496 strain(2) = 1. / ( 2. * GModulus ) * values [ 2 ];
497 strain(3) = 1. / GModulus * values [ 3 ];
498 } else if ( myMode == _PlaneStrainGrad ) {
499 strain(0) = 1. / ( 2. * GModulus ) * values [ 0 ];
500 strain(1) = 1. / ( 2. * GModulus ) * values [ 1 ];
501 strain(2) = 1. / ( 2. * GModulus ) * values [ 2 ];
502 strain(3) = 1. / GModulus * values [ 3 ];
503 strain(4) = values [ 4 ];
504 } else {
505 strain(0) = 1. / ( 2. * GModulus ) * values [ 0 ];
506 strain(1) = 1. / ( 2. * GModulus ) * values [ 1 ];
507 strain(2) = 1. / ( 2. * GModulus ) * values [ 2 ];
508 strain(3) = 1. / GModulus * values [ 3 ];
509 strain(4) = 1. / GModulus * values [ 4 ];
510 strain(5) = 1. / GModulus * values [ 5 ];
511 }
512}
513
514double
515StressVector :: computeStressNorm() const
516{
517 //
518 // This function computes the tensorial Norm of the stress in engineering notation
519 //
520 #ifdef _USE_EIGEN
521 const StressVector& values=*this;
522 #endif
523 MaterialMode myMode = giveStressStrainMode();
524 if ( myMode == _1dMat ) {
525 // 1d problem
526 return fabs(values [ 0 ]);
527 } else if ( myMode == _PlaneStress ) {
528 // 2d problem: plane stress
529 return std::sqrt(values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] + 2. * values [ 2 ] * values [ 2 ]);
530 } else if ( myMode == _PlaneStrain ) {
531 // plane strain or axisymmetry
532 return std::sqrt(values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] + values [ 2 ] * values [ 2 ] +
533 2. * values [ 3 ] * values [ 3 ]);
534 } else if ( myMode == _PlaneStrainGrad ) {
535 // plane strain or axisymmetry
536 return std::sqrt(values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] + values [ 2 ] * values [ 2 ] +
537 2. * values [ 3 ] * values [ 3 ]);
538 } else {
539 // 3d problem
540 return std::sqrt(values [ 0 ] * values [ 0 ] + values [ 1 ] * values [ 1 ] + values [ 2 ] * values [ 2 ] +
541 2. * values [ 3 ] * values [ 3 ] + 2. * values [ 4 ] * values [ 4 ] + 2. * values [ 5 ] * values [ 5 ]);
542 }
543}
544
545
546void
547StressVector :: giveTranformationMtrx(FloatMatrix &answer,
548 const FloatMatrix &base,
549 int transpose) const
550//
551// returns transformation matrix for 3d - stress to another system of axes,
552// given by base.
553// In base (FloatMatrix[3,3]) there are on each column stored vectors of
554// coordinate system to which we do transformation.
555//
556// If transpose == 1 we transpose base matrix before transforming
557//
558{
559 FloatMatrix t;
560 answer.resize(6, 6);
561 answer.zero();
562
563 if ( transpose == 1 ) {
564 t.beTranspositionOf(base);
565 } else {
566 t = base;
567 }
568
569 answer.at(1, 1) = t.at(1, 1) * t.at(1, 1);
570 answer.at(1, 2) = t.at(2, 1) * t.at(2, 1);
571 answer.at(1, 3) = t.at(3, 1) * t.at(3, 1);
572 answer.at(1, 4) = 2.0 * t.at(2, 1) * t.at(3, 1);
573 answer.at(1, 5) = 2.0 * t.at(1, 1) * t.at(3, 1);
574 answer.at(1, 6) = 2.0 * t.at(1, 1) * t.at(2, 1);
575
576 answer.at(2, 1) = t.at(1, 2) * t.at(1, 2);
577 answer.at(2, 2) = t.at(2, 2) * t.at(2, 2);
578 answer.at(2, 3) = t.at(3, 2) * t.at(3, 2);
579 answer.at(2, 4) = 2.0 * t.at(2, 2) * t.at(3, 2);
580 answer.at(2, 5) = 2.0 * t.at(1, 2) * t.at(3, 2);
581 answer.at(2, 6) = 2.0 * t.at(1, 2) * t.at(2, 2);
582
583 answer.at(3, 1) = t.at(1, 3) * t.at(1, 3);
584 answer.at(3, 2) = t.at(2, 3) * t.at(2, 3);
585 answer.at(3, 3) = t.at(3, 3) * t.at(3, 3);
586 answer.at(3, 4) = 2.0 * t.at(2, 3) * t.at(3, 3);
587 answer.at(3, 5) = 2.0 * t.at(1, 3) * t.at(3, 3);
588 answer.at(3, 6) = 2.0 * t.at(1, 3) * t.at(2, 3);
589
590 answer.at(4, 1) = t.at(1, 2) * t.at(1, 3);
591 answer.at(4, 2) = t.at(2, 2) * t.at(2, 3);
592 answer.at(4, 3) = t.at(3, 2) * t.at(3, 3);
593 answer.at(4, 4) = ( t.at(2, 2) * t.at(3, 3) + t.at(3, 2) * t.at(2, 3) );
594 answer.at(4, 5) = ( t.at(1, 2) * t.at(3, 3) + t.at(3, 2) * t.at(1, 3) );
595 answer.at(4, 6) = ( t.at(1, 2) * t.at(2, 3) + t.at(2, 2) * t.at(1, 3) );
596
597 answer.at(5, 1) = t.at(1, 1) * t.at(1, 3);
598 answer.at(5, 2) = t.at(2, 1) * t.at(2, 3);
599 answer.at(5, 3) = t.at(3, 1) * t.at(3, 3);
600 answer.at(5, 4) = ( t.at(2, 1) * t.at(3, 3) + t.at(3, 1) * t.at(2, 3) );
601 answer.at(5, 5) = ( t.at(1, 1) * t.at(3, 3) + t.at(3, 1) * t.at(1, 3) );
602 answer.at(5, 6) = ( t.at(1, 1) * t.at(2, 3) + t.at(2, 1) * t.at(1, 3) );
603
604 answer.at(6, 1) = t.at(1, 1) * t.at(1, 2);
605 answer.at(6, 2) = t.at(2, 1) * t.at(2, 2);
606 answer.at(6, 3) = t.at(3, 1) * t.at(3, 2);
607 answer.at(6, 4) = ( t.at(2, 1) * t.at(3, 2) + t.at(3, 1) * t.at(2, 2) );
608 answer.at(6, 5) = ( t.at(1, 1) * t.at(3, 2) + t.at(3, 1) * t.at(1, 2) );
609 answer.at(6, 6) = ( t.at(1, 1) * t.at(2, 2) + t.at(2, 1) * t.at(1, 2) );
610}
611} // end namespace oofem
Index size() const
Definition floatarray.h:162
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
std::vector< double > values
Stored values.
Definition floatarray.h:102
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
bool jaco_(FloatArray &eval, FloatMatrix &v, int nf)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
StressStrainBaseVector(MaterialMode)
Constructor. Creates zero value stress/strain vector for given material mode.
StressStrainMatMode mode
Stress strain mode.
MaterialMode giveStressStrainMode() const
Returns the material mode of receiver.
void convertToFullForm(FloatArray &fullform) const
void computeDeviatoricVolumetricSplit(StressVector &dev, double &vol) const
double computeFirstInvariant() const
double computeSecondCoordinate() const
double computeSecondInvariant() const
void applyDeviatoricElasticCompliance(StrainVector &strain, const double EModulus, const double nu) const
StressVector(MaterialMode)
Constructor. Creates zero value stress/strain vector for given material mode.
double computeThirdCoordinate() const
double computeFirstCoordinate() const
double computeThirdInvariant() const
#define OOFEM_ERROR(...)
Definition error.h:79
FloatMatrixF< M, N > transpose(const FloatMatrixF< N, M > &mat)
Constructs transposed matrix.
void cubic3r(double a, double b, double c, double d, double *r1, double *r2, double *r3, int *num)
Definition mathfem.C:155

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