50StrainVector :: computeDeviatoricVolumetricSplit(
StrainVector &dev,
double &vol)
const
54 if ( myMode == _1dMat ) {
59 }
else if ( myMode == _PlaneStress ) {
69 vol = ( this->
at(1) + this->
at(2) + this->
at(3) ) / 3.0;
77StrainVector :: computeDeviatoricVolumetricSum(
StrainVector &answer,
const double vol)
const
81 if ( myMode == _1dMat ) {
86 }
else if ( myMode == _PlaneStress ) {
96 for (
int i = 0; i < 3; i++ ) {
103StrainVector :: computePrincipalValues(
FloatArray &answer)
const
113 if ( myMode == _1dMat ) {
115 }
else if ( myMode == _PlaneStress ) {
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 );
127 double I1 = 0.0, I2 = 0.0, I3 = 0.0, s1, s2, s3;
132 for (
int i = 1; i <=
size; i++ ) {
133 if ( fabs( s.
at(i) ) > 1.e-20 ) {
140 if ( nonzeroFlag == 0 ) {
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) );
156 cubic3r( (
double ) -1., I1, -I2, I3, & s1, & s2, & s3, &num);
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) );
182StrainVector :: computeMaxPrincipalDir(
FloatArray &answer)
const
193 for (
int i = 1; i <= n; i++ ) {
194 answer.at(i) = princdir.
at(i, 1);
217 if ( myMode == _1dMat ) {
223 }
else if ( myMode == _PlaneStress ) {
228 for (
int i = 1; i <=
size; i++ ) {
229 if ( fabs( this->
at(i) ) > 1.e-20 ) {
234 if ( nonzeroFlag == 0 ) {
239 ss.
at(1, 1) = this->
at(1);
240 ss.
at(2, 2) = this->
at(2);
242 ss.
at(1, 2) = ss.
at(2, 1) = 0.5 * this->
at(3);
251 for (
int i = 1; i <=
size; i++ ) {
252 if ( fabs( s.
at(i) ) > 1.e-20 ) {
257 if ( nonzeroFlag == 0 ) {
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);
271 ss.Jacobi(& answer, & dir, & i);
273 ss.
jaco_(answer, dir, 10);
277 if ( myMode == _PlaneStress ) {
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) ) {
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) );
295StrainVector :: printYourself()
const
300 printf(
"StrainVector (MaterialMode %d)\n",
mode);
301 for (
double x :
values ) {
302 printf(
"%10.3e ", x);
310StrainVector :: computeVolumeChange()
const
319 if ( myMode == _1dMat ) {
322 }
else if ( myMode == _PlaneStress ) {
332StrainVector :: computeStrainNorm()
const
341 if ( myMode == _1dMat ) {
344 }
else if ( myMode == _PlaneStress ) {
347 }
else if ( myMode == _PlaneStrain ) {
359StrainVector :: applyElasticStiffness(
StressVector &stress,
const double EModulus,
const double nu)
const
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 ) {
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 ] );
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 ] );
403 const double EModulus,
404 const double nu)
const
415 const double GModulus)
const
425 if ( myMode == _1dMat ) {
427 }
else if ( myMode == _PlaneStress ) {
429 }
else if ( myMode == _PlaneStrain ) {
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 ) {
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 ];
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 ];
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);
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);
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);
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) );
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) );
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) );
std::vector< double > values
Stored values.
Index giveSize() const
Returns the size of receiver.
void zero()
Zeroes all coefficients of receiver.
void resize(Index rows, Index cols)
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
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)