Go to the documentation of this file.
35 #define VOIGTFROMDIM( x ) ( ( ( x * x ) + x ) >> 1 )
41 namespace ContinuumMechanics::VoigtNotation {
49 return (
VoigtSize)( ( ( x * x ) + x ) >> 1 );
118 template < enum VoigtSize voigtSize >
122 return ( Eigen::Matrix< double, 1, 1 >() << Voigt3D( 0 ) ).finished();
128 throw std::invalid_argument(
MakeString() << __PRETTY_FUNCTION__ <<
": invalid dimension specified" );
151 template < enum VoigtSize voigtSize >
161 throw std::invalid_argument(
MakeString() << __PRETTY_FUNCTION__ <<
": invalid dimension specified" );
183 template <
typename T >
184 Eigen::Matrix< T, 3, 3 >
voigtToStrain(
const Eigen::Matrix< T, 6, 1 >& voigt )
186 Eigen::Matrix< T, 3, 3 > strain;
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];
212 template <
typename T >
213 Eigen::Matrix< T, 3, 3 >
voigtToStress(
const Eigen::Matrix< T, 6, 1 >& voigt )
215 Eigen::Matrix< T, 3, 3 > stress;
217 stress << voigt[0], voigt[3], voigt[4],
218 voigt[3], voigt[1], voigt[5],
219 voigt[4], voigt[5], voigt[2];
260 template <
typename T =
double >
261 Eigen::Matrix< T, 6, 1 >
stressToVoigt(
const Eigen::Matrix< T, 3, 3 >& stressTensor )
263 Eigen::Matrix< T, 6, 1 > stress;
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 );
275 Eigen::Matrix< double, 6, 6 >
stiffnessToVoigt(
const Eigen::Tensor< double, 4 >& C );
277 template <
int nDim >
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 )
288 throw std::invalid_argument(
MakeString() << __PRETTY_FUNCTION__ <<
": invalid dimension specified" );
291 template <
int nDim >
293 const Eigen::Matrix< double, nDim, nDim >& strain )
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 ),
302 else if constexpr (
nDim == 3 )
305 throw std::invalid_argument(
MakeString() << __PRETTY_FUNCTION__ <<
": invalid dimension specified" );
308 namespace Invariants {
376 template <
typename T >
401 template <
typename T >
402 T
I1(
const Eigen::Matrix< T, 6, 1 >& stress )
404 return stress( 0 ) + stress( 1 ) + stress( 2 );
420 template <
typename T >
421 T
I2(
const Eigen::Matrix< T, 6, 1 >& stress )
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 ) -
441 template <
typename T >
442 T
I3(
const Eigen::Matrix< T, 6, 1 >& stress )
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 );
459 template <
typename T >
460 T
J2(
const Eigen::Matrix< T, 6, 1 >& stress )
462 const T I1_ =
I1( stress );
463 const T I2_ =
I2( stress );
464 const T res = ( 1. / 3 ) * I1_ * I1_ - I2_;
482 template <
typename T >
483 T
J3(
const Eigen::Matrix< T, 6, 1 >& stress )
485 const T I1_ =
I1( stress );
486 const T I2_ =
I2( stress );
487 const T I3_ =
I3( stress );
489 return ( 2. / 27 ) * pow( I1_, 3 ) - ( 1. / 3 ) * I1_ * I2_ + I3_;
499 const Eigen::Matrix< double, 6, 1 >& S );
503 namespace Derivatives {
516 template <
typename T >
517 Eigen::Matrix< T, 6, 1 >
dRho_dStress( T rho,
const Eigen::Matrix< T, 6, 1 >& stress )
520 return Eigen::Matrix< T, 6, 1 >::Zero();
522 Eigen::Matrix< T, 6, 1 > s =
IDev * stress;
524 return T( 1. / rho ) *
P.array() * s.array();
645 namespace Transformations {
Eigen::Vector3d principalStresses(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:103
constexpr int voigtSize
Definition: MarmotFiniteElement.h:113
const Marmot::Vector6d PInv
Definition: MarmotVoigt.cpp:17
double I1Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:245
constexpr VoigtSize voigtSizeFromDimension(int x)
Definition: MarmotVoigt.h:47
double vonMisesEquivalentStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:224
RowVector6d dDeltaEpvneg_dE(const Marmot::Vector6d &dEp, const Matrix6d &CelInv, const Matrix6d &Cep)
Definition: MarmotVoigt.cpp:476
double normStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:232
Eigen::Vector3d sortedPrincipalStrains(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:109
Eigen::Matrix< T, 6, 1 > stressToVoigt(const Eigen::Matrix< T, 3, 3 > &stressTensor)
Definition: MarmotVoigt.h:261
Marmot::Matrix36 dSortedStrainPrincipal_dStrain(const Marmot::Vector6d &dEp)
Definition: MarmotVoigt.cpp:431
T normStrain(const Eigen::Matrix< T, 6, 1 > &strain)
Definition: MarmotVoigt.h:377
Eigen::Matrix< double, 6, 6 > Matrix6d
Definition: MarmotTypedefs.h:35
Eigen::Matrix< T, 6, 1 > dRho_dStress(T rho, const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:517
Marmot::Matrix36 dStressPrincipals_dStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:390
Vector6d dStressMean_dStress()
Definition: MarmotVoigt.cpp:281
double dThetaStrain_dJ2Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:335
@ Axial
Definition: MarmotVoigt.h:45
Marmot::Vector6d make3DVoigt(const Eigen::Matrix< double, voigtSize, 1 > &Voigt)
Definition: MarmotVoigt.h:152
constexpr int nDim
Definition: MarmotFiniteElement.h:112
T I1(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:402
T J2(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:460
Marmot::Vector6d dJ3Strain_dStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:376
Eigen::Vector3d voigtToPlaneVoigt(const Marmot::Vector6d &voigt)
Definition: MarmotVoigt.cpp:74
Eigen::Matrix< double, 3, 1 > Vector3d
Definition: MarmotTypedefs.h:42
double makeReal(const double &value)
Definition: MarmotMath.cpp:37
@ ThreeD
Definition: MarmotVoigt.h:45
Marmot::Vector6d strainToVoigt(const Eigen::Matrix3d &strainTensor)
#define VOIGTFROMDIM(x)
Definition: MarmotVoigt.h:35
double dTheta_dJ3(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:319
T I2(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:421
const Marmot::Vector6d I
Definition: MarmotVoigt.cpp:19
Eigen::Vector3d dStrainVolumetricNegative_dStrainPrincipal(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:410
Matrix6d dEp_dE(const Matrix6d &CelInv, const Matrix6d &Cep)
Definition: MarmotVoigt.cpp:421
Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > voigtFromStrainMatrix(const Eigen::Matrix< double, nDim, nDim > &strain)
Definition: MarmotVoigt.h:292
Marmot::Vector6d dTheta_dStress(double theta, const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:286
Eigen::Matrix< T, 3, 3 > voigtToStrain(const Eigen::Matrix< T, 6, 1 > &voigt)
Definition: MarmotVoigt.h:184
Eigen::Matrix< double, 3, 3 > Matrix3d
Definition: MarmotTypedefs.h:40
Eigen::Matrix< double, voigtSize, 1 > reduce3DVoigt(const Marmot::Vector6d &Voigt3D)
Definition: MarmotVoigt.h:119
T J3(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:483
This file includes functions needed for calculations with stress and strain tensors written in voigt ...
Definition: MarmotTesting.h:32
double StrainVolumetricNegative(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:237
double dTheta_dJ2(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:303
const Marmot::Vector6d IHyd
Definition: MarmotVoigt.cpp:20
RowVector6d dDeltaEpv_dE(const Matrix6d &CelInv, const Matrix6d &Cep)
Definition: MarmotVoigt.cpp:426
const Marmot::Vector6d P
Definition: MarmotVoigt.cpp:16
double J2Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:265
VoigtSize
Definition: MarmotVoigt.h:45
Eigen::Matrix< double, 6, 1 > Vector6d
Definition: MarmotTypedefs.h:43
Marmot::Vector6d dJ2Strain_dStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:371
double vonMisesEquivalentStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:219
Eigen::Matrix< double, 1, 6 > RowVector6d
Definition: MarmotTypedefs.h:48
Eigen::Matrix< T, 3, 3 > voigtToStress(const Eigen::Matrix< T, 6, 1 > &voigt)
Definition: MarmotVoigt.h:213
Eigen::Vector3d principalStrains(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:97
Eigen::Matrix3d principalStressesDirections(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:121
T I3(const Eigen::Matrix< T, 6, 1 > &stress)
Definition: MarmotVoigt.h:442
Marmot::Vector6d planeVoigtToVoigt(const Eigen::Vector3d &voigtPlane)
double J3Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:271
std::pair< Eigen::Vector3d, Eigen::Matrix< double, 3, 6 > > principalValuesAndDerivatives(const Eigen::Matrix< double, 6, 1 > &S)
Definition: MarmotVoigt.cpp:129
Marmot::Vector6d dJ2_dStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:359
@ OneD
Definition: MarmotVoigt.h:45
Marmot::Vector6d dThetaStrain_dStrain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:383
Eigen::Matrix< double, nDim, nDim > stressMatrixFromVoigt(const Eigen::Matrix< double, VOIGTFROMDIM(nDim), 1 > &Voigt)
Definition: MarmotVoigt.h:278
@ TwoD
Definition: MarmotVoigt.h:45
Eigen::Matrix< double, 6, 6 > stiffnessToVoigt(const Eigen::Tensor< double, 4 > &C)
Definition: MarmotVoigt.cpp:46
Marmot::Vector6d dJ3_dStress(const Marmot::Vector6d &stress)
Definition: MarmotVoigt.cpp:364
double I3Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:259
double I2Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:250
const Matrix6d IDev
Definition: MarmotVoigt.cpp:22
double dThetaStrain_dJ3Strain(const Marmot::Vector6d &strain)
Definition: MarmotVoigt.cpp:347
Eigen::Matrix< double, 3, 6 > Matrix36d
Definition: MarmotTypedefs.h:53
Definition: MarmotJournal.h:32
Eigen::Matrix< double, 3, 6 > Matrix36
Definition: MarmotTypedefs.h:54