Alamo
ThermoElastic.H
Go to the documentation of this file.
1 #ifndef INTEGRATOR_THERMOELASTIC_H
2 #define INTEGRATOR_THERMOELASTIC_H
3 
5 #include "Integrator/Mechanics.H"
7 #include "Numeric/Stencil.H"
8 
9 namespace Integrator
10 {
11 class ThermoElastic :
12  virtual public HeatConduction,
13  virtual public Mechanics<Model::Solid::Affine::Isotropic>
14 {
15 public:
17  HeatConduction(3),
18  Mechanics<Model::Solid::Affine::Isotropic>()
19  { }
21  {Parse(*this,pp);}
22  static void Parse(ThermoElastic &value, IO::ParmParse &pp)
23  {
24  pp_queryclass("hc",static_cast<HeatConduction*>(&value));
26 
27  pp_queryarr("alpha",value.alpha); // Diffusion coefficient
28  }
29 
30 protected:
31  void Initialize(int lev) override
32  {
35  }
36 
37  void UpdateModel(int a_step, Set::Scalar a_time) override
38  {
40 
41  //Set::Scalar alpha[2];
42  //alpha[0] = 0.001; alpha[1] = 0.002;
43 
44  for (int lev = 0; lev <= finest_level; ++lev)
45  {
46  for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
47  {
48  amrex::Box bx = mfi.grownnodaltilebox();
49  amrex::Array4<Model::Solid::Affine::Isotropic> const &model = model_mf[lev]->array(mfi);
50  amrex::Array4<const Set::Scalar> const &eta = eta_mf[lev]->array(mfi);
51  amrex::Array4<const Set::Scalar> const &temp = temp_old_mf[lev]->array(mfi);
52  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
53  Set::Matrix F0 = Set::Matrix::Zero();
55  for (int n = 0; n < eta.nComp(); n++)
56  {
57  F0 += (eta(i,j,k,n) * alpha[n]) * tempavg * Set::Matrix::Identity();
58  }
59  model(i, j, k).F0 = F0;
60  });
61  }
62 
63  Util::RealFillBoundary(*model_mf[lev], geom[lev]);
64  }
65 
66  }
67 
68  void TimeStepBegin(Set::Scalar a_time, int a_step) override
69  {
70  HeatConduction::TimeStepBegin(a_time, a_step);
72  }
73 
74  void Advance(int a_lev, amrex::Real a_time, amrex::Real a_dt) override
75  {
76  HeatConduction::Advance(a_lev, a_time, a_dt);
78  }
79 
80  void TagCellsForRefinement(int a_lev, amrex::TagBoxArray& a_tags, Set::Scalar a_time, int a_ngrow) override
81  {
82  HeatConduction::TagCellsForRefinement(a_lev, a_tags, a_time, a_ngrow);
84  }
85 
86  std::vector<Set::Scalar> alpha;
87 };
88 } // namespace Integrator
89 #endif
Numeric::Interpolate::CellToNodeAverage
static AMREX_FORCE_INLINE T CellToNodeAverage(const amrex::Array4< const T > &f, const int &i, const int &j, const int &k, const int &m, std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Definition: Stencil.H:1256
Util::RealFillBoundary
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T >> &a_mf, const amrex::Geometry &)
Definition: Util.H:307
Integrator::Base::Mechanics::TimeStepBegin
virtual void TimeStepBegin(Set::Scalar a_time, int a_step) override
Definition: Mechanics.H:160
Integrator::Mechanics
Definition: Mechanics.H:101
Integrator::HeatConduction::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int)
Definition: HeatConduction.H:148
Integrator::ThermoElastic::ThermoElastic
ThermoElastic(IO::ParmParse &pp)
Definition: ThermoElastic.H:20
Set::Isotropic
@ Isotropic
Definition: Base.H:197
Mechanics.H
Integrator::ThermoElastic::Initialize
void Initialize(int lev) override
Definition: ThermoElastic.H:31
Integrator::ThermoElastic::UpdateModel
void UpdateModel(int a_step, Set::Scalar a_time) override
Definition: ThermoElastic.H:37
Integrator::Mechanics< Model::Solid::Affine::Isotropic >::eta_mf
Set::Field< Set::Scalar > eta_mf
Definition: Mechanics.H:359
Integrator::Base::Mechanics< Model::Solid::Affine::Isotropic >::model_mf
Set::Field< Model::Solid::Affine::Isotropic > model_mf
Definition: Mechanics.H:453
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
Integrator::ThermoElastic::Parse
static void Parse(ThermoElastic &value, IO::ParmParse &pp)
Definition: ThermoElastic.H:22
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
Integrator::ThermoElastic::TimeStepBegin
void TimeStepBegin(Set::Scalar a_time, int a_step) override
Definition: ThermoElastic.H:68
Integrator::Mechanics::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Definition: Mechanics.H:314
Integrator::HeatConduction::Advance
void Advance(int lev, Set::Scalar, Set::Scalar dt)
Definition: HeatConduction.H:115
Integrator::Mechanics::UpdateModel
virtual void UpdateModel(int a_step, Set::Scalar time) override
Definition: Mechanics.H:219
Integrator::ThermoElastic::alpha
std::vector< Set::Scalar > alpha
Definition: ThermoElastic.H:86
Integrator::Integrator::TimeStepBegin
virtual void TimeStepBegin(Set::Scalar, int)
Run another system calculation (e.g. implicit solve) before integration step.
Definition: Integrator.H:189
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
Integrator::HeatConduction::Initialize
void Initialize(int lev)
Definition: HeatConduction.H:109
Integrator::ThermoElastic::ThermoElastic
ThermoElastic()
Definition: ThermoElastic.H:16
IO::ParmParse
Definition: ParmParse.H:110
Integrator::HeatConduction
Definition: HeatConduction.H:41
Integrator::Base::Mechanics::Advance
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition: Mechanics.H:238
Isotropic.H
Integrator::HeatConduction::temp_old_mf
Set::Field< Set::Scalar > temp_old_mf
Definition: HeatConduction.H:179
Integrator::ThermoElastic::Advance
void Advance(int a_lev, amrex::Real a_time, amrex::Real a_dt) override
Definition: ThermoElastic.H:74
Integrator::ThermoElastic
Definition: ThermoElastic.H:11
Integrator::Mechanics::Initialize
void Initialize(int lev) override
Definition: Mechanics.H:208
Model
Definition: Constant.H:12
HeatConduction.H
Integrator::ThermoElastic::TagCellsForRefinement
void TagCellsForRefinement(int a_lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Definition: ThermoElastic.H:80
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103