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