15#ifndef INTEGRATOR_HEATCONDUCTION_H
16#define INTEGRATOR_HEATCONDUCTION_H
19#include "AMReX_Array4.H"
20#include "AMReX_GpuComplex.H"
21#include "AMReX_MFIter.H"
23#include <AMReX_ParallelDescriptor.H>
24#include <AMReX_ParmParse.H>
40#include "Numeric/Stencil.H"
47 static constexpr const char*
name =
"heatconduction";
98 value.RegisterNewFab(value.temp_mf, value.bc, value.number_of_components, value.number_of_ghost_cells,
"Temp",
true);
99 if (value.method ==
"realspace")
101 value.RegisterNewFab(value.temp_old_mf, value.bc, value.number_of_components, value.number_of_ghost_cells,
"Temp_old",
false);
110 ic->Initialize(lev, temp_mf);
111 if (method ==
"realspace") ic->Initialize(lev, temp_old_mf);
118 if (method ==
"spectral")
120 AdvanceSpectral(lev,time,dt);
126 std::swap(*temp_mf[lev], *temp_old_mf[lev]);
132 for (amrex::MFIter mfi(*temp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
135 const amrex::Box &bx = mfi.tilebox();
142 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
147 temp(i, j, k) = temp_old(i, j, k) + dt * alpha *
Numeric::Laplacian(temp_old, i, j, k, 0, DX);
159 amrex::Box
const & domain = this->geom[lev].Domain();
161 amrex::FFT::R2C my_fft(this->geom[lev].Domain());
162 auto const& [cba, cdm] = my_fft.getSpectralDataLayout();
163 amrex::FabArray<amrex::BaseFab<amrex::GpuComplex<Set::Scalar> > > Temp_hat(cba, cdm, 1, 0);
164 my_fft.forward(*temp_mf[lev], Temp_hat);
174 Set::Scalar scaling = 1.0 / geom[lev].Domain().d_numPts();
176 for (amrex::MFIter mfi(Temp_hat, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
178 const amrex::Box &bx = mfi.tilebox();
179 amrex::Array4<amrex::GpuComplex<Set::Scalar>>
const & T_hat = Temp_hat.array(mfi);
180 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int m,
int n,
int p) {
184 Set::Scalar k2 = (n < domain.length(1)/2 ? n * pi_Ly : (n - domain.length(1)) * pi_Ly);,
185 Set::Scalar k3 = (p < domain.length(2)/2 ? p * pi_Lz : (p - domain.length(2)) * pi_Lz););
187 Set::Scalar lap = AMREX_D_TERM(k1 * k1, + k2 * k2, + k3*k3);
190 T_hat(m,n,p) *= factor*scaling;
194 my_fft.backward(Temp_hat, *temp_mf[lev]);
208 if (method==
"spectral")
return;
213 Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
216 for (amrex::MFIter mfi(*temp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
219 const amrex::Box& bx = mfi.tilebox();
220 amrex::Array4<char>
const& tags = a_tags.array(mfi);
221 amrex::Array4<Set::Scalar>
const& temp = (*temp_mf[lev]).array(mfi);
224 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
231 if (grad.lpNorm<2>() * dr > refinement_threshold)
232 tags(i, j, k) = amrex::TagBox::SET;
249 const int number_of_components = 1;
250 const int number_of_ghost_cells = 2;
#define pp_query_required(...)
#define pp_query_default(...)
Pure abstract IC object from which all other IC objects inherit.
void select_default(std::string name, PTRTYPE *&ic_eta, Args &&... args)
int query_validate(std::string name, int &value, std::vector< int > possibleintvals, std::string file="", std::string func="", int line=-1)
const int number_of_ghost_cells
Set::Field< Set::Scalar > temp_mf
static void Parse(HeatConduction &value, IO::ParmParse &pp)
BC::BC< Set::Scalar > * bc
HeatConduction(int a_nghost=2)
void Advance(int lev, Set::Scalar time, Set::Scalar dt)
IC::IC< Set::Scalar > * ic
HeatConduction(IO::ParmParse &pp)
void AdvanceSpectral(int, Set::Scalar, Set::Scalar)
virtual ~HeatConduction()
Set::Scalar refinement_threshold
Set::Field< Set::Scalar > temp_old_mf
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int)
static constexpr const char * name
Collection of numerical integrator objects.
AMREX_FORCE_INLINE Set::Scalar Laplacian(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > &stencil=DefaultType)
AMREX_FORCE_INLINE Set::Vector Gradient(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
static const Set::Scalar Pi
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
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)