Alamo
Cubic.H
Go to the documentation of this file.
1#ifndef MODEL_SOLID_FINITE_PSEUDOAFFINE_CUBIC_H_
2#define MODEL_SOLID_FINITE_PSEUDOAFFINE_CUBIC_H_
3
4#include "IO/ParmParse.H"
6
7namespace Model
8{
9namespace Solid
10{
11namespace Finite
12{
13namespace PseudoAffine
14{
16{
17public:
18 static constexpr const char *name = "cubic";
19
20 Cubic() {};
21 Cubic(PseudoLinear::Cubic base) : PseudoLinear::Cubic(base) {};
22 virtual ~Cubic() {};
23
24 Set::Scalar W(const Set::Matrix& F) const override
25 {
26 Set::Matrix F0inv = F0.inverse();
27 return PseudoLinear::Cubic::W(F * F0inv);
28 }
29 Set::Matrix DW(const Set::Matrix& F) const override
30 {
31 Set::Matrix F0inv = F0.inverse();
32 return PseudoLinear::Cubic::DW(F * F0inv) * F0inv.transpose();
33 }
35 {
36 Set::Matrix F0inv = F0.inverse();
37 auto DDW = PseudoLinear::Cubic::DDW(F * F0inv);
39 for (int i = 0; i < AMREX_SPACEDIM; i++)
40 for (int J = 0; J < AMREX_SPACEDIM; J++)
41 for (int k = 0; k < AMREX_SPACEDIM; k++)
42 for (int L = 0; L < AMREX_SPACEDIM; L++)
43 {
44 ret(i,J,k,L) = 0.0;
45 for (int Q = 0; Q < AMREX_SPACEDIM; Q++)
46 for (int R = 0; R < AMREX_SPACEDIM; R++)
47 ret(i,J,k,L) += DDW(i,Q,k,R)*F0inv(J,Q)*F0inv(L,R);
48 }
49 return ret;
50 }
51
52public:
53 Set::Matrix F0 = Set::Matrix::Identity();
54
55public:
56 static Cubic Zero()
57 {
59 ret.F0 = Set::Matrix::Zero();
60 return ret;
61 }
62 static Cubic Random()
63 {
65 ret.F0 = Set::Matrix::Random();
66 return ret;
67 }
68
69 AMREX_FORCE_INLINE
70 static Cubic Combine(const std::vector<Cubic> &models, const std::vector<Set::Scalar> &eta, int order)
71 {
72 Cubic ret = Cubic::Zero();
73 if (order == 1)
74 {
75 Set::Scalar etasum = 0.;
76 for (unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n];
77 for (unsigned int n = 0 ; n < models.size(); n++)
78 {
79 ret.C11 += models[n].C11 * (eta[n] / etasum);
80 ret.C12 += models[n].C12 * (eta[n] / etasum);
81 ret.C44 += models[n].C44 * (eta[n] / etasum);
82 ret.q += models[n].q * (eta[n] / etasum);
83 ret.F0 += models[n].F0 * (eta[n] / etasum);
84 }
85 return ret;
86 }
87 else if (order == 2)
88 {
89 Set::Scalar etasum = 0.;
90 for (unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n]*eta[n];
91 for (unsigned int n = 0 ; n < models.size(); n++)
92 {
93 ret.C11 += models[n].C11 * (eta[n]*eta[n] / etasum);
94 ret.C12 += models[n].C12 * (eta[n]*eta[n] / etasum);
95 ret.C44 += models[n].C44 * (eta[n]*eta[n] / etasum);
96 ret.q += models[n].q * (eta[n]*eta[n] / etasum);
97 ret.F0 += models[n].F0 * (eta[n]*eta[n] / etasum);
98 }
99 return ret;
100 }
101
102 Util::Exception(INFO,"invalid value used for order, ",order);
103 return ret; // should never happen
104 }
105
106
107 static void Parse(Cubic& value, IO::ParmParse& pp)
108 {
110 if (pp.contains("eps0") && pp.contains("F0"))
111 {
112 Util::Abort("Cannot specify both F0 and eps0");
113 }
114 else if (pp.contains("F0"))
115 {
117 pp_queryarr("F0", F0); // Large-strain eigendeformation (Identity = no deformation)
118 value.F0 = F0;
119 }
120 else if (pp.contains("eps0"))
121 {
122 Set::Matrix eps0;
123 pp_queryarr("eps0",eps0); // Small-strain eigendeformation (Zero = no deformation)
124 value.F0 = eps0 + Set::Matrix::Identity();
125 }
126 else
127 {
128 value.F0 = Set::Matrix::Identity();
129 }
130 Util::Assert(INFO,TEST(fabs(value.F0.determinant()) > 1E-8 ),"F0 must be non-singular");
131 }
132
133#define OP_CLASS Cubic
134#define OP_VARS X(C11) X(C12) X(C44) X(q) X(F0)
136};
138
139}
140}
141}
142}
143
144
145template<>
151
152template<>
155{
156 if (i == 0) return name + "_F0xx";
157 if (i == 1) return name + "_F0xy";
158 if (i == 2) return name + "_F0yx";
159 if (i == 3) return name + "_F0yy";
160 return name;
161}
162
163template<>
166{
167 for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
168 {
169 const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
170 if (bx.ok())
171 {
172 amrex::Array4<const Model::Solid::Finite::PseudoAffine::Cubic> const& src = ((*this)[a_lev])->array(mfi);
173 amrex::Array4<Set::Scalar> const& dst = a_dst.array(mfi);
174 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
175 {
176 dst(i, j, k, a_dstcomp + 0) = src(i, j, k).F0(0, 0);
177 dst(i, j, k, a_dstcomp + 1) = src(i, j, k).F0(0, 1);
178 dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(1, 0);
179 dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(1, 1);
180 });
181 }
182 }
183}
184
185
186#endif
#define pp_queryarr(...)
Definition ParmParse.H:103
#define TEST(x)
Definition Util.H:21
#define ALAMO_SINGLE_DEFINITION
Definition Util.H:25
#define INFO
Definition Util.H:20
bool contains(std::string name)
Definition ParmParse.H:154
static constexpr const char * name
Definition Cubic.H:18
Cubic(PseudoLinear::Cubic base)
Definition Cubic.H:21
static AMREX_FORCE_INLINE Cubic Combine(const std::vector< Cubic > &models, const std::vector< Set::Scalar > &eta, int order)
Definition Cubic.H:70
static void Parse(Cubic &value, IO::ParmParse &pp)
Definition Cubic.H:107
Set::Matrix DW(const Set::Matrix &F) const override
Definition Cubic.H:29
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Major > DDW(const Set::Matrix &F) const override
Definition Cubic.H:34
Set::Scalar W(const Set::Matrix &F) const override
Definition Cubic.H:24
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
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
std::string Name(int) const
Definition Set.H:73
void Copy(int, amrex::MultiFab &, int, int) const
Definition Set.H:68
int NComp() const
Definition Set.H:72
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
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:70
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:205