1#ifndef ICNS_DIFFUSION_H
2#define ICNS_DIFFUSION_H
12#include "AMReX_REAL.H"
14using namespace amrex::literals;
22 PDEFields& fields,
const bool has_overset,
const bool mesh_mapping)
29 this->m_pdefields.field, amrex::Orientation::high));
30 this->m_applier->setDomainBC(
32 this->m_pdefields.field, amrex::Orientation::low),
34 this->m_pdefields.field, amrex::Orientation::high));
37 template <
typename Scheme>
42 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
45 auto& divtau = this->
m_pdefields.diff_term.state(tau_state);
48 mlmg.apply(divtau.vec_ptrs(), this->m_pdefields.field.vec_ptrs());
57 const bool has_overset,
58 const bool mesh_mapping,
59 const std::string& prefix =
"diffusion")
61 ,
m_density(fields.repo.get_field(
"density"))
65 amrex::LPInfo isolve =
m_options.lpinfo();
68 iapply.setMaxCoarseningLevel(0);
79 mesh.Geom(0, mesh.finestLevel()),
80 mesh.boxArray(0, mesh.finestLevel()),
81 mesh.DistributionMap(0, mesh.finestLevel()), isolve,
83 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
86 mesh.Geom(0, mesh.finestLevel()),
87 mesh.boxArray(0, mesh.finestLevel()),
88 mesh.DistributionMap(0, mesh.finestLevel()), iapply,
90 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
96 mesh.Geom(0, mesh.finestLevel()),
97 mesh.boxArray(0, mesh.finestLevel()),
98 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve,
100 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
103 mesh.Geom(0, mesh.finestLevel()),
104 mesh.boxArray(0, mesh.finestLevel()),
105 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply,
107 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
118 template <
typename Scheme>
121 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
125 const int nlevels = repo.num_active_levels();
126 const auto& geom = repo.mesh().Geom();
129 auto& divtau =
m_pdefields.diff_term.state(tau_state);
130 const auto& density =
m_density.state(fstate);
132 Field const* mesh_detJ =
135 std::unique_ptr<ScratchField> rho_times_detJ =
140 const amrex::Real alpha = 0.0_rt;
141 const amrex::Real beta = -1.0_rt;
144 for (
int lev = 0; lev < nlevels; ++lev) {
149 (*rho_times_detJ)(lev).setVal(0.0_rt);
150 amrex::MultiFab::AddProduct(
151 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
160 geom[lev], viscosity(lev));
168 mlmg.apply(divtau.vec_ptrs(),
m_pdefields.field.vec_ptrs());
171 for (
int lev = 0; lev < nlevels; ++lev) {
172 const auto& divtau_arrs = divtau(lev).arrays();
173 const auto& rho_arrs = density(lev).const_arrays();
176 divtau(lev), [=] AMREX_GPU_DEVICE(
177 int nbx,
int i,
int j,
int k)
noexcept {
178 amrex::Real rhoinv = 1.0_rt / rho_arrs[nbx](i, j, k);
179 divtau_arrs[nbx](i, j, k, 0) *= rhoinv;
180 divtau_arrs[nbx](i, j, k, 1) *= rhoinv;
181 divtau_arrs[nbx](i, j, k, 2) *= rhoinv;
184 amrex::Gpu::streamSynchronize();
192 const auto& geom = repo.mesh().Geom();
194 const auto& density =
m_density.state(fstate);
195 const int nlevels = repo.num_active_levels();
196 const int ndim = field.num_comp();
197 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
199 Field const* mesh_detJ =
202 std::unique_ptr<ScratchField> rho_times_detJ =
207 const amrex::Real alpha = 1.0_rt;
208 const amrex::Real beta = dt;
211 for (
int lev = 0; lev < nlevels; ++lev) {
216 (*rho_times_detJ)(lev).setVal(0.0_rt);
217 amrex::MultiFab::AddProduct(
218 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
227 geom[lev], viscosity(lev));
235 for (
int lev = 0; lev < nlevels; ++lev) {
236 auto& rhs = (*rhs_ptr)(lev);
238 const auto& rhs_arrs = rhs.arrays();
239 const auto& fld_arrs = field(lev).const_arrays();
240 const auto& rho_arrs = density(lev).const_arrays();
243 rhs, amrex::IntVect(0), ndim,
244 [=] AMREX_GPU_DEVICE(
245 int nbx,
int i,
int j,
int k,
int n)
noexcept {
246 rhs_arrs[nbx](i, j, k, n) =
247 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
250 amrex::Gpu::streamSynchronize();
255 m_pdefields.field.vec_ptrs(), rhs_ptr->vec_const_ptrs(),
276 const bool has_overset,
277 const bool mesh_mapping,
278 const std::string& prefix =
"diffusion")
280 ,
m_density(fields.repo.get_field(
"density"))
284 amrex::LPInfo isolve =
m_options.lpinfo();
285 amrex::LPInfo iapply;
287 iapply.setMaxCoarseningLevel(0);
295 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
298 mesh.Geom(0, mesh.finestLevel()),
299 mesh.boxArray(0, mesh.finestLevel()),
300 mesh.DistributionMap(0, mesh.finestLevel()), isolve);
302 mesh.Geom(0, mesh.finestLevel()),
303 mesh.boxArray(0, mesh.finestLevel()),
304 mesh.DistributionMap(0, mesh.finestLevel()), iapply);
309 mesh.Geom(0, mesh.finestLevel()),
310 mesh.boxArray(0, mesh.finestLevel()),
311 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve);
313 mesh.Geom(0, mesh.finestLevel()),
314 mesh.boxArray(0, mesh.finestLevel()),
315 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply);
326 template <
typename Scheme>
329 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
333 const int nlevels = repo.num_active_levels();
334 const auto& geom = repo.mesh().Geom();
336 auto& divtau =
m_pdefields.diff_term.state(tau_state);
338 const auto& density =
m_density.state(fstate);
340 Field const* mesh_detJ =
343 std::unique_ptr<ScratchField> rho_times_detJ =
348 const amrex::Real alpha = 0.0_rt;
349 const amrex::Real beta = -1.0_rt;
351 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
354 for (
int lev = 0; lev < nlevels; ++lev) {
360 (*rho_times_detJ)(lev).setVal(0.0_rt);
361 amrex::MultiFab::AddProduct(
362 (*rho_times_detJ)(lev), density(lev), 0,
363 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
369 lev, (*rho_times_detJ)(lev));
376 geom[lev], viscosity(lev));
383 lev, amrex::GetArrOfConstPtrs(b));
386 auto divtau_comp = divtau.subview(i);
390 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
394 for (
int lev = 0; lev < nlevels; ++lev) {
395 const auto& divtau_arrs = divtau(lev).arrays();
396 const auto& rho_arrs = density(lev).const_arrays();
399 divtau(lev), [=] AMREX_GPU_DEVICE(
400 int nbx,
int i,
int j,
int k)
noexcept {
401 amrex::Real rhoinv = 1.0_rt / rho_arrs[nbx](i, j, k);
402 divtau_arrs[nbx](i, j, k, 0) *= rhoinv;
403 divtau_arrs[nbx](i, j, k, 1) *= rhoinv;
404 divtau_arrs[nbx](i, j, k, 2) *= rhoinv;
407 amrex::Gpu::streamSynchronize();
416 const auto& density =
m_density.state(fstate);
417 const int nlevels = repo.num_active_levels();
418 const int ndim = field.num_comp();
419 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
421 const auto& geom = repo.mesh().Geom();
423 Field const* mesh_detJ =
426 std::unique_ptr<ScratchField> rho_times_detJ =
431 const amrex::Real alpha = 1.0_rt;
432 const amrex::Real beta = dt;
434 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
437 for (
int lev = 0; lev < nlevels; ++lev) {
440 (*rho_times_detJ)(lev).setVal(0.0_rt);
441 amrex::MultiFab::AddProduct(
442 (*rho_times_detJ)(lev), density(lev), 0,
443 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
458 geom[lev], viscosity(lev));
464 lev, amrex::GetArrOfConstPtrs(b));
469 for (
int lev = 0; lev < nlevels; ++lev) {
470 auto& rhs = (*rhs_ptr)(lev);
472 const auto& rhs_arrs = rhs.arrays();
473 const auto& fld_arrs = field(lev).const_arrays();
474 const auto& rho_arrs = density(lev).const_arrays();
477 rhs, amrex::IntVect(0), ndim,
478 [=] AMREX_GPU_DEVICE(
479 int nbx,
int i,
int j,
int k,
int n)
noexcept {
480 rhs_arrs[nbx](i, j, k, n) =
481 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
484 amrex::Gpu::streamSynchronize();
486 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
488 auto rhs_ptr_comp = rhs_ptr->subview(i);
494 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
498 field.name() + std::to_string(i) +
"_solve", mlmg);
508 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
510 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
517template <
typename Scheme>
526 "DiffusionOp invoked for scalar PDE type");
531 PDEFields& fields,
const bool has_overset,
const bool mesh_mapping)
533 bool use_tensor_op =
true;
535 amrex::ParmParse pp(fields.
field.
name() +
"_diffusion");
536 pp.query(
"use_tensor_operator", use_tensor_op);
541 "Tensor and segregated operators should not be enabled "
547 fields, has_overset, mesh_mapping);
551 std::make_unique<ICNSDiffScalarSegregatedOp>(
552 fields, has_overset, mesh_mapping);
555 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:279
amrex::Vector< const amrex::iMultiFab * > vec_const_ptrs() const noexcept
Definition IntField.cpp:46
PDEFields & m_pdefields
Definition DiffusionOps.H:103
std::unique_ptr< amrex::MLTensorOp > m_applier
Definition DiffusionOps.H:111
virtual void setup_operator(amrex::MLTensorOp &linop, const amrex::Real alpha, const amrex::Real beta, const FieldState fstate)
Definition DiffusionOps.cpp:57
std::unique_ptr< amrex::MLTensorOp > m_solver
Definition DiffusionOps.H:110
DiffSolverIface(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition DiffusionOps.cpp:12
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:119
ICNSDiffScalarOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:55
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:188
bool m_mesh_mapping
Definition icns_diffusion.H:265
std::unique_ptr< amrex::MLABecLaplacian > m_solver_scalar
Definition icns_diffusion.H:267
Field & m_density
Definition icns_diffusion.H:263
MLMGOptions m_options
Definition icns_diffusion.H:264
PDEFields & m_pdefields
Definition icns_diffusion.H:262
std::unique_ptr< amrex::MLABecLaplacian > m_applier_scalar
Definition icns_diffusion.H:268
MLMGOptions m_options
Definition icns_diffusion.H:505
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_solver_scalar
Definition icns_diffusion.H:509
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:327
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:411
PDEFields & m_pdefields
Definition icns_diffusion.H:503
Field & m_density
Definition icns_diffusion.H:504
ICNSDiffScalarSegregatedOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:274
bool m_mesh_mapping
Definition icns_diffusion.H:506
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_applier_scalar
Definition icns_diffusion.H:511
ICNSDiffTensorOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:21
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:38
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:173
Definition AdvOp_Godunov.H:21
Definition console_io.cpp:30
void viscosity_to_uniform_space(amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > &b, const amr_wind::FieldRepo &repo, int lev)
Definition incflo_diffusion.cpp:207
amrex::Vector< amrex::Array< amrex::LinOpBCType, AMREX_SPACEDIM > > get_diffuse_tensor_bc(amr_wind::Field &velocity, amrex::Orientation::Side side) noexcept
Definition incflo_diffusion.cpp:10
amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > average_velocity_eta_to_faces(const amrex::Geometry &geom, amrex::MultiFab const &cc_eta)
Definition incflo_diffusion.cpp:107
Definition MLMGOptions.H:27
DiffusionOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:530
bool use_segregated_op
Definition icns_diffusion.H:528
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:573
std::unique_ptr< ICNSDiffScalarSegregatedOp > m_scalar_segregated_op
Definition icns_diffusion.H:522
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:560
std::unique_ptr< ICNSDiffScalarOp > m_scalar_op
Definition icns_diffusion.H:521
std::unique_ptr< ICNSDiffTensorOp > m_tensor_op
Definition icns_diffusion.H:520
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