MarmotVoigt.h
Go to the documentation of this file.
1 /* ---------------------------------------------------------------------
2  * _
3  * _ __ ___ __ _ _ __ _ __ ___ ___ | |_
4  * | '_ ` _ \ / _` | '__| '_ ` _ \ / _ \| __|
5  * | | | | | | (_| | | | | | | | | (_) | |_
6  * |_| |_| |_|\__,_|_| |_| |_| |_|\___/ \__|
7  *
8  * Unit of Strength of Materials and Structural Analysis
9  * University of Innsbruck,
10  * 2020 - today
11  *
12  * festigkeitslehre@uibk.ac.at
13  *
14  * Matthias Neuner matthias.neuner@uibk.ac.at
15  * Magdalena Schreter magdalena.schreter@uibk.ac.at
16  * Alexander Dummer alexander.dummer@uibk.ac.at
17  *
18  * This file is part of the MAteRialMOdellingToolbox (marmot).
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  * The full text of the license can be found in the file LICENSE.md at
26  * the top level directory of marmot.
27  * ---------------------------------------------------------------------
28  */
29 
30 #pragma once
31 #include "Marmot/MarmotJournal.h"
32 #include "Marmot/MarmotMath.h"
33 #include "Marmot/MarmotTypedefs.h"
34 
35 #define VOIGTFROMDIM( x ) ( ( ( x * x ) + x ) >> 1 )
36 
40 namespace Marmot {
41  namespace ContinuumMechanics::VoigtNotation {
42 
43  /* constexpr int VoigtSize = 6; */
44 
45  enum VoigtSize { OneD = 1, TwoD = 3, ThreeD = 6, Axial = 4 };
46 
47  constexpr VoigtSize voigtSizeFromDimension( int x )
48  {
49  return (VoigtSize)( ( ( x * x ) + x ) >> 1 );
50  }
51 
52  extern const Marmot::Vector6d P;
53  extern const Marmot::Vector6d PInv;
54 
55  extern const Marmot::Vector6d I;
56  extern const Marmot::Vector6d IHyd;
57  extern const Matrix6d IDev;
58 
59  // Plane Stress handling
60 
79 
98 
118  template < enum VoigtSize voigtSize >
119  Eigen::Matrix< double, voigtSize, 1 > reduce3DVoigt( const Marmot::Vector6d& Voigt3D )
120  {
121  if constexpr ( voigtSize == OneD )
122  return ( Eigen::Matrix< double, 1, 1 >() << Voigt3D( 0 ) ).finished();
123  else if constexpr ( voigtSize == TwoD )
124  return voigtToPlaneVoigt( Voigt3D );
125  else if constexpr ( voigtSize == ThreeD )
126  return Voigt3D;
127  else
128  throw std::invalid_argument( MakeString() << __PRETTY_FUNCTION__ << ": invalid dimension specified" );
129  }
130 
151  template < enum VoigtSize voigtSize >
152  Marmot::Vector6d make3DVoigt( const Eigen::Matrix< double, voigtSize, 1 >& Voigt )
153  {
154  if constexpr ( voigtSize == OneD )
155  return ( Marmot::Vector6d() << Voigt( 0 ), 0, 0, 0, 0, 0 ).finished();
156  else if constexpr ( voigtSize == TwoD )
157  return planeVoigtToVoigt( Voigt );
158  else if constexpr ( voigtSize == ThreeD )
159  return Voigt;
160  else
161  throw std::invalid_argument( MakeString() << __PRETTY_FUNCTION__ << ": invalid dimension specified" );
162  }
163 
164  // function prototypes for Marmot::Vector6d handling
165 
183  template < typename T >
184  Eigen::Matrix< T, 3, 3 > voigtToStrain( const Eigen::Matrix< T, 6, 1 >& voigt )
185  {
186  Eigen::Matrix< T, 3, 3 > strain;
187  // clang-format off
188  strain << voigt[0], voigt[3] * 0.5, voigt[4] * 0.5,
189  voigt[3] * 0.5, voigt[1], voigt[5] * 0.5,
190  voigt[4] * 0.5, voigt[5] * 0.5, voigt[2];
191  // clang-format on
192  return strain;
193  }
194 
212  template < typename T >
213  Eigen::Matrix< T, 3, 3 > voigtToStress( const Eigen::Matrix< T, 6, 1 >& voigt )
214  {
215  Eigen::Matrix< T, 3, 3 > stress;
216  // clang-format off
217  stress << voigt[0], voigt[3], voigt[4],
218  voigt[3], voigt[1], voigt[5],
219  voigt[4], voigt[5], voigt[2];
220  // clang-format on
221  return stress;
222  }
223 
242 
260  template < typename T = double >
261  Eigen::Matrix< T, 6, 1 > stressToVoigt( const Eigen::Matrix< T, 3, 3 >& stressTensor )
262  {
263  Eigen::Matrix< T, 6, 1 > stress;
264  // clang-format off
265  stress << stressTensor( 0, 0 ),
266  stressTensor( 1, 1 ),
267  stressTensor( 2, 2 ),
268  stressTensor( 0, 1 ),
269  stressTensor( 0, 2 ),
270  stressTensor( 1, 2 );
271  // clang-format on
272  return stress;
273  }
274 
275  Eigen::Matrix< double, 6, 6 > stiffnessToVoigt( const Eigen::Tensor< double, 4 >& C );
276 
277  template < int nDim >
278  Eigen::Matrix< double, nDim, nDim > stressMatrixFromVoigt(
279  const Eigen::Matrix< double, VOIGTFROMDIM( nDim ), 1 >& Voigt )
280  {
281  if constexpr ( nDim == 1 )
282  return ( Eigen::Matrix< double, nDim, nDim >() << Voigt( 0 ) ).finished();
283  else if constexpr ( nDim == 2 )
284  return ( Eigen::Matrix< double, nDim, nDim >() << Voigt( 0 ), Voigt( 2 ), Voigt( 2 ), Voigt( 1 ) ).finished();
285  else if constexpr ( nDim == 3 )
286  return voigtToStress( Voigt );
287  else
288  throw std::invalid_argument( MakeString() << __PRETTY_FUNCTION__ << ": invalid dimension specified" );
289  }
290 
291  template < int nDim >
292  Eigen::Matrix< double, VOIGTFROMDIM( nDim ), 1 > voigtFromStrainMatrix(
293  const Eigen::Matrix< double, nDim, nDim >& strain )
294  {
295  if constexpr ( nDim == 1 )
296  return ( Eigen::Matrix< double, VOIGTFROMDIM( nDim ), 1 >() << strain( 0, 0 ) ).finished();
297  else if constexpr ( nDim == 2 )
298  return ( Eigen::Matrix< double, VOIGTFROMDIM( nDim ), 1 >() << strain( 0, 0 ),
299  strain( 1, 1 ),
300  2 * strain( 0, 1 ) )
301  .finished();
302  else if constexpr ( nDim == 3 )
303  return strainToVoigt( strain );
304  else
305  throw std::invalid_argument( MakeString() << __PRETTY_FUNCTION__ << ": invalid dimension specified" );
306  }
307 
308  namespace Invariants {
309 
318 
327  // principal strains calculated from haigh westergaard strains ( sorted --> e1 > e2 > e3 )
328 
348  // principal stressDirections calculated by solving eigenvalue problem ( !NOT sorted! )
349 
358 
365  double vonMisesEquivalentStress( const Marmot::Vector6d& stress );
366 
372  double vonMisesEquivalentStrain( const Marmot::Vector6d& strain );
373 
376  template < typename T >
377  T normStrain( const Eigen::Matrix< T, 6, 1 >& strain )
378  {
379  return ContinuumMechanics::VoigtNotation::voigtToStrain( strain ).norm();
380  }
381 
384  double normStress( const Marmot::Vector6d& stress );
385  // Trace of compressive strains
386 
394  double StrainVolumetricNegative( const Marmot::Vector6d& strain );
395 
401  template < typename T >
402  T I1( const Eigen::Matrix< T, 6, 1 >& stress )
403  {
404  return stress( 0 ) + stress( 1 ) + stress( 2 );
405  }
406 
413  double I1Strain( const Marmot::Vector6d& strain ); //
414 
420  template < typename T >
421  T I2( const Eigen::Matrix< T, 6, 1 >& stress )
422  {
423  const Eigen::Matrix< T, 6, 1 >& s = stress;
424  return s( 0 ) * s( 1 ) + s( 1 ) * s( 2 ) + s( 2 ) * s( 0 ) - s( 3 ) * s( 3 ) - s( 4 ) * s( 4 ) -
425  s( 5 ) * s( 5 );
426  }
427 
434  double I2Strain( const Marmot::Vector6d& strain );
435 
441  template < typename T >
442  T I3( const Eigen::Matrix< T, 6, 1 >& stress )
443  {
444  const Eigen::Matrix< T, 6, 1 >& s = stress;
445  return s( 0 ) * s( 1 ) * s( 2 ) + 2. * s( 3 ) * s( 4 ) * s( 5 ) - s( 0 ) * s( 5 ) * s( 5 ) -
446  s( 1 ) * s( 4 ) * s( 4 ) - s( 2 ) * s( 3 ) * s( 3 );
447  }
448 
452  double I3Strain( const Marmot::Vector6d& strain );
453 
459  template < typename T >
460  T J2( const Eigen::Matrix< T, 6, 1 >& stress )
461  {
462  const T I1_ = I1( stress );
463  const T I2_ = I2( stress );
464  const T res = ( 1. / 3 ) * I1_ * I1_ - I2_;
465  return Marmot::Math::makeReal( res ) >= 0 ? res : T( 0.0 );
466  }
467 
475  double J2Strain( const Marmot::Vector6d& strain );
476 
482  template < typename T >
483  T J3( const Eigen::Matrix< T, 6, 1 >& stress )
484  {
485  const T I1_ = I1( stress );
486  const T I2_ = I2( stress );
487  const T I3_ = I3( stress );
488 
489  return ( 2. / 27 ) * pow( I1_, 3 ) - ( 1. / 3 ) * I1_ * I2_ + I3_;
490  }
491 
495  double J3Strain( const Marmot::Vector6d& strain );
496 
497  // principal values in voigt
498  std::pair< Eigen::Vector3d, Eigen::Matrix< double, 3, 6 > > principalValuesAndDerivatives(
499  const Eigen::Matrix< double, 6, 1 >& S );
500 
501  } // namespace Invariants
502 
503  namespace Derivatives {
504 
505  // derivatives of Haigh Westergaard stresses with respect to cauchy stress in eng. notation
506 
516  template < typename T >
517  Eigen::Matrix< T, 6, 1 > dRho_dStress( T rho, const Eigen::Matrix< T, 6, 1 >& stress )
518  {
519  if ( Marmot::Math::makeReal( rho ) <= 1e-16 )
520  return Eigen::Matrix< T, 6, 1 >::Zero();
521 
522  Eigen::Matrix< T, 6, 1 > s = IDev * stress;
523 
524  return T( 1. / rho ) * P.array() * s.array();
525  }
526 
531  Marmot::Vector6d dTheta_dStress( double theta, const Marmot::Vector6d& stress );
532 
537  double dTheta_dJ2( const Marmot::Vector6d& stress );
538 
543  double dTheta_dJ3( const Marmot::Vector6d& stress );
544 
550  double dThetaStrain_dJ2Strain( const Marmot::Vector6d& strain );
551 
557  double dThetaStrain_dJ3Strain( const Marmot::Vector6d& strain );
558 
564 
570 
577 
584 
591 
592  // derivatives of principalStess with respect to stress
593 
599 
600  // derivatives of plastic strains with respect to strains
601 
608 
621  Matrix6d dEp_dE( const Matrix6d& CelInv, const Matrix6d& Cep );
622 
628  RowVector6d dDeltaEpv_dE( const Matrix6d& CelInv, const Matrix6d& Cep );
629 
635 
641  RowVector6d dDeltaEpvneg_dE( const Marmot::Vector6d& dEp, const Matrix6d& CelInv, const Matrix6d& Cep );
642 
643  } // namespace Derivatives
644 
645  namespace Transformations {
646 
651  Matrix6d transformationMatrixStrainVoigt( const Matrix3d& transformedCoordinateSystem );
652 
657  Matrix6d transformationMatrixStressVoigt( const Matrix3d& transformedCoordinateSystem );
658 
667  Matrix36d projectVoigtStressToPlane( const Vector3d& normalVector );
668 
674  Matrix36d projectVoigtStrainToPlane( const Vector3d& normalVector );
675 
684 
685  } // namespace Transformations
686 
687  } // namespace ContinuumMechanics::VoigtNotation
688 } // namespace Marmot
Marmot::ContinuumMechanics::VoigtNotation::Transformations::transformationMatrixStressVoigt
Matrix6d transformationMatrixStressVoigt(const Matrix3d &transformedCoordinateSystem)
Definition: MarmotVoigt.cpp:493
Marmot::ContinuumMechanics::VoigtNotation::Invariants::principalStresses
Eigen::Vector3d principalStresses(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:103
Marmot::FiniteElement::Spatial2D::voigtSize
constexpr int voigtSize
Definition: MarmotFiniteElement.h:113
Marmot::ContinuumMechanics::VoigtNotation::PInv
const Marmot::Vector6d PInv
Definition: MarmotVoigt.cpp:17
Marmot::ContinuumMechanics::VoigtNotation::Invariants::I1Strain
double I1Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:245
Marmot::ContinuumMechanics::VoigtNotation::voigtSizeFromDimension
constexpr VoigtSize voigtSizeFromDimension(int x)
Definition: MarmotVoigt.h:47
Marmot::ContinuumMechanics::VoigtNotation::Invariants::vonMisesEquivalentStrain
double vonMisesEquivalentStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:224
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dDeltaEpvneg_dE
RowVector6d dDeltaEpvneg_dE(const Marmot::Vector6d &dEp, const Matrix6d &CelInv, const Matrix6d &Cep)
Definition: MarmotVoigt.cpp:476
Marmot::ContinuumMechanics::VoigtNotation::Invariants::normStress
double normStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:232
Marmot::ContinuumMechanics::VoigtNotation::Invariants::sortedPrincipalStrains
Eigen::Vector3d sortedPrincipalStrains(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:109
Marmot::ContinuumMechanics::VoigtNotation::Transformations::projectVoigtStressToPlane
Matrix36d projectVoigtStressToPlane(const Vector3d &normalVector)
Definition: MarmotVoigt.cpp:512
Marmot::ContinuumMechanics::VoigtNotation::stressToVoigt
Eigen::Matrix< T, 6, 1 > stressToVoigt(const Eigen::Matrix< T, 3, 3 > &stressTensor)
Definition: MarmotVoigt.h:261
MarmotJournal.h
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dSortedStrainPrincipal_dStrain
Marmot::Matrix36 dSortedStrainPrincipal_dStrain(const Marmot::Vector6d &dEp)
Definition: MarmotVoigt.cpp:431
MarmotTypedefs.h
Marmot::ContinuumMechanics::VoigtNotation::Invariants::normStrain
T normStrain(const Eigen::Matrix< T, 6, 1 > &strain)
Definition: MarmotVoigt.h:377
MarmotMath.h
Marmot::Matrix6d
Eigen::Matrix< double, 6, 6 > Matrix6d
Definition: MarmotTypedefs.h:35
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dRho_dStress
Eigen::Matrix< T, 6, 1 > dRho_dStress(T rho, const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:517
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dStressPrincipals_dStress
Marmot::Matrix36 dStressPrincipals_dStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:390
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dStressMean_dStress
Vector6d dStressMean_dStress()
Definition: MarmotVoigt.cpp:281
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dThetaStrain_dJ2Strain
double dThetaStrain_dJ2Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:335
Marmot::ContinuumMechanics::VoigtNotation::Axial
@ Axial
Definition: MarmotVoigt.h:45
Marmot::ContinuumMechanics::VoigtNotation::make3DVoigt
Marmot::Vector6d make3DVoigt(const Eigen::Matrix< double, voigtSize, 1 > &Voigt)
Definition: MarmotVoigt.h:152
Marmot::ContinuumMechanics::VoigtNotation::Transformations::projectVoigtStrainToPlane
Matrix36d projectVoigtStrainToPlane(const Vector3d &normalVector)
Definition: MarmotVoigt.cpp:525
Marmot::FiniteElement::Spatial2D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:112
Marmot::ContinuumMechanics::VoigtNotation::Invariants::I1
T I1(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:402
Marmot::ContinuumMechanics::VoigtNotation::Invariants::J2
T J2(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:460
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dJ3Strain_dStrain
Marmot::Vector6d dJ3Strain_dStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:376
Marmot::ContinuumMechanics::VoigtNotation::voigtToPlaneVoigt
Eigen::Vector3d voigtToPlaneVoigt(const Marmot::Vector6d &voigt)
Definition: MarmotVoigt.cpp:74
Marmot::ContinuumMechanics::VoigtNotation::Transformations::rotateVoigtStress
Marmot::Vector6d rotateVoigtStress(const Eigen::Matrix3d &Q, const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:533
Marmot::Vector3d
Eigen::Matrix< double, 3, 1 > Vector3d
Definition: MarmotTypedefs.h:42
Marmot::Math::makeReal
double makeReal(const double &value)
Definition: MarmotMath.cpp:37
Marmot::ContinuumMechanics::VoigtNotation::ThreeD
@ ThreeD
Definition: MarmotVoigt.h:45
Marmot::ContinuumMechanics::VoigtNotation::strainToVoigt
Marmot::Vector6d strainToVoigt(const Eigen::Matrix3d &strainTensor)
VOIGTFROMDIM
#define VOIGTFROMDIM(x)
Definition: MarmotVoigt.h:35
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dTheta_dJ3
double dTheta_dJ3(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:319
Marmot::ContinuumMechanics::VoigtNotation::Invariants::I2
T I2(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:421
Marmot::ContinuumMechanics::VoigtNotation::I
const Marmot::Vector6d I
Definition: MarmotVoigt.cpp:19
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dStrainVolumetricNegative_dStrainPrincipal
Eigen::Vector3d dStrainVolumetricNegative_dStrainPrincipal(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:410
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dEp_dE
Matrix6d dEp_dE(const Matrix6d &CelInv, const Matrix6d &Cep)
Definition: MarmotVoigt.cpp:421
Marmot::ContinuumMechanics::VoigtNotation::voigtFromStrainMatrix
Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > voigtFromStrainMatrix(const Eigen::Matrix< double, nDim, nDim > &strain)
Definition: MarmotVoigt.h:292
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dTheta_dStress
Marmot::Vector6d dTheta_dStress(double theta, const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:286
Marmot::ContinuumMechanics::VoigtNotation::voigtToStrain
Eigen::Matrix< T, 3, 3 > voigtToStrain(const Eigen::Matrix< T, 6, 1 > &voigt)
Definition: MarmotVoigt.h:184
Marmot::ContinuumMechanics::VoigtNotation::Transformations::transformationMatrixStrainVoigt
Matrix6d transformationMatrixStrainVoigt(const Matrix3d &transformedCoordinateSystem)
Definition: MarmotVoigt.cpp:484
Marmot::Matrix3d
Eigen::Matrix< double, 3, 3 > Matrix3d
Definition: MarmotTypedefs.h:40
Marmot::ContinuumMechanics::VoigtNotation::reduce3DVoigt
Eigen::Matrix< double, voigtSize, 1 > reduce3DVoigt(const Marmot::Vector6d &Voigt3D)
Definition: MarmotVoigt.h:119
Marmot::ContinuumMechanics::VoigtNotation::Invariants::J3
T J3(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:483
Marmot
This file includes functions needed for calculations with stress and strain tensors written in voigt ...
Definition: MarmotTesting.h:32
Marmot::ContinuumMechanics::VoigtNotation::Invariants::StrainVolumetricNegative
double StrainVolumetricNegative(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:237
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dTheta_dJ2
double dTheta_dJ2(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:303
Marmot::ContinuumMechanics::VoigtNotation::IHyd
const Marmot::Vector6d IHyd
Definition: MarmotVoigt.cpp:20
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dDeltaEpv_dE
RowVector6d dDeltaEpv_dE(const Matrix6d &CelInv, const Matrix6d &Cep)
Definition: MarmotVoigt.cpp:426
Marmot::ContinuumMechanics::VoigtNotation::P
const Marmot::Vector6d P
Definition: MarmotVoigt.cpp:16
Marmot::ContinuumMechanics::VoigtNotation::Invariants::J2Strain
double J2Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:265
Marmot::ContinuumMechanics::VoigtNotation::VoigtSize
VoigtSize
Definition: MarmotVoigt.h:45
Marmot::Vector6d
Eigen::Matrix< double, 6, 1 > Vector6d
Definition: MarmotTypedefs.h:43
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dJ2Strain_dStrain
Marmot::Vector6d dJ2Strain_dStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:371
Marmot::ContinuumMechanics::VoigtNotation::Invariants::vonMisesEquivalentStress
double vonMisesEquivalentStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:219
Marmot::RowVector6d
Eigen::Matrix< double, 1, 6 > RowVector6d
Definition: MarmotTypedefs.h:48
Marmot::ContinuumMechanics::VoigtNotation::voigtToStress
Eigen::Matrix< T, 3, 3 > voigtToStress(const Eigen::Matrix< T, 6, 1 > &voigt)
Definition: MarmotVoigt.h:213
Marmot::ContinuumMechanics::VoigtNotation::Invariants::principalStrains
Eigen::Vector3d principalStrains(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:97
Marmot::ContinuumMechanics::VoigtNotation::Invariants::principalStressesDirections
Eigen::Matrix3d principalStressesDirections(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:121
Marmot::ContinuumMechanics::VoigtNotation::Invariants::I3
T I3(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:442
Marmot::ContinuumMechanics::VoigtNotation::planeVoigtToVoigt
Marmot::Vector6d planeVoigtToVoigt(const Eigen::Vector3d &voigtPlane)
Marmot::ContinuumMechanics::VoigtNotation::Invariants::J3Strain
double J3Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:271
Marmot::ContinuumMechanics::VoigtNotation::Invariants::principalValuesAndDerivatives
std::pair< Eigen::Vector3d, Eigen::Matrix< double, 3, 6 > > principalValuesAndDerivatives(const Eigen::Matrix< double, 6, 1 > &S)
Definition: MarmotVoigt.cpp:129
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dJ2_dStress
Marmot::Vector6d dJ2_dStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:359
Marmot::ContinuumMechanics::VoigtNotation::OneD
@ OneD
Definition: MarmotVoigt.h:45
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dThetaStrain_dStrain
Marmot::Vector6d dThetaStrain_dStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:383
Marmot::ContinuumMechanics::VoigtNotation::stressMatrixFromVoigt
Eigen::Matrix< double, nDim, nDim > stressMatrixFromVoigt(const Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > &Voigt)
Definition: MarmotVoigt.h:278
Marmot::ContinuumMechanics::VoigtNotation::TwoD
@ TwoD
Definition: MarmotVoigt.h:45
Marmot::ContinuumMechanics::VoigtNotation::stiffnessToVoigt
Eigen::Matrix< double, 6, 6 > stiffnessToVoigt(const Eigen::Tensor< double, 4 > &C)
Definition: MarmotVoigt.cpp:46
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dJ3_dStress
Marmot::Vector6d dJ3_dStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:364
Marmot::ContinuumMechanics::VoigtNotation::Invariants::I3Strain
double I3Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:259
Marmot::ContinuumMechanics::VoigtNotation::Invariants::I2Strain
double I2Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:250
Marmot::ContinuumMechanics::VoigtNotation::IDev
const Matrix6d IDev
Definition: MarmotVoigt.cpp:22
Marmot::ContinuumMechanics::VoigtNotation::Derivatives::dThetaStrain_dJ3Strain
double dThetaStrain_dJ3Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:347
Marmot::Matrix36d
Eigen::Matrix< double, 3, 6 > Matrix36d
Definition: MarmotTypedefs.h:53
MakeString
Definition: MarmotJournal.h:32
Marmot::Matrix36
Eigen::Matrix< double, 3, 6 > Matrix36
Definition: MarmotTypedefs.h:54