Line data Source code
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"
15 : #include "Model/Solid/Linear/Hexagonal.H"
16 : #include "IO/ParmParse.H"
17 :
18 : namespace Model
19 : {
20 : namespace Solid
21 : {
22 : namespace Affine
23 : {
24 : class Hexagonal : public Linear::Hexagonal
25 : {
26 : public:
27 : static constexpr const char *name = "hexagonal";
28 :
29 859 : Hexagonal() {};
30 172 : Hexagonal(Linear::Hexagonal base) : Linear::Hexagonal(base) {};
31 1184 : virtual ~Hexagonal() {};
32 :
33 : void
34 : 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())
35 : {
36 : Linear::Hexagonal::Define(C11,C12,C13,C33,C44,phi1,Phi,phi2);
37 : F0 = a_F0;
38 : }
39 : void
40 : 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())
41 : {
42 : Linear::Hexagonal::Define(C11,C12,C13,C33,C44,R);
43 : F0 = a_F0;
44 : }
45 260 : Set::Scalar W(const Set::Matrix & gradu) const override
46 : {
47 260 : return Linear::Hexagonal::W(gradu - F0);
48 : }
49 26960 : Set::Matrix DW(const Set::Matrix & gradu) const override
50 : {
51 53920 : return Linear::Hexagonal::DW(gradu - F0);
52 : }
53 12520 : Set::Matrix4<AMREX_SPACEDIM,Set::Sym::MajorMinor> DDW(const Set::Matrix & gradu) const override
54 : {
55 12520 : return Linear::Hexagonal::DDW(gradu - F0);
56 : }
57 0 : virtual void Print(std::ostream &out) const override
58 : {
59 0 : out << ddw;
60 0 : }
61 : AMREX_FORCE_INLINE
62 : void SetF0 (Set::Matrix &a_F0) {F0 = a_F0;}
63 :
64 : public:
65 : Set::Matrix F0 = Set::Matrix::Zero();
66 : static const KinematicVariable kinvar = KinematicVariable::gradu;
67 :
68 128 : static Hexagonal Zero()
69 : {
70 128 : Hexagonal ret = Linear::Hexagonal::Zero();
71 128 : ret.F0 = Set::Matrix::Zero();
72 128 : return ret;
73 0 : }
74 44 : static Hexagonal Random()
75 : {
76 44 : return Random(Util::Random(), Util::Random(), Util::Random(),Util::Random(),Util::Random());
77 : }
78 44 : static Hexagonal Random(Set::Scalar C11, Set::Scalar C12, Set::Scalar C13, Set::Scalar C33, Set::Scalar C44)
79 : {
80 44 : Hexagonal ret = Linear::Hexagonal::Random(C11,C12,C13,C33,C44);
81 44 : ret.F0 = Set::Matrix::Random();
82 44 : return ret;
83 0 : }
84 1 : static void Parse(Hexagonal & value, IO::ParmParse & pp)
85 : {
86 4 : pp.queryclass<Linear::Hexagonal>(value);
87 5 : pp_queryarr("F0",value.F0); // Eigenstrain matrix. Can be defined in 2D or 3D.
88 1 : }
89 :
90 : AMREX_FORCE_INLINE
91 : static Hexagonal Combine(const std::vector<Hexagonal> &models, const std::vector<Set::Scalar> &eta, int order)
92 : {
93 0 : Hexagonal ret;
94 0 : ret.ddw = Set::Matrix4<AMREX_SPACEDIM,Set::Sym::MajorMinor>::Zero();
95 0 : ret.F0 = Set::Matrix::Zero();
96 0 : if (order == 1)
97 : {
98 0 : Set::Scalar etasum = 0.;
99 0 : for (unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n];
100 0 : for (unsigned int n = 0 ; n < models.size(); n++)
101 : {
102 0 : ret.ddw += models[n].ddw * (eta[n] / etasum);
103 0 : ret.F0 += models[n].F0 * (eta[n] / etasum);
104 : }
105 0 : return ret;
106 : }
107 0 : Util::Exception(INFO,"Order value of ", order, " not supported yet");
108 0 : return ret;
109 0 : }
110 :
111 : #define OP_CLASS Hexagonal
112 : #define OP_VARS X(ddw) X(F0)
113 : #include "Model/Solid/InClassOperators.H"
114 : };
115 : #include "Model/Solid/ExtClassOperators.H"
116 :
117 : }
118 : }
119 : }
120 :
121 : template<>
122 : ALAMO_SINGLE_DEFINITION
123 44 : int Set::Field<Model::Solid::Affine::Hexagonal>::NComp() const
124 : {
125 44 : return 1;
126 : }
127 :
128 : template<>
129 : ALAMO_SINGLE_DEFINITION
130 11 : std::string Set::Field<Model::Solid::Affine::Hexagonal>::Name(int i) const
131 : {
132 11 : if (i==0) return name + ".Cxxxx";
133 0 : return name;
134 : }
135 :
136 : template<>
137 : ALAMO_SINGLE_DEFINITION
138 11 : void Set::Field<Model::Solid::Affine::Hexagonal>::Copy(int a_lev, amrex::MultiFab &a_dst, int a_dstcomp, int a_nghost) const
139 : {
140 22 : for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
141 : {
142 11 : const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
143 11 : if (bx.ok())
144 : {
145 11 : amrex::Array4<const Model::Solid::Affine::Hexagonal> const & src = ((*this)[a_lev])->array(mfi);
146 11 : amrex::Array4<Set::Scalar> const & dst = a_dst.array(mfi);
147 11 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
148 : {
149 4125 : dst(i,j,k,a_dstcomp + 0) = src(i,j,k).ddw(0,0,0,0);
150 1375 : });
151 : }
152 11 : }
153 11 : }
154 :
155 :
156 : #endif
|