/home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/icns/icns_diffusion.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/icns/icns_diffusion.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
icns_diffusion.H
Go to the documentation of this file.
1#ifndef ICNS_DIFFUSION_H
2#define ICNS_DIFFUSION_H
3
4#include <memory>
5
12#include "AMReX_REAL.H"
13
14using namespace amrex::literals;
15
16namespace amr_wind::pde {
17
18class ICNSDiffTensorOp : public DiffSolverIface<amrex::MLTensorOp>
19{
20public:
22 PDEFields& fields, const bool has_overset, const bool mesh_mapping)
23 : DiffSolverIface<amrex::MLTensorOp>(fields, has_overset, mesh_mapping)
24 {
25 this->m_solver->setDomainBC(
27 this->m_pdefields.field, amrex::Orientation::low),
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));
35 }
36
37 template <typename Scheme>
38 void compute_diff_term(const FieldState fstate)
39 {
40 this->setup_operator(*this->m_applier, 0.0_rt, -1.0_rt, fstate);
41
42 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
44 : fstate;
45 auto& divtau = this->m_pdefields.diff_term.state(tau_state);
46
47 amrex::MLMG mlmg(*this->m_applier);
48 mlmg.apply(divtau.vec_ptrs(), this->m_pdefields.field.vec_ptrs());
49 }
50};
51
53{
54public:
56 PDEFields& fields,
57 const bool has_overset,
58 const bool mesh_mapping,
59 const std::string& prefix = "diffusion")
60 : m_pdefields(fields)
61 , m_density(fields.repo.get_field("density"))
62 , m_options(prefix, m_pdefields.field.name() + "_" + prefix)
63 , m_mesh_mapping(mesh_mapping)
64 {
65 amrex::LPInfo isolve = m_options.lpinfo();
66 amrex::LPInfo iapply;
67
68 iapply.setMaxCoarseningLevel(0);
69 const auto& mesh = m_pdefields.repo.mesh();
70
71 const auto& bclo = diffusion::get_diffuse_tensor_bc(
72 m_pdefields.field, amrex::Orientation::low);
73 const auto& bchi = diffusion::get_diffuse_tensor_bc(
74 m_pdefields.field, amrex::Orientation::high);
75
76 if (!has_overset) {
77
78 m_solver_scalar = std::make_unique<amrex::MLABecLaplacian>(
79 mesh.Geom(0, mesh.finestLevel()),
80 mesh.boxArray(0, mesh.finestLevel()),
81 mesh.DistributionMap(0, mesh.finestLevel()), isolve,
82 amrex::Vector<
83 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
84 AMREX_SPACEDIM);
85 m_applier_scalar = std::make_unique<amrex::MLABecLaplacian>(
86 mesh.Geom(0, mesh.finestLevel()),
87 mesh.boxArray(0, mesh.finestLevel()),
88 mesh.DistributionMap(0, mesh.finestLevel()), iapply,
89 amrex::Vector<
90 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
91 AMREX_SPACEDIM);
92 } else {
93 auto imask =
94 fields.repo.get_int_field("mask_cell").vec_const_ptrs();
95 m_solver_scalar = std::make_unique<amrex::MLABecLaplacian>(
96 mesh.Geom(0, mesh.finestLevel()),
97 mesh.boxArray(0, mesh.finestLevel()),
98 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve,
99 amrex::Vector<
100 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
101 AMREX_SPACEDIM);
102 m_applier_scalar = std::make_unique<amrex::MLABecLaplacian>(
103 mesh.Geom(0, mesh.finestLevel()),
104 mesh.boxArray(0, mesh.finestLevel()),
105 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply,
106 amrex::Vector<
107 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
108 AMREX_SPACEDIM);
109 }
110
111 m_solver_scalar->setMaxOrder(m_options.max_order);
112 m_applier_scalar->setMaxOrder(m_options.max_order);
113
114 m_solver_scalar->setDomainBC(bclo, bchi);
115 m_applier_scalar->setDomainBC(bclo, bchi);
116 }
117
118 template <typename Scheme>
119 void compute_diff_term(const FieldState fstate)
120 {
121 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
123 : fstate;
124 const auto& repo = m_pdefields.repo;
125 const int nlevels = repo.num_active_levels();
126 const auto& geom = repo.mesh().Geom();
127
128 bool diff_for_RHS(fstate == amr_wind::FieldState::New);
129 auto& divtau = m_pdefields.diff_term.state(tau_state);
130 const auto& density = m_density.state(fstate);
131 const auto& viscosity = m_pdefields.mueff;
132 Field const* mesh_detJ =
133 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
134 : nullptr;
135 std::unique_ptr<ScratchField> rho_times_detJ =
136 m_mesh_mapping ? repo.create_scratch_field(
137 1, m_density.num_grow()[0], FieldLoc::CELL)
138 : nullptr;
139
140 const amrex::Real alpha = 0.0_rt;
141 const amrex::Real beta = -1.0_rt;
142 m_applier_scalar->setScalars(alpha, beta);
143
144 for (int lev = 0; lev < nlevels; ++lev) {
145 m_applier_scalar->setLevelBC(lev, &m_pdefields.field(lev));
146
147 // A coeffs
148 if (m_mesh_mapping) {
149 (*rho_times_detJ)(lev).setVal(0.0_rt);
150 amrex::MultiFab::AddProduct(
151 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
152 0, 0, 1, m_density.num_grow()[0]);
153 m_applier_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
154 } else {
155 m_applier_scalar->setACoeffs(lev, density(lev));
156 }
157
158 // B coeffs
160 geom[lev], viscosity(lev));
161 if (m_mesh_mapping) {
163 }
164 m_applier_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
165 }
166
167 amrex::MLMG mlmg(*m_applier_scalar);
168 mlmg.apply(divtau.vec_ptrs(), m_pdefields.field.vec_ptrs());
169
170 if (!diff_for_RHS) {
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();
174
175 amrex::ParallelFor(
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;
182 });
183 }
184 amrex::Gpu::streamSynchronize();
185 }
186 }
187
188 void linsys_solve(const amrex::Real dt)
189 {
190 const FieldState fstate = FieldState::New;
191 auto& repo = m_pdefields.repo;
192 const auto& geom = repo.mesh().Geom();
193 const auto& field = m_pdefields.field;
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);
198 const auto& viscosity = m_pdefields.mueff;
199 Field const* mesh_detJ =
200 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
201 : nullptr;
202 std::unique_ptr<ScratchField> rho_times_detJ =
203 m_mesh_mapping ? repo.create_scratch_field(
204 1, m_density.num_grow()[0], FieldLoc::CELL)
205 : nullptr;
206
207 const amrex::Real alpha = 1.0_rt;
208 const amrex::Real beta = dt;
209 m_solver_scalar->setScalars(alpha, beta);
210
211 for (int lev = 0; lev < nlevels; ++lev) {
212 m_solver_scalar->setLevelBC(lev, &m_pdefields.field(lev));
213
214 // A coeffs
215 if (m_mesh_mapping) {
216 (*rho_times_detJ)(lev).setVal(0.0_rt);
217 amrex::MultiFab::AddProduct(
218 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
219 0, 0, 1, m_density.num_grow()[0]);
220 m_solver_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
221 } else {
222 m_solver_scalar->setACoeffs(lev, density(lev));
223 }
224
225 // B coeffs
227 geom[lev], viscosity(lev));
228 if (m_mesh_mapping) {
230 }
231 m_solver_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
232 }
233
234 // Always multiply with rho since there is no diffusion term for density
235 for (int lev = 0; lev < nlevels; ++lev) {
236 auto& rhs = (*rhs_ptr)(lev);
237
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();
241
242 amrex::ParallelFor(
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);
248 });
249 }
250 amrex::Gpu::streamSynchronize();
251
252 amrex::MLMG mlmg(*m_solver_scalar);
253 m_options(mlmg);
254 mlmg.solve(
255 m_pdefields.field.vec_ptrs(), rhs_ptr->vec_const_ptrs(),
256 m_options.rel_tol, m_options.abs_tol);
257
258 io::print_mlmg_info(field.name() + "_multicomponent_solve", mlmg);
259 }
260
261protected:
265 bool m_mesh_mapping{false};
266
267 std::unique_ptr<amrex::MLABecLaplacian> m_solver_scalar;
268 std::unique_ptr<amrex::MLABecLaplacian> m_applier_scalar;
269};
270
272{
273public:
275 PDEFields& fields,
276 const bool has_overset,
277 const bool mesh_mapping,
278 const std::string& prefix = "diffusion")
279 : m_pdefields(fields)
280 , m_density(fields.repo.get_field("density"))
281 , m_options(prefix, m_pdefields.field.name() + "_" + prefix)
282 , m_mesh_mapping(mesh_mapping)
283 {
284 amrex::LPInfo isolve = m_options.lpinfo();
285 amrex::LPInfo iapply;
286
287 iapply.setMaxCoarseningLevel(0);
288 const auto& mesh = m_pdefields.repo.mesh();
289
290 const auto& bclo = diffusion::get_diffuse_tensor_bc(
291 m_pdefields.field, amrex::Orientation::low);
292 const auto& bchi = diffusion::get_diffuse_tensor_bc(
293 m_pdefields.field, amrex::Orientation::high);
294
295 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
296 if (!has_overset) {
297 m_solver_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
298 mesh.Geom(0, mesh.finestLevel()),
299 mesh.boxArray(0, mesh.finestLevel()),
300 mesh.DistributionMap(0, mesh.finestLevel()), isolve);
301 m_applier_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
302 mesh.Geom(0, mesh.finestLevel()),
303 mesh.boxArray(0, mesh.finestLevel()),
304 mesh.DistributionMap(0, mesh.finestLevel()), iapply);
305 } else {
306 auto imask =
307 fields.repo.get_int_field("mask_cell").vec_const_ptrs();
308 m_solver_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
309 mesh.Geom(0, mesh.finestLevel()),
310 mesh.boxArray(0, mesh.finestLevel()),
311 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve);
312 m_applier_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
313 mesh.Geom(0, mesh.finestLevel()),
314 mesh.boxArray(0, mesh.finestLevel()),
315 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply);
316 }
317
318 m_solver_scalar[i]->setMaxOrder(m_options.max_order);
319 m_applier_scalar[i]->setMaxOrder(m_options.max_order);
320
321 m_solver_scalar[i]->setDomainBC(bclo[i], bchi[i]);
322 m_applier_scalar[i]->setDomainBC(bclo[i], bchi[i]);
323 }
324 }
325
326 template <typename Scheme>
327 void compute_diff_term(const FieldState fstate)
328 {
329 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
331 : fstate;
332 const auto& repo = m_pdefields.repo;
333 const int nlevels = repo.num_active_levels();
334 const auto& geom = repo.mesh().Geom();
335
336 auto& divtau = m_pdefields.diff_term.state(tau_state);
337 bool diff_for_RHS(fstate == amr_wind::FieldState::New);
338 const auto& density = m_density.state(fstate);
339 const auto& viscosity = m_pdefields.mueff;
340 Field const* mesh_detJ =
341 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
342 : nullptr;
343 std::unique_ptr<ScratchField> rho_times_detJ =
344 m_mesh_mapping ? repo.create_scratch_field(
345 1, m_density.num_grow()[0], FieldLoc::CELL)
346 : nullptr;
347
348 const amrex::Real alpha = 0.0_rt;
349 const amrex::Real beta = -1.0_rt;
350
351 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
352 m_applier_scalar[i]->setScalars(alpha, beta);
353
354 for (int lev = 0; lev < nlevels; ++lev) {
355
356 m_applier_scalar[i]->setLevelBC(
357 lev, &m_pdefields.field.subview(i)(lev));
358
359 if (m_mesh_mapping) {
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]);
364 }
365
366 // A coeffs
367 if (m_mesh_mapping) {
368 m_applier_scalar[i]->setACoeffs(
369 lev, (*rho_times_detJ)(lev));
370 } else {
371 m_applier_scalar[i]->setACoeffs(lev, density(lev));
372 }
373
374 // B coeffs
376 geom[lev], viscosity(lev));
377
378 if (m_mesh_mapping) {
380 }
381
382 m_applier_scalar[i]->setBCoeffs(
383 lev, amrex::GetArrOfConstPtrs(b));
384 }
385
386 auto divtau_comp = divtau.subview(i);
387 auto vel_comp = m_pdefields.field.subview(i);
388
389 amrex::MLMG mlmg(*m_applier_scalar[i]);
390 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
391 }
392
393 if (!diff_for_RHS) {
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();
397
398 amrex::ParallelFor(
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;
405 });
406 }
407 amrex::Gpu::streamSynchronize();
408 }
409 }
410
411 void linsys_solve(const amrex::Real dt)
412 {
413 const FieldState fstate = FieldState::New;
414 auto& repo = m_pdefields.repo;
415 const auto& field = m_pdefields.field;
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);
420 const auto& viscosity = m_pdefields.mueff;
421 const auto& geom = repo.mesh().Geom();
422
423 Field const* mesh_detJ =
424 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
425 : nullptr;
426 std::unique_ptr<ScratchField> rho_times_detJ =
427 m_mesh_mapping ? repo.create_scratch_field(
428 1, m_density.num_grow()[0], FieldLoc::CELL)
429 : nullptr;
430
431 const amrex::Real alpha = 1.0_rt;
432 const amrex::Real beta = dt;
433
434 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
435 m_solver_scalar[i]->setScalars(alpha, beta);
436
437 for (int lev = 0; lev < nlevels; ++lev) {
438
439 if (m_mesh_mapping) {
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]);
444 }
445
446 m_solver_scalar[i]->setLevelBC(
447 lev, &m_pdefields.field.subview(i)(lev));
448
449 // A coeffs
450 if (m_mesh_mapping) {
451 m_solver_scalar[i]->setACoeffs(lev, (*rho_times_detJ)(lev));
452 } else {
453 m_solver_scalar[i]->setACoeffs(lev, density(lev));
454 }
455
456 // B coeffs
458 geom[lev], viscosity(lev));
459 if (m_mesh_mapping) {
461 }
462
463 m_solver_scalar[i]->setBCoeffs(
464 lev, amrex::GetArrOfConstPtrs(b));
465 }
466 }
467
468 // Always multiply with rho since there is no diffusion term for density
469 for (int lev = 0; lev < nlevels; ++lev) {
470 auto& rhs = (*rhs_ptr)(lev);
471
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();
475
476 amrex::ParallelFor(
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);
482 });
483 }
484 amrex::Gpu::streamSynchronize();
485
486 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
487 auto vel_comp = m_pdefields.field.subview(i);
488 auto rhs_ptr_comp = rhs_ptr->subview(i);
489
490 amrex::MLMG mlmg(*m_solver_scalar[i]);
491 m_options(mlmg);
492
493 mlmg.solve(
494 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
495 m_options.rel_tol, m_options.abs_tol);
496
498 field.name() + std::to_string(i) + "_solve", mlmg);
499 }
500 }
501
502protected:
506 bool m_mesh_mapping{false};
507
508 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
510 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
512};
513
517template <typename Scheme>
518struct DiffusionOp<ICNS, Scheme>
519{
520 std::unique_ptr<ICNSDiffTensorOp> m_tensor_op;
521 std::unique_ptr<ICNSDiffScalarOp> m_scalar_op;
522 std::unique_ptr<ICNSDiffScalarSegregatedOp> m_scalar_segregated_op;
523
524 static_assert(
525 ICNS::ndim == AMREX_SPACEDIM,
526 "DiffusionOp invoked for scalar PDE type");
527
528 bool use_segregated_op = false;
529
531 PDEFields& fields, const bool has_overset, const bool mesh_mapping)
532 {
533 bool use_tensor_op = true;
534
535 amrex::ParmParse pp(fields.field.name() + "_diffusion");
536 pp.query("use_tensor_operator", use_tensor_op);
537 pp.query("use_segregated_operator", use_segregated_op);
538
539 if (use_tensor_op && use_segregated_op) {
540 amrex::Abort(
541 "Tensor and segregated operators should not be enabled "
542 "simultaneously.");
543 }
544
545 if (use_tensor_op) {
546 m_tensor_op = std::make_unique<ICNSDiffTensorOp>(
547 fields, has_overset, mesh_mapping);
548 } else {
549 if (use_segregated_op) {
551 std::make_unique<ICNSDiffScalarSegregatedOp>(
552 fields, has_overset, mesh_mapping);
553 } else {
554 m_scalar_op = std::make_unique<ICNSDiffScalarOp>(
555 fields, has_overset, mesh_mapping);
556 }
557 }
558 }
559
560 void compute_diff_term(const FieldState fstate)
561 {
562 if (m_tensor_op) {
563 m_tensor_op->compute_diff_term<Scheme>(fstate);
564 } else {
565 if (use_segregated_op) {
566 m_scalar_segregated_op->compute_diff_term<Scheme>(fstate);
567 } else {
568 m_scalar_op->compute_diff_term<Scheme>(fstate);
569 }
570 }
571 }
572
573 void linsys_solve(const amrex::Real dt)
574 {
575 if (m_tensor_op) {
576 m_tensor_op->linsys_solve(dt);
577 } else {
578 if (use_segregated_op) {
579 m_scalar_segregated_op->linsys_solve(dt);
580 } else {
581 m_scalar_op->linsys_solve(dt);
582 }
583 }
584 }
585};
586
587} // namespace amr_wind::pde
588
589#endif /* ICNS_DIFFUSION_H */
Definition Field.H:116
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
Definition icns.H:34
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