Line data Source code
1 : //
2 : // This implements a basic heat conduction method in Alamo.
3 : // The partial differential equation to be solved is
4 : //
5 : // .. math::
6 : //
7 : // \frac{\partial T}{\partial t} = \alpha\,\Delta T
8 : //
9 : // where :math:`T` is temperature, :math:`t` is time, and :math:`alpha`
10 : // is the thermal diffusivity.
11 : // Integration is performed explicitly in time using forward Euler, and
12 : // differentiation is performed using the finite difference method.
13 : //
14 :
15 : #ifndef INTEGRATOR_HEATCONDUCTION_H // Include guards
16 : #define INTEGRATOR_HEATCONDUCTION_H //
17 :
18 : // Standard library includes
19 : #include <iostream>
20 : #include <fstream>
21 : #include <iomanip>
22 :
23 : // AMReX Includes
24 : #include "AMReX.H"
25 : #include "AMReX_ParallelDescriptor.H"
26 : #include "AMReX_ParmParse.H"
27 :
28 : // Alamo Includes
29 : #include "IO/ParmParse.H"
30 : #include "Integrator/Integrator.H"
31 : #include "BC/Constant.H"
32 : #include "BC/Expression.H"
33 : #include "IC/IC.H"
34 : #include "IC/Sphere.H"
35 : #include "IC/Constant.H"
36 : #include "IC/Expression.H"
37 : #include "Numeric/Stencil.H"
38 :
39 : namespace Integrator
40 : {
41 : class HeatConduction : virtual public Integrator
42 : {
43 : public:
44 : // Empty constructor
45 5 : HeatConduction(int a_nghost = 2) :
46 : Integrator(),
47 5 : number_of_ghost_cells(a_nghost)
48 5 : {}
49 :
50 : // Constructor that triggers parse
51 3 : HeatConduction(IO::ParmParse& pp) : HeatConduction()
52 : {
53 3 : Parse(*this, pp);
54 3 : }
55 :
56 6 : virtual ~HeatConduction()
57 3 : {
58 3 : delete ic;
59 3 : delete bc;
60 6 : }
61 :
62 : // The Parse function initializes the HeatConduction object using
63 : // a parser, pp.
64 : // Note that this is a static function, which means it does not have
65 : // direct access to member variables. Instead, it initializes the variables
66 : // inside the argument, "value", and so all references to member items are
67 : // prefixed by "value."
68 5 : static void Parse(HeatConduction& value, IO::ParmParse& pp)
69 : {
70 : // Diffusion coefficient :math:`\alpha`.
71 : // *This is an example of a required input variable -
72 : // - program will terminate unless it is provided.*
73 30 : pp_query_required("heat.alpha", value.alpha);
74 :
75 : // Criterion for mesh refinement.
76 : // *This is an example of a default input variable.
77 : // The default value is provided here, not in the
78 : // definition of the variable.*
79 25 : pp_query_default("heat.refinement_threshold", value.refinement_threshold, 0.01);
80 :
81 : // Initial condition type.
82 10 : pp.select_default<IC::Constant,IC::Sphere,IC::Expression>("ic",value.ic,value.geom);
83 :
84 : // Select BC object for temperature
85 10 : pp.select_default<BC::Constant,BC::Expression>("bc.temp",value.bc,1);
86 :
87 : // Register the temperature and old temperature fields.
88 : // temp_mf and temp_old_mf are defined near the bottom of this Header file.
89 15 : value.RegisterNewFab(value.temp_mf, value.bc, value.number_of_components, value.number_of_ghost_cells, "Temp", true);
90 15 : value.RegisterNewFab(value.temp_old_mf, value.bc, value.number_of_components, value.number_of_ghost_cells, "Temp_old", false);
91 5 : }
92 :
93 : protected:
94 :
95 : // Use the ic object to initialize the temperature field
96 30 : void Initialize(int lev)
97 : {
98 30 : ic->Initialize(lev, temp_old_mf);
99 30 : }
100 :
101 : // Integrate the heat equation
102 12122 : void Advance(int lev, Set::Scalar /*time*/, Set::Scalar dt)
103 : {
104 : // Swap the old temp fab and the new temp fab so we use
105 : // the new one.
106 12122 : std::swap(*temp_mf[lev], *temp_old_mf[lev]);
107 :
108 : // Get the cell size corresponding to this level
109 12122 : const Set::Scalar* DX = geom[lev].CellSize();
110 :
111 : // Iterate over all of the patches on this level
112 439384 : for (amrex::MFIter mfi(*temp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
113 : {
114 : // Get the box (index dimensions) for this patch
115 427263 : const amrex::Box& bx = mfi.tilebox();
116 :
117 : // Get an array-accessible handle to the data on this patch.
118 : // amrex::Array4<const Set::Scalar> const& temp_old = (*temp_old_mf[lev]).array(mfi);
119 : // amrex::Array4<Set::Scalar> const& temp = temp_mf.Patch(lev,mfi);//(*temp_mf[lev]).array(mfi);
120 427263 : Set::Patch<const Set::Scalar> temp_old = temp_old_mf.Patch(lev,mfi);
121 427263 : Set::Patch<Set::Scalar> temp = temp_mf.Patch(lev,mfi);
122 :
123 : // Iterate over the grid on this patch
124 427263 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
125 : {
126 : // Do the physics!
127 : // Note that Numeric::Laplacian is an inlined function so there is no overhead.
128 : // You can calculate the derivatives yourself if you want.
129 225697527 : temp(i, j, k) = temp_old(i, j, k) + dt * alpha * Numeric::Laplacian(temp_old, i, j, k, 0, DX);
130 75232509 : });
131 12121 : }
132 12121 : }
133 :
134 : // Tag cells for mesh refinement based on temperature gradient
135 685 : void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar /*time*/, int /*ngrow*/)
136 : {
137 : // Get cell dimensions as done above.
138 685 : const Set::Scalar* DX = geom[lev].CellSize();
139 : // Calculate the diagonal.
140 685 : Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
141 :
142 : // Iterate over the patches on this level
143 8511 : for (amrex::MFIter mfi(*temp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
144 : {
145 : // Get the box and handles as done above.
146 7826 : const amrex::Box& bx = mfi.tilebox();
147 7826 : amrex::Array4<char> const& tags = a_tags.array(mfi);
148 7826 : amrex::Array4<Set::Scalar> const& temp = (*temp_mf[lev]).array(mfi);
149 :
150 : // Iterate over the grid as done above.
151 7826 : amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
152 : {
153 : // Calculate the temperature gradient.
154 909192 : Set::Vector grad = Numeric::Gradient(temp, i, j, k, 0, DX);
155 :
156 : // Is the gradient * cell_size too big? If so, then
157 : // mark this cell as needing refinement.
158 909192 : if (grad.lpNorm<2>() * dr > refinement_threshold)
159 744600 : tags(i, j, k) = amrex::TagBox::SET;
160 909192 : });
161 685 : }
162 685 : }
163 :
164 : protected:
165 : Set::Field<Set::Scalar> temp_mf; // Temperature field variable (current timestep)
166 : Set::Field<Set::Scalar> temp_old_mf; // Temperature field variable (previous timestep)
167 :
168 : private:
169 :
170 : //
171 : // Definition of parameters set only at instantiation by
172 : // constructors.
173 : //
174 : const int number_of_components = 1; // Number of components
175 : const int number_of_ghost_cells = 2; // Number of ghost cells
176 :
177 : //
178 : // Definition of user-determined variables.
179 : //
180 : // Instantiate all variables to NAN if possible.
181 : // Default values may be set in Parse using query_default.
182 : //
183 :
184 : Set::Scalar alpha = NAN; // Thermal diffusivity
185 : Set::Scalar refinement_threshold = NAN ; // Criterion for cell refinement
186 :
187 : //
188 : // Definition of user-determined pointer variables.
189 : //
190 : // These should be set to nullptr. Make sure that they are deleted
191 : // in the ~HeatConduction destructor.
192 : //
193 :
194 : IC::IC<Set::Scalar>* ic = nullptr; // Object used to initialize temperature field
195 : BC::BC<Set::Scalar>* bc = nullptr; // Object used to update temp field boundary ghost cells
196 : };
197 : } // namespace Integrator
198 : #endif
|