Alamo
Cubic.H
Go to the documentation of this file.
1#ifndef MODEL_SOLID_FINITE_PSEUDOLINEAR_CUBIC_H_
2#define MODEL_SOLID_FINITE_PSEUDOLINEAR_CUBIC_H_
3
4#include "IO/ParmParse.H"
5#include "Model/Solid/Solid.H"
6
7namespace Model
8{
9namespace Solid
10{
11namespace Finite
12{
13namespace PseudoLinear
14{
15class Cubic : public Solid<Set::Sym::Major>
16{
17public:
18 Cubic() {};
19 //Cubic(Solid<Set::Sym::Major> base) : Solid<Set::Sym::Major>(base) {};
20 virtual ~Cubic() {};
21
22 Set::Scalar W(const Set::Matrix & F) const override
23 {
24 q.normalize();
27 Set::Matrix E = 0.5 * (F.transpose() * F - Set::Matrix::Identity());
28 return 0.5 * (E * (CC*E).transpose()).trace();
29 }
30 Set::Matrix DW(const Set::Matrix & F) const override
31 {
32 q.normalize();
35 Set::Matrix C = F.transpose() * F;
36 Set::Matrix dw = Set::Matrix::Zero();
37 for (int i = 0; i < AMREX_SPACEDIM; i++)
38 for (int J = 0; J < AMREX_SPACEDIM; J++)
39 for (int A = 0; A < AMREX_SPACEDIM; A++)
40 for (int R = 0; R < AMREX_SPACEDIM; R++)
41 for (int S = 0; S < AMREX_SPACEDIM; S++)
42 {
43 dw(i,J) += 0.5 * F(i,A)*CC(J,A,R,S) * (C(R,S) - ((R==S) ? 1.0 : 0.0));
44 }
45 return dw;
46 }
48 {
49 q.normalize();
53 Set::Matrix C = F.transpose() * F;
54 Set::Matrix I = Set::Matrix::Identity();
55
56 for (int i = 0; i < AMREX_SPACEDIM; i++)
57 for (int J = 0; J < AMREX_SPACEDIM; J++)
58 for (int k = 0; k < AMREX_SPACEDIM; k++)
59 for (int L = 0; L < AMREX_SPACEDIM; L++)
60 {
61 ddw(i,J,k,L) = 0.0;
62 for (int A = 0; A < AMREX_SPACEDIM; A++)
63 for (int B = 0; B < AMREX_SPACEDIM; B++)
64 {
65 ddw(i,J,k,L) += 0.5*I(i,k)*CC(J,L,A,B) * (C(A,B) - I(A,B));
66 ddw(i,J,k,L) += F(i,A)*F(k,B)*CC(A,J,B,L);
67 if (std::isnan(ddw(i,J,k,L))) Util::Abort(INFO);
68 }
69 }
70
71 return ddw;
72 }
73 virtual void Print(std::ostream &out) const override
74 {
75 //out << "CC = " << CC;
76 out << "C11=" << C11 << " C12=" << C12 << " C44=" << C44 <<
77 " q=" << q.w() << " + " << q.x() << "i + " << q.y() << "j + " << q.z() << "k";
78 }
79
80public:
81 //Set::Matrix4<AMREX_SPACEDIM,Set::Sym::MajorMinor> CC;
83 Set::Scalar C11=1.68, C12=1.21, C44=0.75;
84 mutable Set::Quaternion q = Set::Quaternion(1.0, 0.0, 0.0, 0.0);
85
86public:
87 static Cubic Zero()
88 {
89 Cubic ret;
90 ret.C11 = 0.0;
91 ret.C12 = 0.0;
92 ret.C44 = 0.0;
93 ret.q = Set::Quaternion(0.0,0.0,0.0,0.0);
94 //ret.CC = Set::Matrix4<AMREX_SPACEDIM,Set::Sym::MajorMinor>::Zero();
95 return ret;
96 }
97 static Cubic Random()
98 {
100 }
102 {
103 Cubic ret;
104 ret.C11 = C11;
105 ret.C12 = C12;
106 ret.C44 = C44;
107 ret.q = Set::Quaternion::UnitRandom();
108 //Set::Scalar phi1 = 2.0*Set::Constant::Pi * Util::Random();
109 //Set::Scalar Phi = Set::Constant::Pi * Util::Random();
110 //Set::Scalar phi2 = 2.0*Set::Constant::Pi * Util::Random();
111 //ret.CC = Set::Matrix4<AMREX_SPACEDIM,Set::Sym::MajorMinor>::Cubic(C11,C12,C44,phi1,Phi,phi2);
112 //if (ret.CC.contains_nan()) Util::Abort(INFO);
113 return ret;
114 }
115 static void Parse(Cubic & value, IO::ParmParse & pp)
116 {
117 pp_query_default("C11",value.C11, 1.68); // Elastic constant
118 pp_query_default("C12",value.C12, 1.21); // Elastic constant
119 pp_query_default("C44",value.C44, 0.75); // Elastic constant
120
121 if (pp.contains("random"))
122 {
123 value = Cubic::Random(value.C11,value.C12,value.C44);
124 return;
125 }
126
127 Set::Scalar phi1 = 0.0, Phi = 0.0, phi2 = 0.0;
128 pp_query_default("phi1",phi1,0.0); // Bunge Euler angle :math:`\phi_1`
129 pp_query_default("Phi",Phi,0.0); // Bunge Euler angle :math:`\Phi`
130 pp_query_default("phi2",phi2,0.0); // Bunge Euler angle :math:`\phi_2`
131 Eigen::Matrix3d R;
132 R = Eigen::AngleAxisd(phi2, Eigen::Vector3d::UnitX()) *
133 Eigen::AngleAxisd(Phi, Eigen::Vector3d::UnitZ()) *
134 Eigen::AngleAxisd(phi1, Eigen::Vector3d::UnitX());
135 value.q = Set::Quaternion(R);
136 }
137
138#define OP_CLASS Cubic
139#define OP_VARS X(C11) X(C12) X(C44) X(q)
141};
143
144}
145}
146}
147}
148
149#endif
#define pp_query_default(...)
Definition ParmParse.H:100
#define INFO
Definition Util.H:20
bool contains(std::string name)
Definition ParmParse.H:154
Set::Scalar W(const Set::Matrix &F) const override
Definition Cubic.H:22
Set::Matrix DW(const Set::Matrix &F) const override
Definition Cubic.H:30
virtual void Print(std::ostream &out) const override
Definition Cubic.H:73
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Major > DDW(const Set::Matrix &F) const override
Definition Cubic.H:47
static void Parse(Cubic &value, IO::ParmParse &pp)
Definition Cubic.H:115
static constexpr KinematicVariable kinvar
Definition Cubic.H:82
static Cubic Random(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44)
Definition Cubic.H:101
KinematicVariable
Definition Solid.H:26
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:23
void Abort(const char *msg)
Definition Util.cpp:170
Set::Scalar Random()
Definition Set.cpp:9