1#ifndef MODEL_SOLID_FINITE_PSEUDOAFFINE_CUBIC_H_
2#define MODEL_SOLID_FINITE_PSEUDOAFFINE_CUBIC_H_
18 static constexpr const char *
name =
"cubic";
39 for (
int i = 0; i < AMREX_SPACEDIM; i++)
40 for (
int J = 0; J < AMREX_SPACEDIM; J++)
41 for (
int k = 0; k < AMREX_SPACEDIM; k++)
42 for (
int L = 0; L < AMREX_SPACEDIM; L++)
45 for (
int Q = 0; Q < AMREX_SPACEDIM; Q++)
46 for (
int R = 0; R < AMREX_SPACEDIM; R++)
47 ret(i,J,k,L) +=
DDW(i,Q,k,R)*F0inv(J,Q)*F0inv(L,R);
59 ret.
F0 = Set::Matrix::Zero();
65 ret.
F0 = Set::Matrix::Random();
70 static Cubic Combine(
const std::vector<Cubic> &models,
const std::vector<Set::Scalar> &eta,
int order)
76 for (
unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n];
77 for (
unsigned int n = 0 ; n < models.size(); n++)
79 ret.
C11 += models[n].C11 * (eta[n] / etasum);
80 ret.
C12 += models[n].C12 * (eta[n] / etasum);
81 ret.
C44 += models[n].C44 * (eta[n] / etasum);
82 ret.
q += models[n].q * (eta[n] / etasum);
83 ret.
F0 += models[n].F0 * (eta[n] / etasum);
90 for (
unsigned int n = 0 ; n < models.size(); n++) etasum += eta[n]*eta[n];
91 for (
unsigned int n = 0 ; n < models.size(); n++)
93 ret.
C11 += models[n].C11 * (eta[n]*eta[n] / etasum);
94 ret.
C12 += models[n].C12 * (eta[n]*eta[n] / etasum);
95 ret.
C44 += models[n].C44 * (eta[n]*eta[n] / etasum);
96 ret.
q += models[n].q * (eta[n]*eta[n] / etasum);
97 ret.
F0 += models[n].F0 * (eta[n]*eta[n] / etasum);
124 value.
F0 = eps0 + Set::Matrix::Identity();
128 value.
F0 = Set::Matrix::Identity();
133#define OP_CLASS Cubic
134#define OP_VARS X(C11) X(C12) X(C44) X(q) X(F0)
156 if (
i == 0)
return name +
"_F0xx";
157 if (
i == 1)
return name +
"_F0xy";
158 if (
i == 2)
return name +
"_F0yx";
159 if (
i == 3)
return name +
"_F0yy";
167 for (amrex::MFIter
mfi(
a_dst, amrex::TilingIfNotGPU());
mfi.isValid(); ++
mfi)
169 const amrex::Box&
bx =
mfi.growntilebox(amrex::IntVect(
a_nghost));
172 amrex::Array4<const Model::Solid::Finite::PseudoAffine::Cubic>
const&
src = ((*this)[
a_lev])->
array(
mfi);
173 amrex::Array4<Set::Scalar>
const&
dst =
a_dst.array(
mfi);
#define ALAMO_SINGLE_DEFINITION
bool contains(std::string name)
static constexpr const char * name
Cubic(PseudoLinear::Cubic base)
static AMREX_FORCE_INLINE Cubic Combine(const std::vector< Cubic > &models, const std::vector< Set::Scalar > &eta, int order)
static void Parse(Cubic &value, IO::ParmParse &pp)
Set::Matrix DW(const Set::Matrix &F) const override
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Major > DDW(const Set::Matrix &F) const override
Set::Scalar W(const Set::Matrix &F) const override
Set::Scalar W(const Set::Matrix &F) const override
Set::Matrix DW(const Set::Matrix &F) const override
Set::Matrix4< AMREX_SPACEDIM, Set::Sym::Major > DDW(const Set::Matrix &F) const override
static void Parse(Cubic &value, IO::ParmParse &pp)
std::string Name(int) const
void Copy(int, amrex::MultiFab &, int, int) const
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
void Abort(const char *msg)
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
void Exception(std::string file, std::string func, int line, Args const &... args)