Line data Source code
1 : #ifndef INTEGRATOR_THERMOELASTIC_H
2 : #define INTEGRATOR_THERMOELASTIC_H
3 :
4 : #include "Model/Solid/Affine/Isotropic.H"
5 : #include "Integrator/Mechanics.H"
6 : #include "Integrator/HeatConduction.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:
16 2 : ThermoElastic():
17 : HeatConduction(3),
18 2 : Mechanics<Model::Solid::Affine::Isotropic>()
19 2 : { }
20 :
21 0 : ~ThermoElastic() = default;
22 :
23 2 : ThermoElastic(IO::ParmParse &pp) : ThermoElastic()
24 2 : {Parse(*this,pp);}
25 2 : static void Parse(ThermoElastic &value, IO::ParmParse &pp)
26 : {
27 12 : pp.queryclass<HeatConduction>("hc",value);
28 12 : pp.queryclass<Mechanics<Model::Solid::Affine::Isotropic>>("el",value);
29 :
30 10 : pp_queryarr("alpha",value.alpha); // Diffusion coefficient
31 2 : }
32 :
33 : protected:
34 16 : void Initialize(int lev) override
35 : {
36 16 : HeatConduction::Initialize(lev);
37 16 : Mechanics<Model::Solid::Affine::Isotropic>::Initialize(lev);
38 16 : }
39 :
40 305 : void UpdateModel(int a_step, Set::Scalar a_time) override
41 : {
42 305 : Mechanics<Model::Solid::Affine::Isotropic>::UpdateModel(a_step, a_time);
43 :
44 : //Set::Scalar alpha[2];
45 : //alpha[0] = 0.001; alpha[1] = 0.002;
46 :
47 1525 : for (int lev = 0; lev <= finest_level; ++lev)
48 : {
49 51713 : for (MFIter mfi(*model_mf[lev], false); mfi.isValid(); ++mfi)
50 : {
51 50493 : amrex::Box bx = mfi.grownnodaltilebox();
52 50493 : amrex::Array4<Model::Solid::Affine::Isotropic> const &model = model_mf[lev]->array(mfi);
53 50493 : amrex::Array4<const Set::Scalar> const &eta = eta_mf[lev]->array(mfi);
54 50493 : amrex::Array4<const Set::Scalar> const &temp = temp_old_mf[lev]->array(mfi);
55 50493 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
56 17238917 : Set::Matrix F0 = Set::Matrix::Zero();
57 17238917 : Set::Scalar tempavg = Numeric::Interpolate::CellToNodeAverage(temp,i,j,k,0);
58 103433502 : for (int n = 0; n < eta.nComp(); n++)
59 : {
60 68955668 : F0 += (eta(i,j,k,n) * alpha[n]) * tempavg * Set::Matrix::Identity();
61 : }
62 34477834 : model(i, j, k).F0 = F0;
63 17238917 : });
64 1220 : }
65 :
66 1220 : Util::RealFillBoundary(*model_mf[lev], geom[lev]);
67 : }
68 :
69 305 : }
70 :
71 305 : void TimeStepBegin(Set::Scalar a_time, int a_step) override
72 : {
73 305 : HeatConduction::TimeStepBegin(a_time, a_step);
74 305 : Mechanics<Model::Solid::Affine::Isotropic>::TimeStepBegin(a_time, a_step);
75 305 : }
76 :
77 4562 : void Advance(int a_lev, amrex::Real a_time, amrex::Real a_dt) override
78 : {
79 4562 : HeatConduction::Advance(a_lev, a_time, a_dt);
80 4561 : Mechanics<Model::Solid::Affine::Isotropic>::Advance(a_lev, a_time, a_dt);
81 4560 : }
82 :
83 102 : void TagCellsForRefinement(int a_lev, amrex::TagBoxArray& a_tags, Set::Scalar a_time, int a_ngrow) override
84 : {
85 102 : HeatConduction::TagCellsForRefinement(a_lev, a_tags, a_time, a_ngrow);
86 102 : Mechanics<Model::Solid::Affine::Isotropic>::TagCellsForRefinement(a_lev, a_tags, a_time, a_ngrow);
87 102 : }
88 :
89 : std::vector<Set::Scalar> alpha;
90 : };
91 : } // namespace Integrator
92 : #endif
|