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