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();
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_rt;
138 const amrex::Real beta = -1.0_rt;
141 for (
int lev = 0; lev < nlevels; ++lev) {
146 (*rho_times_detJ)(lev).setVal(0.0_rt);
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());
172 const auto& geom = repo.mesh().Geom();
174 const auto& density =
m_density.state(fstate);
175 const int nlevels = repo.num_active_levels();
176 const int ndim = field.num_comp();
177 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
179 Field const* mesh_detJ =
182 std::unique_ptr<ScratchField> rho_times_detJ =
187 const amrex::Real alpha = 1.0_rt;
188 const amrex::Real beta = dt;
191 for (
int lev = 0; lev < nlevels; ++lev) {
196 (*rho_times_detJ)(lev).setVal(0.0_rt);
197 amrex::MultiFab::AddProduct(
198 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
207 geom[lev], viscosity(lev));
215 for (
int lev = 0; lev < nlevels; ++lev) {
216 auto& rhs = (*rhs_ptr)(lev);
218 const auto& rhs_arrs = rhs.arrays();
219 const auto& fld_arrs = field(lev).const_arrays();
220 const auto& rho_arrs = density(lev).const_arrays();
223 rhs, amrex::IntVect(0), ndim,
224 [=] AMREX_GPU_DEVICE(
int nbx,
int i,
int j,
int k,
int n) {
225 rhs_arrs[nbx](i, j, k, n) =
226 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
229 amrex::Gpu::streamSynchronize();
234 m_pdefields.field.vec_ptrs(), rhs_ptr->vec_const_ptrs(),
255 const bool has_overset,
256 const bool mesh_mapping,
257 const std::string& prefix =
"diffusion")
259 ,
m_density(fields.repo.get_field(
"density"))
263 amrex::LPInfo isolve =
m_options.lpinfo();
264 amrex::LPInfo iapply;
266 iapply.setMaxCoarseningLevel(0);
274 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
277 mesh.Geom(0, mesh.finestLevel()),
278 mesh.boxArray(0, mesh.finestLevel()),
279 mesh.DistributionMap(0, mesh.finestLevel()), isolve);
281 mesh.Geom(0, mesh.finestLevel()),
282 mesh.boxArray(0, mesh.finestLevel()),
283 mesh.DistributionMap(0, mesh.finestLevel()), iapply);
288 mesh.Geom(0, mesh.finestLevel()),
289 mesh.boxArray(0, mesh.finestLevel()),
290 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve);
292 mesh.Geom(0, mesh.finestLevel()),
293 mesh.boxArray(0, mesh.finestLevel()),
294 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply);
305 template <
typename Scheme>
311 const int nlevels = repo.num_active_levels();
312 const auto& geom = repo.mesh().Geom();
314 auto& divtau =
m_pdefields.diff_term.state(tau_state);
315 const auto& density =
m_density.state(fstate);
317 Field const* mesh_detJ =
320 std::unique_ptr<ScratchField> rho_times_detJ =
325 const amrex::Real alpha = 0.0_rt;
326 const amrex::Real beta = -1.0_rt;
328 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
331 for (
int lev = 0; lev < nlevels; ++lev) {
337 (*rho_times_detJ)(lev).setVal(0.0_rt);
338 amrex::MultiFab::AddProduct(
339 (*rho_times_detJ)(lev), density(lev), 0,
340 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
346 lev, (*rho_times_detJ)(lev));
353 geom[lev], viscosity(lev));
360 lev, amrex::GetArrOfConstPtrs(b));
363 auto divtau_comp = divtau.subview(i);
367 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
376 const auto& density =
m_density.state(fstate);
377 const int nlevels = repo.num_active_levels();
378 const int ndim = field.num_comp();
379 auto rhs_ptr = repo.create_scratch_field(
"rhs", field.num_comp(), 0);
381 const auto& geom = repo.mesh().Geom();
383 Field const* mesh_detJ =
386 std::unique_ptr<ScratchField> rho_times_detJ =
391 const amrex::Real alpha = 1.0_rt;
392 const amrex::Real beta = dt;
394 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
397 for (
int lev = 0; lev < nlevels; ++lev) {
400 (*rho_times_detJ)(lev).setVal(0.0_rt);
401 amrex::MultiFab::AddProduct(
402 (*rho_times_detJ)(lev), density(lev), 0,
403 (*mesh_detJ)(lev), 0, 0, 1,
m_density.num_grow()[0]);
418 geom[lev], viscosity(lev));
424 lev, amrex::GetArrOfConstPtrs(b));
429 for (
int lev = 0; lev < nlevels; ++lev) {
430 auto& rhs = (*rhs_ptr)(lev);
432 const auto& rhs_arrs = rhs.arrays();
433 const auto& fld_arrs = field(lev).const_arrays();
434 const auto& rho_arrs = density(lev).const_arrays();
437 rhs, amrex::IntVect(0), ndim,
438 [=] AMREX_GPU_DEVICE(
int nbx,
int i,
int j,
int k,
int n) {
439 rhs_arrs[nbx](i, j, k, n) =
440 rho_arrs[nbx](i, j, k) * fld_arrs[nbx](i, j, k, n);
443 amrex::Gpu::streamSynchronize();
445 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
447 auto rhs_ptr_comp = rhs_ptr->subview(i);
453 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
457 field.name() + std::to_string(i) +
"_solve", mlmg);
467 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
469 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
476template <
typename Scheme>
485 "DiffusionOp invoked for scalar PDE type");
490 PDEFields& fields,
const bool has_overset,
const bool mesh_mapping)
492 bool use_tensor_op =
true;
494 amrex::ParmParse pp(fields.
field.
name() +
"_diffusion");
495 pp.query(
"use_tensor_operator", use_tensor_op);
500 "Tensor and segregated operators should not be enabled "
506 fields, has_overset, mesh_mapping);
510 std::make_unique<ICNSDiffScalarSegregatedOp>(
511 fields, has_overset, mesh_mapping);
514 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
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
std::unique_ptr< amrex::MLTensorOp > m_applier
Definition DiffusionOps.H:108
PDEFields & m_pdefields
Definition DiffusionOps.H:100
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
bool m_mesh_mapping
Definition icns_diffusion.H:244
MLMGOptions m_options
Definition icns_diffusion.H:243
std::unique_ptr< amrex::MLABecLaplacian > m_solver_scalar
Definition icns_diffusion.H:246
ICNSDiffScalarOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:54
PDEFields & m_pdefields
Definition icns_diffusion.H:241
std::unique_ptr< amrex::MLABecLaplacian > m_applier_scalar
Definition icns_diffusion.H:247
Field & m_density
Definition icns_diffusion.H:242
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:168
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:306
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_applier_scalar
Definition icns_diffusion.H:470
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:371
bool m_mesh_mapping
Definition icns_diffusion.H:465
ICNSDiffScalarSegregatedOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping, const std::string &prefix="diffusion")
Definition icns_diffusion.H:253
PDEFields & m_pdefields
Definition icns_diffusion.H:462
Field & m_density
Definition icns_diffusion.H:463
amrex::Array< std::unique_ptr< amrex::MLABecLaplacian >, AMREX_SPACEDIM > m_solver_scalar
Definition icns_diffusion.H:468
MLMGOptions m_options
Definition icns_diffusion.H:464
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
@ CELL
Cell-centered (default)
Definition FieldDescTypes.H:30
@ New
Same as FieldState::NP1.
Definition FieldDescTypes.H:22
Definition console_io.cpp:30
amrex::Vector< amrex::Array< amrex::LinOpBCType, AMREX_SPACEDIM > > get_diffuse_tensor_bc(kynema_sgf::Field &velocity, amrex::Orientation::Side side)
Definition incflo_diffusion.cpp:11
void viscosity_to_uniform_space(amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > &b, const kynema_sgf::FieldRepo &repo, int lev)
Definition incflo_diffusion.cpp:208
amrex::Array< amrex::MultiFab, AMREX_SPACEDIM > average_velocity_eta_to_faces(const amrex::Geometry &geom, amrex::MultiFab const &cc_eta)
Definition incflo_diffusion.cpp:108
void print_mlmg_info(const std::string &solve_name, const amrex::MLMG &mlmg)
Definition console_io.cpp:170
Definition AdvOp_Godunov.H:22
Definition MLMGOptions.H:27
std::unique_ptr< ICNSDiffTensorOp > m_tensor_op
Definition icns_diffusion.H:479
void compute_diff_term(const FieldState fstate)
Definition icns_diffusion.H:519
DiffusionOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:489
std::unique_ptr< ICNSDiffScalarSegregatedOp > m_scalar_segregated_op
Definition icns_diffusion.H:481
void linsys_solve(const amrex::Real dt)
Definition icns_diffusion.H:532
std::unique_ptr< ICNSDiffScalarOp > m_scalar_op
Definition icns_diffusion.H:480
bool use_segregated_op
Definition icns_diffusion.H:487
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