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>
44 auto& divtau = this->
m_pdefields.diff_term.state(tau_state);
47 mlmg.apply(divtau.vec_ptrs(), this->m_pdefields.field.vec_ptrs());
56 const bool has_overset,
57 const bool mesh_mapping,
58 const std::string& prefix =
"diffusion")
60 ,
m_density(fields.repo.get_field(
"density"))
64 amrex::LPInfo isolve =
m_options.lpinfo();
67 iapply.setMaxCoarseningLevel(0);
78 mesh.Geom(0, mesh.finestLevel()),
79 mesh.boxArray(0, mesh.finestLevel()),
80 mesh.DistributionMap(0, mesh.finestLevel()), isolve,
82 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
85 mesh.Geom(0, mesh.finestLevel()),
86 mesh.boxArray(0, mesh.finestLevel()),
87 mesh.DistributionMap(0, mesh.finestLevel()), iapply,
89 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
95 mesh.Geom(0, mesh.finestLevel()),
96 mesh.boxArray(0, mesh.finestLevel()),
97 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve,
99 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
102 mesh.Geom(0, mesh.finestLevel()),
103 mesh.boxArray(0, mesh.finestLevel()),
104 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply,
106 amrex::FabFactory<amrex::MLABecLaplacian::FAB>
const*>(),
117 template <
typename Scheme>
123 const int nlevels = repo.num_active_levels();
124 const auto& geom = repo.mesh().Geom();
127 auto& divtau =
m_pdefields.diff_term.state(tau_state);
128 const auto& density =
m_density.state(fstate);
130 Field const* mesh_detJ =
133 std::unique_ptr<ScratchField> rho_times_detJ =
138 const amrex::Real alpha = 0.0_rt;
139 const amrex::Real beta = -1.0_rt;
142 for (
int lev = 0; lev < nlevels; ++lev) {
147 (*rho_times_detJ)(lev).setVal(0.0_rt);
148 amrex::MultiFab::AddProduct(
149 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
158 geom[lev], viscosity(lev));
166 mlmg.apply(divtau.vec_ptrs(),
m_pdefields.field.vec_ptrs());
169 for (
int lev = 0; lev < nlevels; ++lev) {
170 const auto& divtau_arrs = divtau(lev).arrays();
171 const auto& rho_arrs = density(lev).const_arrays();
175 [=] AMREX_GPU_DEVICE(
int nbx,
int i,
int j,
int k) {
176 amrex::Real rhoinv = 1.0_rt / rho_arrs[nbx](i, j, k);
177 divtau_arrs[nbx](i, j, k, 0) *= rhoinv;
178 divtau_arrs[nbx](i, j, k, 1) *= rhoinv;
179 divtau_arrs[nbx](i, j, k, 2) *= rhoinv;
182 amrex::Gpu::streamSynchronize();
190 const auto& geom = repo.mesh().Geom();
192 const auto& density =
m_density.state(fstate);
193 const int nlevels = repo.num_active_levels();
194 const int ndim = field.num_comp();
195 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
197 Field const* mesh_detJ =
200 std::unique_ptr<ScratchField> rho_times_detJ =
205 const amrex::Real alpha = 1.0_rt;
206 const amrex::Real beta = dt;
209 for (
int lev = 0; lev < nlevels; ++lev) {
214 (*rho_times_detJ)(lev).setVal(0.0_rt);
215 amrex::MultiFab::AddProduct(
216 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
225 geom[lev], viscosity(lev));
233 for (
int lev = 0; lev < nlevels; ++lev) {
234 auto& rhs = (*rhs_ptr)(lev);
236 const auto& rhs_arrs = rhs.arrays();
237 const auto& fld_arrs = field(lev).const_arrays();
238 const auto& rho_arrs = density(lev).const_arrays();
241 rhs, amrex::IntVect(0), ndim,
242 [=] AMREX_GPU_DEVICE(
int nbx,
int i,
int j,
int k,
int n) {
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>
329 const int nlevels = repo.num_active_levels();
330 const auto& geom = repo.mesh().Geom();
332 auto& divtau =
m_pdefields.diff_term.state(tau_state);
334 const auto& density =
m_density.state(fstate);
336 Field const* mesh_detJ =
339 std::unique_ptr<ScratchField> rho_times_detJ =
344 const amrex::Real alpha = 0.0_rt;
345 const amrex::Real beta = -1.0_rt;
347 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
350 for (
int lev = 0; lev < nlevels; ++lev) {
356 (*rho_times_detJ)(lev).setVal(0.0_rt);
357 amrex::MultiFab::AddProduct(
358 (*rho_times_detJ)(lev), density(lev), 0,
359 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
365 lev, (*rho_times_detJ)(lev));
372 geom[lev], viscosity(lev));
379 lev, amrex::GetArrOfConstPtrs(b));
382 auto divtau_comp = divtau.subview(i);
386 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
390 for (
int lev = 0; lev < nlevels; ++lev) {
391 const auto& divtau_arrs = divtau(lev).arrays();
392 const auto& rho_arrs = density(lev).const_arrays();
396 [=] AMREX_GPU_DEVICE(
int nbx,
int i,
int j,
int k) {
397 amrex::Real rhoinv = 1.0_rt / rho_arrs[nbx](i, j, k);
398 divtau_arrs[nbx](i, j, k, 0) *= rhoinv;
399 divtau_arrs[nbx](i, j, k, 1) *= rhoinv;
400 divtau_arrs[nbx](i, j, k, 2) *= rhoinv;
403 amrex::Gpu::streamSynchronize();
412 const auto& density =
m_density.state(fstate);
413 const int nlevels = repo.num_active_levels();
414 const int ndim = field.num_comp();
415 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
417 const auto& geom = repo.mesh().Geom();
419 Field const* mesh_detJ =
422 std::unique_ptr<ScratchField> rho_times_detJ =
427 const amrex::Real alpha = 1.0_rt;
428 const amrex::Real beta = dt;
430 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
433 for (
int lev = 0; lev < nlevels; ++lev) {
436 (*rho_times_detJ)(lev).setVal(0.0_rt);
437 amrex::MultiFab::AddProduct(
438 (*rho_times_detJ)(lev), density(lev), 0,
439 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
454 geom[lev], viscosity(lev));
460 lev, amrex::GetArrOfConstPtrs(b));
465 for (
int lev = 0; lev < nlevels; ++lev) {
466 auto& rhs = (*rhs_ptr)(lev);
468 const auto& rhs_arrs = rhs.arrays();
469 const auto& fld_arrs = field(lev).const_arrays();
470 const auto& rho_arrs = density(lev).const_arrays();
473 rhs, amrex::IntVect(0), ndim,
474 [=] AMREX_GPU_DEVICE(
int nbx,
int i,
int j,
int k,
int n) {
475 rhs_arrs[nbx](i, j, k, n) =
476 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
479 amrex::Gpu::streamSynchronize();
481 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
483 auto rhs_ptr_comp = rhs_ptr->subview(i);
489 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
493 field.name() + std::to_string(i) +
"_solve", mlmg);
503 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
505 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
512template <
typename Scheme>
521 "DiffusionOp invoked for scalar PDE type");
526 PDEFields& fields,
const bool has_overset,
const bool mesh_mapping)
528 bool use_tensor_op =
true;
530 amrex::ParmParse pp(fields.
field.
name() +
"_diffusion");
531 pp.query(
"use_tensor_operator", use_tensor_op);
536 "Tensor and segregated operators should not be enabled "
542 fields, has_overset, mesh_mapping);
546 std::make_unique<ICNSDiffScalarSegregatedOp>(
547 fields, has_overset, mesh_mapping);
550 fields, has_overset, mesh_mapping);
const std::string & name() const
Name of this field (including state information)
Definition Field.H:121
IntField & get_int_field(const std::string &name, 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
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, amrex::Real alpha, amrex::Real beta, FieldState fstate)
Definition DiffusionOps.cpp:57
std::unique_ptr< amrex::MLTensorOp > m_solver
Definition DiffusionOps.H:107
DiffSolverIface(PDEFields &fields, bool has_overset, bool mesh_mapping, const std::string &prefix="diffusion")
Definition DiffusionOps.cpp:12
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:118
ICNSDiffScalarOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:54
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:186
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:500
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_solver_scalar
Definition icns_diffusion.H:504
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:324
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:407
PDEFields & m_pdefields
Definition icns_diffusion.H:498
Field & m_density
Definition icns_diffusion.H:499
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:501
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_applier_scalar
Definition icns_diffusion.H:506
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:16
@ New
Same as FieldState::NP1.
Definition FieldDescTypes.H:22
@ CELL
Cell-centered (default)
Definition FieldDescTypes.H:30
void print_mlmg_info(const std::string &solve_name, const amrex::MLMG &mlmg)
Definition console_io.cpp:170
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:206
amrex::Vector< amrex::Array< amrex::LinOpBCType, AMREX_SPACEDIM > > get_diffuse_tensor_bc(amr_wind::Field &velocity, amrex::Orientation::Side side)
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:106
Definition MLMGOptions.H:27
DiffusionOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:525
bool use_segregated_op
Definition icns_diffusion.H:523
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:568
std::unique_ptr< ICNSDiffScalarSegregatedOp > m_scalar_segregated_op
Definition icns_diffusion.H:517
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:555
std::unique_ptr< ICNSDiffScalarOp > m_scalar_op
Definition icns_diffusion.H:516
std::unique_ptr< ICNSDiffTensorOp > m_tensor_op
Definition icns_diffusion.H:515
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