Line data Source code
1 : #ifndef MODEL_SOLID_FINITE_NEOHOOKEANPREDEFORMED_H_
2 : #define MODEL_SOLID_FINITE_NEOHOOKEANPREDEFORMED_H_
3 :
4 : #include "IO/ParmParse.H"
5 : #include "Model/Solid/Finite/NeoHookean.H"
6 : #include "Unit/Unit.H"
7 :
8 : namespace Model
9 : {
10 : namespace Solid
11 : {
12 : namespace Finite
13 : {
14 : class NeoHookeanPredeformed : public NeoHookean
15 : {
16 : public:
17 10 : NeoHookeanPredeformed() = default;
18 126 : NeoHookeanPredeformed(NeoHookean base) : NeoHookean(base) {};
19 : ~NeoHookeanPredeformed() = default;
20 :
21 480 : Set::Scalar W(const Set::Matrix& F) const
22 : {
23 480 : return NeoHookean::W(F * F0.inverse());
24 : }
25 1960 : Set::Matrix DW(const Set::Matrix& F) const
26 : {
27 3920 : return NeoHookean::DW(F * F0.inverse());
28 : }
29 20 : Set::Matrix4<AMREX_SPACEDIM, Set::Sym::Major> DDW(const Set::Matrix& F) const
30 : {
31 20 : return NeoHookean::DDW(F * F0.inverse());
32 : }
33 : void Print(std::ostream& out) const
34 : {
35 : out << "mu = " << mu << " kappa = " << kappa << " F0 = " << F0;
36 : }
37 :
38 : public:
39 : Set::Matrix F0 = Set::Matrix::Identity();
40 :
41 : public:
42 2 : static NeoHookeanPredeformed Zero()
43 : {
44 2 : NeoHookeanPredeformed ret = NeoHookean::Zero();
45 2 : ret.F0 = Set::Matrix::Zero();
46 2 : return ret;
47 : }
48 62 : static NeoHookean Random()
49 : {
50 62 : NeoHookeanPredeformed ret = NeoHookean::Random();
51 62 : ret.F0 = Set::Matrix::Random();
52 62 : return ret;
53 : }
54 0 : static void Parse(NeoHookeanPredeformed& value, IO::ParmParse& pp)
55 : {
56 0 : pp.queryclass<NeoHookean>(value);
57 :
58 0 : if (pp.contains("eps0") && pp.contains("F0"))
59 : {
60 0 : Util::Abort("Cannot specify both F0 and eps0");
61 : }
62 0 : else if (pp.contains("F0"))
63 : {
64 0 : Set::Matrix F0;
65 : // Large-deformation eigendeformation (Identity = no deformation)
66 0 : pp_queryarr_default("F0", F0, Set::Matrix::Identity());
67 0 : value.F0 = F0;
68 : }
69 0 : else if (pp.contains("eps0"))
70 : {
71 0 : Set::Matrix eps0;
72 : // Small-deformation eigendeformation (Zero = no deformation)
73 0 : pp_queryarr_default("eps0",eps0, Set::Matrix::Zero());
74 0 : value.F0 = eps0 + Set::Matrix::Identity();
75 : }
76 : else
77 : {
78 0 : value.F0 = Set::Matrix::Identity();
79 : }
80 0 : Util::Assert(INFO,TEST(fabs(value.F0.determinant()) > 1E-8 ),"F0 must be non-singular");
81 0 : }
82 :
83 : #define OP_CLASS NeoHookeanPredeformed
84 : #define OP_VARS X(kappa) X(mu) X(F0)
85 : #include "Model/Solid/InClassOperators.H"
86 : };
87 : #include "Model/Solid/ExtClassOperators.H"
88 :
89 : }
90 : }
91 : }
92 :
93 :
94 :
95 : template<>
96 : ALAMO_SINGLE_DEFINITION
97 0 : int Set::Field<Model::Solid::Finite::NeoHookeanPredeformed>::NComp() const
98 : {
99 0 : return 2 + AMREX_SPACEDIM * AMREX_SPACEDIM;
100 : }
101 :
102 : template<>
103 : ALAMO_SINGLE_DEFINITION
104 0 : std::string Set::Field<Model::Solid::Finite::NeoHookeanPredeformed>::Name(int i) const
105 : {
106 0 : if (i == 0) return name + "_mu";
107 0 : if (i == 1) return name + "_kappa";
108 : #if AMREX_SPACEDIM==2
109 0 : if (i == 2) return name + "_F0xx";
110 0 : if (i == 3) return name + "_F0xy";
111 0 : if (i == 4) return name + "_F0yx";
112 0 : if (i == 5) return name + "_F0yy";
113 : #elif AMREX_SPACEDIM==3
114 : //Util::Abort(INFO, "Not implemented yet");
115 0 : if (i == 2) return name + "_F0xx";
116 0 : if (i == 3) return name + "_F0xy";
117 0 : if (i == 4) return name + "_F0yx";
118 0 : if (i == 5) return name + "_F0yy";
119 : #endif
120 0 : return name;
121 : }
122 :
123 : template<>
124 : ALAMO_SINGLE_DEFINITION
125 0 : void Set::Field<Model::Solid::Finite::NeoHookeanPredeformed>::Copy(int a_lev, amrex::MultiFab& a_dst, int a_dstcomp, int a_nghost) const
126 : {
127 0 : for (amrex::MFIter mfi(a_dst, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
128 : {
129 0 : const amrex::Box& bx = mfi.growntilebox(amrex::IntVect(a_nghost));
130 0 : if (bx.ok())
131 : {
132 0 : amrex::Array4<const Model::Solid::Finite::NeoHookeanPredeformed> const& src = ((*this)[a_lev])->array(mfi);
133 0 : amrex::Array4<Set::Scalar> const& dst = a_dst.array(mfi);
134 0 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
135 : {
136 0 : dst(i, j, k, a_dstcomp + 0) = src(i, j, k).mu;
137 0 : dst(i, j, k, a_dstcomp + 1) = src(i, j, k).kappa;
138 : #if AMREX_SPACEDIM==2
139 0 : dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(0, 0);
140 0 : dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(0, 1);
141 0 : dst(i, j, k, a_dstcomp + 4) = src(i, j, k).F0(1, 0);
142 0 : dst(i, j, k, a_dstcomp + 5) = src(i, j, k).F0(1, 1);
143 : #elif AMREX_SPACEDIM==3
144 0 : dst(i, j, k, a_dstcomp + 2) = src(i, j, k).F0(0, 0);
145 0 : dst(i, j, k, a_dstcomp + 3) = src(i, j, k).F0(0, 1);
146 0 : dst(i, j, k, a_dstcomp + 4) = src(i, j, k).F0(1, 0);
147 0 : dst(i, j, k, a_dstcomp + 5) = src(i, j, k).F0(1, 1);
148 : //Util::Abort(INFO, "Not implemented");
149 : #endif
150 :
151 : //dst(i, j, k, a_dstcomp + 0) = src(i, j, k).ddw(0, 0, 0, 0);
152 0 : });
153 : }
154 0 : }
155 0 : }
156 :
157 :
158 : #endif
|