1#ifndef ICNS_DIFFUSION_H
2#define ICNS_DIFFUSION_H
19 PDEFields& fields,
const bool has_overset,
const bool mesh_mapping)
26 this->m_pdefields.field, amrex::Orientation::high));
27 this->m_applier->setDomainBC(
29 this->m_pdefields.field, amrex::Orientation::low),
31 this->m_pdefields.field, amrex::Orientation::high));
34 template <
typename Scheme>
39 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
42 auto& divtau = this->
m_pdefields.diff_term.state(tau_state);
45 mlmg.apply(divtau.vec_ptrs(), this->m_pdefields.field.vec_ptrs());
54 const bool has_overset,
55 const bool mesh_mapping,
56 const std::string& prefix =
"diffusion")
58 ,
m_density(fields.repo.get_field(
"density"))
62 amrex::LPInfo isolve =
m_options.lpinfo();
65 iapply.setMaxCoarseningLevel(0);
76 mesh.Geom(0, mesh.finestLevel()),
77 mesh.boxArray(0, mesh.finestLevel()),
78 mesh.DistributionMap(0, mesh.finestLevel()), isolve,
80 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
83 mesh.Geom(0, mesh.finestLevel()),
84 mesh.boxArray(0, mesh.finestLevel()),
85 mesh.DistributionMap(0, mesh.finestLevel()), iapply,
87 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
93 mesh.Geom(0, mesh.finestLevel()),
94 mesh.boxArray(0, mesh.finestLevel()),
95 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve,
97 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
100 mesh.Geom(0, mesh.finestLevel()),
101 mesh.boxArray(0, mesh.finestLevel()),
102 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply,
104 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
115 template <
typename Scheme>
118 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
122 const int nlevels = repo.num_active_levels();
123 const auto& geom = repo.mesh().Geom();
126 auto& divtau =
m_pdefields.diff_term.state(tau_state);
127 const auto& density =
m_density.state(fstate);
129 Field const* mesh_detJ =
132 std::unique_ptr<ScratchField> rho_times_detJ =
137 const amrex::Real alpha = 0.0;
138 const amrex::Real beta = -1.0;
141 for (
int lev = 0; lev < nlevels; ++lev) {
146 (*rho_times_detJ)(lev).setVal(0.0);
147 amrex::MultiFab::AddProduct(
148 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
157 geom[lev], viscosity(lev));
165 mlmg.apply(divtau.vec_ptrs(),
m_pdefields.field.vec_ptrs());
168 for (
int lev = 0; lev < nlevels; ++lev) {
169 const auto& divtau_arrs = divtau(lev).arrays();
170 const auto& rho_arrs = density(lev).const_arrays();
173 divtau(lev), [=] AMREX_GPU_DEVICE(
174 int nbx,
int i,
int j,
int k)
noexcept {
175 amrex::Real rhoinv = 1.0 / rho_arrs[nbx](i, j, k);
176 divtau_arrs[nbx](i, j, k, 0) *= rhoinv;
177 divtau_arrs[nbx](i, j, k, 1) *= rhoinv;
178 divtau_arrs[nbx](i, j, k, 2) *= rhoinv;
181 amrex::Gpu::streamSynchronize();
189 const auto& geom = repo.mesh().Geom();
191 const auto& density =
m_density.state(fstate);
192 const int nlevels = repo.num_active_levels();
193 const int ndim = field.num_comp();
194 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
196 Field const* mesh_detJ =
199 std::unique_ptr<ScratchField> rho_times_detJ =
204 const amrex::Real alpha = 1.0;
205 const amrex::Real beta = dt;
208 for (
int lev = 0; lev < nlevels; ++lev) {
213 (*rho_times_detJ)(lev).setVal(0.0);
214 amrex::MultiFab::AddProduct(
215 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
224 geom[lev], viscosity(lev));
232 for (
int lev = 0; lev < nlevels; ++lev) {
233 auto& rhs = (*rhs_ptr)(lev);
235 const auto& rhs_arrs = rhs.arrays();
236 const auto& fld_arrs = field(lev).const_arrays();
237 const auto& rho_arrs = density(lev).const_arrays();
240 rhs, amrex::IntVect(0), ndim,
241 [=] AMREX_GPU_DEVICE(
242 int nbx,
int i,
int j,
int k,
int n)
noexcept {
243 rhs_arrs[nbx](i, j, k, n) =
244 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
247 amrex::Gpu::streamSynchronize();
252 m_pdefields.field.vec_ptrs(), rhs_ptr->vec_const_ptrs(),
273 const bool has_overset,
274 const bool mesh_mapping,
275 const std::string& prefix =
"diffusion")
277 ,
m_density(fields.repo.get_field(
"density"))
281 amrex::LPInfo isolve =
m_options.lpinfo();
282 amrex::LPInfo iapply;
284 iapply.setMaxCoarseningLevel(0);
292 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
295 mesh.Geom(0, mesh.finestLevel()),
296 mesh.boxArray(0, mesh.finestLevel()),
297 mesh.DistributionMap(0, mesh.finestLevel()), isolve);
299 mesh.Geom(0, mesh.finestLevel()),
300 mesh.boxArray(0, mesh.finestLevel()),
301 mesh.DistributionMap(0, mesh.finestLevel()), iapply);
306 mesh.Geom(0, mesh.finestLevel()),
307 mesh.boxArray(0, mesh.finestLevel()),
308 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve);
310 mesh.Geom(0, mesh.finestLevel()),
311 mesh.boxArray(0, mesh.finestLevel()),
312 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply);
323 template <
typename Scheme>
326 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
330 const int nlevels = repo.num_active_levels();
331 const auto& geom = repo.mesh().Geom();
333 auto& divtau =
m_pdefields.diff_term.state(tau_state);
335 const auto& density =
m_density.state(fstate);
337 Field const* mesh_detJ =
340 std::unique_ptr<ScratchField> rho_times_detJ =
345 const amrex::Real alpha = 0.0;
346 const amrex::Real beta = -1.0;
348 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
351 for (
int lev = 0; lev < nlevels; ++lev) {
357 (*rho_times_detJ)(lev).setVal(0.0);
358 amrex::MultiFab::AddProduct(
359 (*rho_times_detJ)(lev), density(lev), 0,
360 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
366 lev, (*rho_times_detJ)(lev));
373 geom[lev], viscosity(lev));
380 lev, amrex::GetArrOfConstPtrs(b));
383 auto divtau_comp = divtau.subview(i);
387 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
391 for (
int lev = 0; lev < nlevels; ++lev) {
392 const auto& divtau_arrs = divtau(lev).arrays();
393 const auto& rho_arrs = density(lev).const_arrays();
396 divtau(lev), [=] AMREX_GPU_DEVICE(
397 int nbx,
int i,
int j,
int k)
noexcept {
398 amrex::Real rhoinv = 1.0 / rho_arrs[nbx](i, j, k);
399 divtau_arrs[nbx](i, j, k, 0) *= rhoinv;
400 divtau_arrs[nbx](i, j, k, 1) *= rhoinv;
401 divtau_arrs[nbx](i, j, k, 2) *= rhoinv;
404 amrex::Gpu::streamSynchronize();
413 const auto& density =
m_density.state(fstate);
414 const int nlevels = repo.num_active_levels();
415 const int ndim = field.num_comp();
416 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
418 const auto& geom = repo.mesh().Geom();
420 Field const* mesh_detJ =
423 std::unique_ptr<ScratchField> rho_times_detJ =
428 const amrex::Real alpha = 1.0;
429 const amrex::Real beta = dt;
431 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
434 for (
int lev = 0; lev < nlevels; ++lev) {
437 (*rho_times_detJ)(lev).setVal(0.0);
438 amrex::MultiFab::AddProduct(
439 (*rho_times_detJ)(lev), density(lev), 0,
440 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
455 geom[lev], viscosity(lev));
461 lev, amrex::GetArrOfConstPtrs(b));
466 for (
int lev = 0; lev < nlevels; ++lev) {
467 auto& rhs = (*rhs_ptr)(lev);
469 const auto& rhs_arrs = rhs.arrays();
470 const auto& fld_arrs = field(lev).const_arrays();
471 const auto& rho_arrs = density(lev).const_arrays();
474 rhs, amrex::IntVect(0), ndim,
475 [=] AMREX_GPU_DEVICE(
476 int nbx,
int i,
int j,
int k,
int n)
noexcept {
477 rhs_arrs[nbx](i, j, k, n) =
478 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
481 amrex::Gpu::streamSynchronize();
483 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
485 auto rhs_ptr_comp = rhs_ptr->subview(i);
491 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
495 field.name() + std::to_string(i) +
"_solve", mlmg);
505 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
507 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
514template <
typename Scheme>
523 "DiffusionOp invoked for scalar PDE type");
528 PDEFields& fields,
const bool has_overset,
const bool mesh_mapping)
530 bool use_tensor_op =
true;
532 amrex::ParmParse pp(fields.
field.
name() +
"_diffusion");
533 pp.query(
"use_tensor_operator", use_tensor_op);
538 "Tensor and segregated operators should not be enabled "
544 fields, has_overset, mesh_mapping);
548 std::make_unique<ICNSDiffScalarSegregatedOp>(
549 fields, has_overset, mesh_mapping);
552 fields, has_overset, mesh_mapping);
const std::string & name() const
Name of this field (including state information)
Definition Field.H:125
IntField & get_int_field(const std::string &name, const FieldState fstate=FieldState::New) const
Return a reference to an integer field.
Definition FieldRepo.cpp:276
amrex::Vector< const amrex::iMultiFab * > vec_const_ptrs() const noexcept
Definition IntField.cpp:46
PDEFields & m_pdefields
Definition DiffusionOps.H:100
std::unique_ptr< amrex::MLTensorOp > m_applier
Definition DiffusionOps.H:108
virtual void setup_operator(amrex::MLTensorOp &linop, const amrex::Real alpha, const amrex::Real beta, const FieldState fstate)
Definition DiffusionOps.cpp:54
std::unique_ptr< amrex::MLTensorOp > m_solver
Definition DiffusionOps.H:107
DiffSolverIface(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition DiffusionOps.cpp:9
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:116
ICNSDiffScalarOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:52
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:185
bool m_mesh_mapping
Definition icns_diffusion.H:262
std::unique_ptr< amrex::MLABecLaplacian > m_solver_scalar
Definition icns_diffusion.H:264
Field & m_density
Definition icns_diffusion.H:260
MLMGOptions m_options
Definition icns_diffusion.H:261
PDEFields & m_pdefields
Definition icns_diffusion.H:259
std::unique_ptr< amrex::MLABecLaplacian > m_applier_scalar
Definition icns_diffusion.H:265
MLMGOptions m_options
Definition icns_diffusion.H:502
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_solver_scalar
Definition icns_diffusion.H:506
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:324
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:408
PDEFields & m_pdefields
Definition icns_diffusion.H:500
Field & m_density
Definition icns_diffusion.H:501
ICNSDiffScalarSegregatedOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:271
bool m_mesh_mapping
Definition icns_diffusion.H:503
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_applier_scalar
Definition icns_diffusion.H:508
ICNSDiffTensorOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:18
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:35
FieldState
Definition FieldDescTypes.H:14
@ CELL
Cell-centered (default)
Definition FieldDescTypes.H:28
@ New
Same as FieldState::NP1.
Definition FieldDescTypes.H:20
void print_mlmg_info(const std::string &solve_name, const amrex::MLMG &mlmg)
Definition console_io.cpp:170
Definition AdvOp_Godunov.H:18
Definition console_io.cpp:27
void viscosity_to_uniform_space(amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > &b, const amr_wind::FieldRepo &repo, int lev)
Definition incflo_diffusion.cpp:205
amrex::Vector< amrex::Array< amrex::LinOpBCType, AMREX_SPACEDIM > > get_diffuse_tensor_bc(amr_wind::Field &velocity, amrex::Orientation::Side side) noexcept
amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > average_velocity_eta_to_faces(const amrex::Geometry &geom, amrex::MultiFab const &cc_eta)
Definition MLMGOptions.H:25
DiffusionOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:527
bool use_segregated_op
Definition icns_diffusion.H:525
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:570
std::unique_ptr< ICNSDiffScalarSegregatedOp > m_scalar_segregated_op
Definition icns_diffusion.H:519
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:557
std::unique_ptr< ICNSDiffScalarOp > m_scalar_op
Definition icns_diffusion.H:518
std::unique_ptr< ICNSDiffTensorOp > m_tensor_op
Definition icns_diffusion.H:517
static constexpr int ndim
Definition icns.H:40
Definition PDEFields.H:27
FieldRepo & repo
Reference to the field repository instance.
Definition PDEFields.H:31
Field & field
Solution variable (e.g., velocity, temperature)
Definition PDEFields.H:34