Alamo
Hexagonal.H
Go to the documentation of this file.
1 //
2 // This model implements an affine hexagonal material:
3 //
4 // .. math::
5 //
6 // W_{aff}(\nabla\mathbf{u}) = W_{lin}(\nabla\mathbf{u} - \mathbf{F}_0)
7 //
8 // where :math:`W_{aff}` is the present model, :math:`\mathbf{F}_0` is the eigenstrain,
9 // and :math:`W_{lin}` is the cubic model defined in :ref:`Model::Solid::Linear::Hexagonal`.
10 //
11 #ifndef MODEL_SOLID_AFFINE_HEXAGONAL_H_
12 #define MODEL_SOLID_AFFINE_HEXAGONAL_H_
13 
14 #include "Model/Solid/Solid.H"
16 #include "IO/ParmParse.H"
17 
18 namespace Model
19 {
20 namespace Solid
21 {
22 namespace Affine
23 {
25 {
26 public:
27 
28  Hexagonal() {};
29  Hexagonal(Linear::Hexagonal base) : Linear::Hexagonal(base) {};
30  virtual ~Hexagonal() {};
31 
32  void
33  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,Set::Matrix a_F0=Set::Matrix::Zero())
34  {
35  Linear::Hexagonal::Define(C11,C12,C13,C33,C44,phi1,Phi,phi2);
36  F0 = a_F0;
37  }
38  void
39  Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44, Eigen::Matrix3d R, Set::Matrix a_F0=Set::Matrix::Zero())
40  {
41  Linear::Hexagonal::Define(C11,C12,C13,C33,C44,R);
42  F0 = a_F0;
43  }
44  Set::Scalar W(const Set::Matrix & gradu) const override
45  {
46  return Linear::Hexagonal::W(gradu - F0);
47  }
48  Set::Matrix DW(const Set::Matrix & gradu) const override
49  {
50  return Linear::Hexagonal::DW(gradu - F0);
51  }
53  {
55  }
56  virtual void Print(std::ostream &out) const override
57  {
58  out << ddw;
59  }
60  AMREX_FORCE_INLINE
61  void SetF0 (Set::Matrix &a_F0) {F0 = a_F0;}
62 
63 public:
64  Set::Matrix F0 = Set::Matrix::Zero();
66 
67  static Hexagonal Zero()
68  {
70  ret.F0 = Set::Matrix::Zero();
71  return ret;
72  }
73  static Hexagonal Random()
74  {
76  }
78  {
79  Hexagonal ret = Linear::Hexagonal::Random(C11,C12,C13,C33,C44);
80  ret.F0 = Set::Matrix::Random();
81  return ret;
82  }
83  static void Parse(Hexagonal & value, IO::ParmParse & pp)
84  {
85  Linear::Hexagonal::Parse(value,pp);
86  pp_queryarr("F0",value.F0); // Eigenstrain matrix. Can be defined in 2D or 3D.
87  }
88 
89  AMREX_FORCE_INLINE
90  static Hexagonal Combine(const std::vector<Hexagonal> &models, const std::vector<Set::Scalar> &eta, int order)
91  {
92  Hexagonal ret;
94  ret.F0 = Set::Matrix::Zero();
95  if (order == 1)
96  {
97  Set::Scalar etasum = 0.;
98  for (unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n];
99  for (unsigned int n = 0 ; n < models.size(); n++)
100  {
101  ret.ddw += models[n].ddw * (eta[n] / etasum);
102  ret.F0 += models[n].F0 * (eta[n] / etasum);
103  }
104  return ret;
105  }
106  Util::Exception(INFO,"Order value of ", order, " not supported yet");
107  return ret;
108  }
109 
110  #define OP_CLASS Hexagonal
111  #define OP_VARS X(ddw) X(F0)
113 };
115 
116 }
117 }
118 }
119 
120 template<>
123 {
124  return 1;
125 }
126 
127 template<>
130 {
131  if (i==0) return name + ".Cxxxx";
132  return name;
133 }
134 
135 template<>
137 void Set::Field<Model::Solid::Affine::Hexagonal>::Copy(int a_lev, amrex::MultiFab &a_dst, int a_dstcomp, int a_nghost) const
138 {
139  for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
140  {
141  const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
142  if (bx.ok())
143  {
144  amrex::Array4<const Model::Solid::Affine::Hexagonal> const & src = ((*this)[a_lev])->array(mfi);
145  amrex::Array4<Set::Scalar> const & dst = a_dst.array(mfi);
146  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
147  {
148  dst(i,j,k,a_dstcomp + 0) = src(i,j,k).ddw(0,0,0,0);
149  });
150  }
151  }
152 }
153 
154 
155 #endif
Model::Solid::gradu
@ gradu
Definition: Solid.H:26
Model::Solid::Affine::Hexagonal::kinvar
static const KinematicVariable kinvar
Definition: Hexagonal.H:65
Model::Solid::Affine::Hexagonal::DW
Set::Matrix DW(const Set::Matrix &gradu) const override
Definition: Hexagonal.H:48
Model::Solid::Affine::Hexagonal::Random
static Hexagonal Random(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44)
Definition: Hexagonal.H:77
Model::Solid::Affine::Hexagonal::W
Set::Scalar W(const Set::Matrix &gradu) const override
Definition: Hexagonal.H:44
Model::Solid::Linear::Hexagonal::DDW
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > DDW(const Set::Matrix &) const override
Definition: Hexagonal.H:120
Model::Solid::Affine::Hexagonal::Random
static Hexagonal Random()
Definition: Hexagonal.H:73
ALAMO_SINGLE_DEFINITION
#define ALAMO_SINGLE_DEFINITION
Definition: Util.H:25
Model::Solid::Affine::Hexagonal::SetF0
AMREX_FORCE_INLINE void SetF0(Set::Matrix &a_F0)
Definition: Hexagonal.H:61
Model::Solid::Linear::Hexagonal
Definition: Hexagonal.H:49
Set::Field::Copy
void Copy(int, amrex::MultiFab &, int, int) const
Definition: Set.H:68
ParmParse.H
Model::Solid::Affine::Hexagonal::Zero
static Hexagonal Zero()
Definition: Hexagonal.H:67
Model::Solid::Affine::Hexagonal::Parse
static void Parse(Hexagonal &value, IO::ParmParse &pp)
Definition: Hexagonal.H:83
Model::Solid::Affine::Hexagonal::Define
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, Set::Matrix a_F0=Set::Matrix::Zero())
Definition: Hexagonal.H:33
Model::Solid::Linear::Hexagonal::Parse
static void Parse(Hexagonal &value, IO::ParmParse &pp)
Definition: Hexagonal.H:168
Model::Solid::Linear::Hexagonal::Define
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)
Definition: Hexagonal.H:58
Set::Field::Name
std::string Name(int) const
Definition: Set.H:73
Hexagonal.H
Util::Random
Set::Scalar Random()
Definition: Set.cpp:9
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
Model::Solid::Affine::Hexagonal::Print
virtual void Print(std::ostream &out) const override
Definition: Hexagonal.H:56
Model::Solid::Affine::Hexagonal::Hexagonal
Hexagonal(Linear::Hexagonal base)
Definition: Hexagonal.H:29
Solid.H
Model::Solid::Affine::Hexagonal
Definition: Hexagonal.H:24
Model::Solid::Affine::Hexagonal::DDW
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > DDW(const Set::Matrix &gradu) const override
Definition: Hexagonal.H:52
InClassOperators.H
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
ExtClassOperators.H
Model::Solid::Affine::Hexagonal::Combine
static AMREX_FORCE_INLINE Hexagonal Combine(const std::vector< Hexagonal > &models, const std::vector< Set::Scalar > &eta, int order)
Definition: Hexagonal.H:90
Util::Exception
void Exception(std::string file, std::string func, int line, Args const &... args)
Definition: Util.H:194
Model::Solid::Linear::Hexagonal::Random
static Hexagonal Random()
Definition: Hexagonal.H:154
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor >
Model::Solid::Affine::Hexagonal::Hexagonal
Hexagonal()
Definition: Hexagonal.H:28
Model::Solid::KinematicVariable
KinematicVariable
Definition: Solid.H:26
IO::ParmParse
Definition: ParmParse.H:110
Model::Solid::Affine::Hexagonal::~Hexagonal
virtual ~Hexagonal()
Definition: Hexagonal.H:30
Model::Solid::Affine::Hexagonal::F0
Set::Matrix F0
Definition: Hexagonal.H:64
Model::Solid::Linear::Hexagonal::DW
Set::Matrix DW(const Set::Matrix &gradu) const override
Definition: Hexagonal.H:116
Model::Solid::Linear::Hexagonal::W
Set::Scalar W(const Set::Matrix &gradu) const override
Definition: Hexagonal.H:112
Model::Solid::Affine::Hexagonal::Define
void Define(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44, Eigen::Matrix3d R, Set::Matrix a_F0=Set::Matrix::Zero())
Definition: Hexagonal.H:39
Set::Field::NComp
int NComp() const
Definition: Set.H:72
Model::Solid::Linear::Hexagonal::ddw
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::MajorMinor > ddw
Definition: Hexagonal.H:130
Model::Solid::Linear::Hexagonal::Zero
static Hexagonal Zero()
Definition: Hexagonal.H:147
INFO
#define INFO
Definition: Util.H:20
Model
Definition: Constant.H:12
Set::Matrix3d
Eigen::Matrix3d Matrix3d
Definition: Base.H:24
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103