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