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