Alamo
Cubic.H
Go to the documentation of this file.
1//
2// This class extends :ref:`Model::Solid::Linear::Cubic`
3// by adding an eigenstrain :math:`\mathbf{F}_0`.
4//
5
6#ifndef MODEL_SOLID_AFFINE_CUBIC_H_
7#define MODEL_SOLID_AFFINE_CUBIC_H_
8
9#include "Model/Solid/Solid.H"
11#include "IO/ParmParse.H"
12
13namespace Model
14{
15namespace Solid
16{
17namespace Affine
18{
19class Cubic : public Linear::Cubic
20{
21public:
22 static constexpr const char *name = "cubic";
23
24 Cubic() = default;
25 //Cubic(Linear::Cubic base) : Linear::Cubic(base) {};
26 ~Cubic() = default;
27
28 void
29 Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44, Set::Scalar phi1, Set::Scalar Phi, Set::Scalar phi2,Set::Matrix a_F0=Set::Matrix::Zero())
30 {
31 Linear::Cubic::Define(C11,C12,C44,phi1,Phi,phi2);
32 F0 = a_F0;
33 }
34 void
35 Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44, Eigen::Matrix3d R, Set::Matrix a_F0=Set::Matrix::Zero())
36 {
37 Linear::Cubic::Define(C11,C12,C44,R);
38 F0 = a_F0;
39 }
41 {
42 return Linear::Cubic::W(gradu - F0);
43 }
45 {
46 return Linear::Cubic::DW(gradu - F0);
47 }
52 void Print(std::ostream &out) const
53 {
54 out << ddw;
55 }
56 AMREX_FORCE_INLINE
57 void SetF0 (Set::Matrix &a_F0) {F0 = a_F0;}
58
59public:
60 Set::Matrix F0 = Set::Matrix::Zero();
62
63 static Cubic Zero()
64 {
65 Cubic ret;
67 ret.F0 = Set::Matrix::Zero();
68 return ret;
69 }
70 static Cubic Random()
71 {
73 }
75 {
76 Cubic ret;
77 ret.ddw = Linear::Cubic::Random(C11,C12,C44).ddw;
78 ret.F0 = Set::Matrix::Random();
79 return ret;
80 }
81 // This class extends :ref:`Model::Solid::Linear::Cubic` by adding
82 // an eigenstrain. (See the Linear::Cubic class for more inputs for this model)
83 static void Parse(Cubic & value, IO::ParmParse & pp)
84 {
85 pp.queryclass<Linear::Cubic>(value);
86 if (pp.contains("F0"))
87 pp_queryarr("F0",value.F0); // Eigenstrain
88 }
89
90 AMREX_FORCE_INLINE
91 static Cubic Combine(const std::vector<Cubic> &models, const std::vector<Set::Scalar> &eta, int order)
92 {
93 Cubic ret;
95 ret.F0 = Set::Matrix::Zero();
96 if (order == 1)
97 {
98 Set::Scalar etasum = 0.;
99 for (unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n];
100 for (unsigned int n = 0 ; n < models.size(); n++)
101 {
102 ret.ddw += models[n].ddw * (eta[n] / etasum);
103 ret.F0 += models[n].F0 * (eta[n] / etasum);
104 }
105 return ret;
106 }
107 else if (order == 2)
108 {
109 Set::Scalar etasum = 0.;
110 for (unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n]*eta[n];
111 for (unsigned int n = 0 ; n < models.size(); n++)
112 {
113 ret.ddw += models[n].ddw * (eta[n]*eta[n] / etasum);
114 ret.F0 += models[n].F0 * (eta[n]*eta[n] / etasum);
115 }
116 return ret;
117 }
118
119 Util::Exception(INFO,"invalid value used for order, ",order);
120 return ret; // should never happen
121 }
122
123 #define OP_CLASS Cubic
124 #define OP_VARS X(ddw) X(F0)
126};
128
129}
130}
131}
132
133template<>
136{
137 return 1;
138}
139
140template<>
143{
144 if (i==0) return name + ".Cxxxx";
145 return name;
146}
147
148template<>
151{
152 for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
153 {
154 const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
155 if (bx.ok())
156 {
157 amrex::Array4<const Model::Solid::Affine::Cubic> const & src = ((*this)[a_lev])->array(mfi);
158 amrex::Array4<Set::Scalar> const & dst = a_dst.array(mfi);
159 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
160 {
161 dst(i,j,k,a_dstcomp + 0) = src(i,j,k).ddw(0,0,0,0);
162 });
163}
164 }
165}
166
167#endif
168
#define pp_queryarr(...)
Definition ParmParse.H:105
#define ALAMO_SINGLE_DEFINITION
Definition Util.H:28
#define INFO
Definition Util.H:23
bool contains(std::string name)
Definition ParmParse.H:173
void queryclass(std::string name, T *value)
Definition ParmParse.H:960
AMREX_FORCE_INLINE void SetF0(Set::Matrix &a_F0)
Definition Cubic.H:57
void Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44, Eigen::Matrix3d R, Set::Matrix a_F0=Set::Matrix::Zero())
Definition Cubic.H:35
static AMREX_FORCE_INLINE Cubic Combine(const std::vector< Cubic > &models, const std::vector< Set::Scalar > &eta, int order)
Definition Cubic.H:91
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > DDW(const Set::Matrix &gradu) const
Definition Cubic.H:48
static Cubic Zero()
Definition Cubic.H:63
void Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44, Set::Scalar phi1, Set::Scalar Phi, Set::Scalar phi2, Set::Matrix a_F0=Set::Matrix::Zero())
Definition Cubic.H:29
Set::Matrix DW(const Set::Matrix &gradu) const
Definition Cubic.H:44
static void Parse(Cubic &value, IO::ParmParse &pp)
Definition Cubic.H:83
static Cubic Random()
Definition Cubic.H:70
static Cubic Random(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44)
Definition Cubic.H:74
void Print(std::ostream &out) const
Definition Cubic.H:52
Set::Scalar W(const Set::Matrix &gradu) const
Definition Cubic.H:40
static const KinematicVariable kinvar
Definition Cubic.H:61
static constexpr const char * name
Definition Cubic.H:22
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > ddw
Definition Cubic.H:54
static Cubic Zero()
Definition Cubic.H:71
Set::Matrix DW(const Set::Matrix &gradu) const
Definition Cubic.H:40
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > DDW(const Set::Matrix &) const
Definition Cubic.H:44
static Cubic Random()
Definition Cubic.H:78
void Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C44, Set::Scalar phi1, Set::Scalar Phi, Set::Scalar phi2)
Definition Cubic.H:27
Set::Scalar W(const Set::Matrix &gradu) const
Definition Cubic.H:36
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
KinematicVariable
Definition Solid.H:34
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:22
Set::Scalar Random()
Definition Set.cpp:34
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:225