OOFEM 3.0
Loading...
Searching...
No Matches
strainvector.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 "strainvector.h"
36#include "stressvector.h"
37#include "mathfem.h"
38#include "floatarray.h"
39#include "floatmatrix.h"
40#include "error.h"
41
42namespace oofem {
43StrainVector :: StrainVector(MaterialMode m) : StressStrainBaseVector(m)
44{ }
45
46StrainVector :: StrainVector(const FloatArray &src, MaterialMode m) : StressStrainBaseVector(src, m)
47{ }
48
49void
50StrainVector :: computeDeviatoricVolumetricSplit(StrainVector &dev, double &vol) const
51{
52 MaterialMode myMode = this->giveStressStrainMode();
53
54 if ( myMode == _1dMat ) {
55 // 1D model
56 OOFEM_ERROR("No Split for 1D!");
57 // dev.resize(1); dev.at(1) = 0.0;
58 // vol = this->at (1);
59 } else if ( myMode == _PlaneStress ) {
60 // plane stress problem
61 OOFEM_ERROR("No Split for plane stress!");
62 // dev = *this;
63 // vol = (this->at(1)+this->at(2))/2.0;
64 // dev.at(1) -= vol;
65 // dev.at(2) -= vol;
66 } else {
67 // 3d problem or plane strain problem
68 dev = * this;
69 vol = ( this->at(1) + this->at(2) + this->at(3) ) / 3.0;
70 dev.at(1) -= vol;
71 dev.at(2) -= vol;
72 dev.at(3) -= vol;
73 }
74}
75
76void
77StrainVector :: computeDeviatoricVolumetricSum(StrainVector &answer, const double vol) const
78{
79 MaterialMode myMode = this->giveStressStrainMode();
80
81 if ( myMode == _1dMat ) {
82 // 1D model
83 OOFEM_ERROR("No sum for 1D!");
84 // dev.resize(1); dev.at(1) = 0.0;
85 // vol = this->at (1);
86 } else if ( myMode == _PlaneStress ) {
87 // plane stress problem
88 OOFEM_ERROR("No sum for plane stress!");
89 // dev = *this;
90 // vol = (this->at(1)+this->at(2))/2.0;
91 // dev.at(1) -= vol;
92 // dev.at(2) -= vol;
93 } else {
94 // 3d, plane strain or axisymmetric problem
95 answer = * this;
96 for ( int i = 0; i < 3; i++ ) {
97 answer(i) += vol;
98 }
99 }
100}
101
102void
103StrainVector :: computePrincipalValues(FloatArray &answer) const
104{
105 //
106 // This function computes sorted principal values of a strain vector.
107 // Engineering notation is used.
108 //
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 + 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++ ) {
133 if ( fabs( s.at(i) ) > 1.e-20 ) {
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 0.25 * ( 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) +
148 0.25 * ( s.at(4) * s.at(5) * s.at(6) - s.at(1) * s.at(4) * s.at(4) -
149 s.at(2) * s.at(5) * s.at(5) - 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 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
182StrainVector :: computeMaxPrincipalDir(FloatArray &answer) const
183//
184// This function computes the principal direction of the receiver
185// associated with the maximum principal value.
186//
187{
188 FloatArray princval;
189 FloatMatrix princdir;
190 this->computePrincipalValDir(princval, princdir);
191 int n = princval.giveSize();
192 answer.resize(n);
193 for ( int i = 1; i <= n; i++ ) {
194 answer.at(i) = princdir.at(i, 1);
195 }
196}
197
198void
199StrainVector :: computePrincipalValDir(FloatArray &answer, FloatMatrix &dir) const
200{
201 //
202 // This function computes principal values & directions of the receiver.
203 //
204 // Return Values:
205 //
206 // answer -> principal strains (ordered from largest to smallest)
207 // dir -> principal strain directions
208 //
209
210 FloatMatrix ss;
211 FloatArray sp;
212 int nval;
213 int nonzeroFlag = 0;
214 int size = this->giveSize();
215 MaterialMode myMode = this->giveStressStrainMode();
216
217 if ( myMode == _1dMat ) {
218 // 1D problem
219 answer = * this;
220 dir.resize(1, 1);
221 dir.at(1, 1) = 1.0;
222 return;
223 } else if ( myMode == _PlaneStress ) {
224 // 2D problem
225 ss.resize(2, 2);
226 answer.resize(2);
227
228 for ( int i = 1; i <= size; i++ ) {
229 if ( fabs( this->at(i) ) > 1.e-20 ) {
230 nonzeroFlag = 1;
231 }
232 }
233
234 if ( nonzeroFlag == 0 ) {
235 answer.zero();
236 return;
237 }
238
239 ss.at(1, 1) = this->at(1);
240 ss.at(2, 2) = this->at(2);
241
242 ss.at(1, 2) = ss.at(2, 1) = 0.5 * this->at(3);
243 } else {
244 // 3D problem
245 ss.resize(3, 3);
246 FloatArray s;
247
248 answer.resize(3);
249
250 this->convertToFullForm(s);
251 for ( int i = 1; i <= size; i++ ) {
252 if ( fabs( s.at(i) ) > 1.e-20 ) {
253 nonzeroFlag = 1;
254 }
255 }
256
257 if ( nonzeroFlag == 0 ) {
258 answer.zero();
259 return;
260 }
261
262 ss.at(1, 1) = s.at(1);
263 ss.at(2, 2) = s.at(2);
264 ss.at(3, 3) = s.at(3);
265 ss.at(1, 2) = ss.at(2, 1) = 0.5 * s.at(6);
266 ss.at(1, 3) = ss.at(3, 1) = 0.5 * s.at(5);
267 ss.at(2, 3) = ss.at(3, 2) = 0.5 * s.at(4);
268 }
269
270#if 0
271 ss.Jacobi(& answer, & dir, & i);
272#else
273 ss.jaco_(answer, dir, 10);
274#endif
275 // sort results (from the largest to the smallest eigenvalue)
276 nval = 3;
277 if ( myMode == _PlaneStress ) {
278 nval = 2;
279 }
280
281 for ( int ii = 1; ii < nval; ii++ ) {
282 for ( int jj = 1; jj < nval; jj++ ) {
283 if ( answer.at(jj + 1) > answer.at(jj) ) {
284 // swap eigenvalues and eigenvectors
285 std :: swap( answer.at(jj + 1), answer.at(jj) );
286 for ( int kk = 1; kk <= nval; kk++ ) {
287 std :: swap( dir.at(kk, jj + 1), dir.at(kk, jj) );
288 }
289 }
290 }
291 }
292}
293
294void
295StrainVector :: printYourself() const
296{
297 #ifdef _USE_EIGEN
298 const StrainVector& values(*this);
299 #endif
300 printf("StrainVector (MaterialMode %d)\n", mode);
301 for ( double x : values ) {
302 printf("%10.3e ", x);
303 }
304
305 printf("\n");
306}
307
308
309double
310StrainVector :: computeVolumeChange() const
311{
312 //
313 // This function computes the change of volume. This is different from the volumetric strain.
314 //
315 #ifdef _USE_EIGEN
316 const StrainVector& values(*this);
317 #endif
318 MaterialMode myMode = giveStressStrainMode();
319 if ( myMode == _1dMat ) {
320 // 1d problem
321 return values [ 0 ];
322 } else if ( myMode == _PlaneStress ) {
323 // 2d problem: plane stress
324 return values [ 0 ] + values [ 1 ];
325 } else {
326 // plane strain, axisymmetry or full 3d
327 return values [ 0 ] + values [ 1 ] + values [ 2 ];
328 }
329}
330
331double
332StrainVector :: computeStrainNorm() const
333{
334 //
335 // This function computes the tensorial Norm of the strain in engineering notation
336 //
337 #ifdef _USE_EIGEN
338 const StrainVector& values(*this);
339 #endif
340 MaterialMode myMode = giveStressStrainMode();
341 if ( myMode == _1dMat ) {
342 // 1d problem
343 return std::sqrt(values [ 0 ] * values [ 0 ]);
344 } else if ( myMode == _PlaneStress ) {
345 // 2d problem: plane stress
346 return std::sqrt( .5 * ( 2. * values [ 0 ] * values [ 0 ] + 2 * values [ 1 ] * values [ 1 ] + values [ 2 ] * values [ 2 ] ) );
347 } else if ( myMode == _PlaneStrain ) {
348 // plane strain or axisymmetry
349 return std::sqrt( .5 * ( 2. * values [ 0 ] * values [ 0 ] + 2. * values [ 1 ] * values [ 1 ] + 2. * values [ 2 ] * values [ 2 ] +
350 values [ 3 ] * values [ 3 ] ) );
351 } else {
352 // 3d problem
353 return std::sqrt( .5 * ( 2. * values [ 0 ] * values [ 0 ] + 2. * values [ 1 ] * values [ 1 ] + 2. * values [ 2 ] * values [ 2 ] +
354 values [ 3 ] * values [ 3 ] + values [ 4 ] * values [ 4 ] + values [ 5 ] * values [ 5 ] ) );
355 }
356}
357
358void
359StrainVector :: applyElasticStiffness(StressVector &stress, const double EModulus, const double nu) const
360{
361 //
362 // This function applies the elastic stiffness to the total strain vector
363 //
364
365 #ifdef _USE_EIGEN
366 const StrainVector& values(*this);
367 #endif
368 MaterialMode myMode = giveStressStrainMode();
369 if ( myMode == _1dMat ) {
370 stress(0) = EModulus * values [ 0 ];
371 } else if ( myMode == _PlaneStress ) {
372 double factor = EModulus / ( 1. - nu * nu );
373 stress(0) = factor * ( values [ 0 ] + nu * values [ 1 ] );
374 stress(1) = factor * ( nu * values [ 0 ] + values [ 1 ] );
375 stress(2) = factor * ( ( ( 1 - nu ) / 2. ) * values [ 2 ] );
376 } else if ( myMode == _PlaneStrain ) {
377 if ( nu >= 0.5 ) {
378 OOFEM_ERROR("nu must be less than 0.5");
379 }
380
381 double factor = EModulus / ( ( 1. + nu ) * ( 1. - 2. * nu ) );
382 stress(0) = factor * ( ( 1. - nu ) * values [ 0 ] + nu * values [ 1 ] );
383 stress(1) = factor * ( nu * values [ 0 ] + ( 1. - nu ) * values [ 1 ] );
384 stress(2) = nu * ( stress(0) + stress(1) );
385 stress(3) = factor * ( ( ( 1. - 2. * nu ) / 2. ) * values [ 3 ] );
386 } else {
387 if ( nu >= 0.5 ) {
388 OOFEM_ERROR("nu must be less than 0.5");
389 }
390
391 double factor = EModulus / ( ( 1. + nu ) * ( 1. - 2. * nu ) );
392 stress(0) = factor * ( ( 1. - nu ) * values [ 0 ] + nu * values [ 1 ] + nu * values [ 2 ] );
393 stress(1) = factor * ( nu * values [ 0 ] + ( 1. - nu ) * values [ 1 ] + nu * values [ 2 ] );
394 stress(2) = factor * ( nu * values [ 0 ] + nu * values [ 1 ] + ( 1. - nu ) * values [ 2 ] );
395 stress(3) = factor * ( ( ( 1. - 2. * nu ) / 2. ) * values [ 3 ] );
396 stress(4) = factor * ( ( ( 1. - 2. * nu ) / 2. ) * values [ 4 ] );
397 stress(5) = factor * ( ( ( 1. - 2. * nu ) / 2. ) * values [ 5 ] );
398 }
399}
400
401void
402StrainVector :: applyDeviatoricElasticStiffness(StressVector &stress,
403 const double EModulus,
404 const double nu) const
405{
406 //
407 // This function applies the elastic stiffness to the deviatoric strain vector
408 //
409
410 applyDeviatoricElasticStiffness( stress, EModulus / ( 2. * ( 1. + nu ) ) );
411}
412
413void
414StrainVector :: applyDeviatoricElasticStiffness(StressVector &stress,
415 const double GModulus) const
416{
417 //
418 // This function applies the elastic stiffness to the deviatoric strain vector
419 //
420
421 #ifdef _USE_EIGEN
422 const StrainVector& values(*this);
423 #endif
424 MaterialMode myMode = giveStressStrainMode();
425 if ( myMode == _1dMat ) {
426 OOFEM_ERROR("No Split for 1D");
427 } else if ( myMode == _PlaneStress ) {
428 OOFEM_ERROR("No Split for Plane Stress");
429 } else if ( myMode == _PlaneStrain ) {
430 if ( stress.giveStressStrainMode() != _PlaneStrain ) {
431 stress.letStressStrainModeBe(_PlaneStrain);
432 }
433
434 stress(0) = 2. * GModulus * values [ 0 ];
435 stress(1) = 2. * GModulus * values [ 1 ];
436 stress(2) = 2. * GModulus * values [ 2 ];
437 stress(3) = GModulus * values [ 3 ];
438 } else if ( myMode == _PlaneStrainGrad ) {
439 if ( stress.giveStressStrainMode() != _PlaneStrainGrad ) {
440 stress.letStressStrainModeBe(_PlaneStrainGrad);
441 }
442
443 stress(0) = 2. * GModulus * values [ 0 ];
444 stress(1) = 2. * GModulus * values [ 1 ];
445 stress(2) = 2. * GModulus * values [ 2 ];
446 stress(3) = GModulus * values [ 3 ];
447 stress(4) = values [ 4 ];
448 } else {
449 if ( stress.giveStressStrainMode() != _3dMat ) {
450 stress.letStressStrainModeBe(_3dMat);
451 }
452
453 stress(0) = 2. * GModulus * values [ 0 ];
454 stress(1) = 2. * GModulus * values [ 1 ];
455 stress(2) = 2. * GModulus * values [ 2 ];
456 stress(3) = GModulus * values [ 3 ];
457 stress(4) = GModulus * values [ 4 ];
458 stress(5) = GModulus * values [ 5 ];
459 }
460}
461
462
463
464void
465StrainVector :: giveTranformationMtrx(FloatMatrix &answer,
466 const FloatMatrix &base,
467 int transpose) const
468//
469// returns transformation matrix for 3d - strains to another system of axes,
470// given by base.
471// In base (FloatMatrix[3,3]) there are on each column stored vectors of
472// coordinate system to which we do transformation.
473//
474// If transpose == 1 we transpose base matrix before transforming
475//
476{
477 FloatMatrix t;
478 answer.resize(6, 6);
479 answer.zero();
480
481 if ( transpose == 1 ) {
482 t.beTranspositionOf(base);
483 } else {
484 t = base;
485 }
486
487 answer.at(1, 1) = t.at(1, 1) * t.at(1, 1);
488 answer.at(1, 2) = t.at(2, 1) * t.at(2, 1);
489 answer.at(1, 3) = t.at(3, 1) * t.at(3, 1);
490 answer.at(1, 4) = t.at(2, 1) * t.at(3, 1);
491 answer.at(1, 5) = t.at(1, 1) * t.at(3, 1);
492 answer.at(1, 6) = t.at(1, 1) * t.at(2, 1);
493
494 answer.at(2, 1) = t.at(1, 2) * t.at(1, 2);
495 answer.at(2, 2) = t.at(2, 2) * t.at(2, 2);
496 answer.at(2, 3) = t.at(3, 2) * t.at(3, 2);
497 answer.at(2, 4) = t.at(2, 2) * t.at(3, 2);
498 answer.at(2, 5) = t.at(1, 2) * t.at(3, 2);
499 answer.at(2, 6) = t.at(1, 2) * t.at(2, 2);
500
501 answer.at(3, 1) = t.at(1, 3) * t.at(1, 3);
502 answer.at(3, 2) = t.at(2, 3) * t.at(2, 3);
503 answer.at(3, 3) = t.at(3, 3) * t.at(3, 3);
504 answer.at(3, 4) = t.at(2, 3) * t.at(3, 3);
505 answer.at(3, 5) = t.at(1, 3) * t.at(3, 3);
506 answer.at(3, 6) = t.at(1, 3) * t.at(2, 3);
507
508 answer.at(4, 1) = 2.0 * t.at(1, 2) * t.at(1, 3);
509 answer.at(4, 2) = 2.0 * t.at(2, 2) * t.at(2, 3);
510 answer.at(4, 3) = 2.0 * t.at(3, 2) * t.at(3, 3);
511 answer.at(4, 4) = ( t.at(2, 2) * t.at(3, 3) + t.at(3, 2) * t.at(2, 3) );
512 answer.at(4, 5) = ( t.at(1, 2) * t.at(3, 3) + t.at(3, 2) * t.at(1, 3) );
513 answer.at(4, 6) = ( t.at(1, 2) * t.at(2, 3) + t.at(2, 2) * t.at(1, 3) );
514
515 answer.at(5, 1) = 2.0 * t.at(1, 1) * t.at(1, 3);
516 answer.at(5, 2) = 2.0 * t.at(2, 1) * t.at(2, 3);
517 answer.at(5, 3) = 2.0 * t.at(3, 1) * t.at(3, 3);
518 answer.at(5, 4) = ( t.at(2, 1) * t.at(3, 3) + t.at(3, 1) * t.at(2, 3) );
519 answer.at(5, 5) = ( t.at(1, 1) * t.at(3, 3) + t.at(3, 1) * t.at(1, 3) );
520 answer.at(5, 6) = ( t.at(1, 1) * t.at(2, 3) + t.at(2, 1) * t.at(1, 3) );
521
522 answer.at(6, 1) = 2.0 * t.at(1, 1) * t.at(1, 2);
523 answer.at(6, 2) = 2.0 * t.at(2, 1) * t.at(2, 2);
524 answer.at(6, 3) = 2.0 * t.at(3, 1) * t.at(3, 2);
525 answer.at(6, 4) = ( t.at(2, 1) * t.at(3, 2) + t.at(3, 1) * t.at(2, 2) );
526 answer.at(6, 5) = ( t.at(1, 1) * t.at(3, 2) + t.at(3, 1) * t.at(1, 2) );
527 answer.at(6, 6) = ( t.at(1, 1) * t.at(2, 2) + t.at(2, 1) * t.at(1, 2) );
528}
529} // 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)
double at(std::size_t i, std::size_t j) const
void applyDeviatoricElasticStiffness(StressVector &stress, const double EModulus, const double nu) const
void computePrincipalValDir(FloatArray &answer, FloatMatrix &dir) const override
StrainVector(MaterialMode)
Constructor. Creates zero value stress/strain vector for given material mode.
void letStressStrainModeBe(const MaterialMode newMode)
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
#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