37#ifndef MODEL_SOLID_LINEAR_HEXAGONAL_H_
38#define MODEL_SOLID_LINEAR_HEXAGONAL_H_
61 m = Eigen::AngleAxisd(phi2, Eigen::Vector3d::UnitX()) *
62 Eigen::AngleAxisd(Phi, Eigen::Vector3d::UnitZ()) *
63 Eigen::AngleAxisd(phi1, Eigen::Vector3d::UnitX());
64 Define(C11,C12,C13,C33,C44,m);
70 amrex::Real Ctmp[3][3][3][3];
73 for(
int i = 0; i < 3; i++)
74 for(
int j = 0; j < 3; j++)
75 for(
int k = 0; k < 3; k++)
76 for(
int l = 0; l < 3; l++)
78 if(i == j && j == k && k == l)
80 if (i==0) Ctmp[i][j][k][l] = C11;
81 else if (i==1) Ctmp[i][j][k][l] = C11;
82 else if (i==2) Ctmp[i][j][k][l] = C33;
84 else if (i==k && j==l)
86 if ((i==0 && j==1) || (i==1 && j==0))
87 Ctmp[i][j][k][l] = 0.5*(C11-C12);
89 Ctmp[i][j][k][l] = C44;
91 else if (i==j && k==l)
93 if ((i==0 && k==1) || (i==1 && k==0)) Ctmp[i][j][k][l] = C12;
94 else if ((i==0 && k==2) || (i==2 && k==0)) Ctmp[i][j][k][l] = C13;
95 else if ((i==1 && k==2) || (i==2 && k==1)) Ctmp[i][j][k][l] = C13;
97 else Ctmp[i][j][k][l] = 0.0;
99 for(
int p = 0; p < AMREX_SPACEDIM; p++)
100 for(
int q = 0; q < AMREX_SPACEDIM; q++)
101 for(
int s = 0; s < AMREX_SPACEDIM; s++)
102 for(
int t = 0;
t < AMREX_SPACEDIM;
t++)
105 for(
int i = 0; i < 3; i++)
106 for(
int j = 0; j < 3; j++)
107 for(
int k = 0; k < 3; k++)
108 for(
int l = 0; l < 3; l++)
109 ddw(p,q,s,
t) += R(p,i)*R(s,k)*Ctmp[i][j][k][l]*R(q,j)*R(
t,l);
124 virtual void Print(std::ostream &out)
const override
134 static Hexagonal Combine(
const std::vector<Hexagonal> &models,
const std::vector<Set::Scalar> &eta)
139 for (
unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n];
140 for (
unsigned int n = 0 ; n < models.size(); n++)
142 ret.
ddw += models[n].ddw * (eta[n] / etasum);
150 ret.
Define(0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0);
164 ret.
Define(C11,C12,C13,C33,C44,phi1,Phi,phi2);
170 Set::Scalar C11 = NAN, C12 = NAN, C13 = NAN, C33 = NAN, C44 = NAN;
187 value.
Define(C11,C12,C13,C33,C44,phi1,Phi,phi2);
190 #define OP_CLASS Hexagonal
191 #define OP_VARS X(ddw)
#define pp_query_required(...)
#define pp_query_default(...)
bool contains(std::string name)
Set::Scalar W(const Set::Matrix &gradu) const override
void Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44, Eigen::Matrix3d R)
static Hexagonal Random()
static const KinematicVariable kinvar
Hexagonal(Solid< Set::Sym::MajorMinor > base)
static AMREX_FORCE_INLINE Hexagonal Combine(const std::vector< Hexagonal > &models, const std::vector< Set::Scalar > &eta)
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > ddw
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > DDW(const Set::Matrix &) const override
virtual void Print(std::ostream &out) const override
void Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44, Set::Scalar phi1, Set::Scalar Phi, Set::Scalar phi2)
static Hexagonal Random(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44)
static void Parse(Hexagonal &value, IO::ParmParse &pp)
Set::Matrix DW(const Set::Matrix &gradu) const override
static const Set::Scalar Pi
A collection of data types and symmetry-reduced data structures.
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix