/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
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
13namespace amr_wind::pde {
14
15class ICNSDiffTensorOp : public DiffSolverIface<amrex::MLTensorOp>
16{
17public:
19 PDEFields& fields, const bool has_overset, const bool mesh_mapping)
20 : DiffSolverIface<amrex::MLTensorOp>(fields, has_overset, mesh_mapping)
21 {
22 this->m_solver->setDomainBC(
24 this->m_pdefields.field, amrex::Orientation::low),
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));
32 }
33
34 template <typename Scheme>
35 void compute_diff_term(const FieldState fstate)
36 {
37 this->setup_operator(*this->m_applier, 0.0, -1.0, fstate);
38
39 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
41 : fstate;
42 auto& divtau = this->m_pdefields.diff_term.state(tau_state);
43
44 amrex::MLMG mlmg(*this->m_applier);
45 mlmg.apply(divtau.vec_ptrs(), this->m_pdefields.field.vec_ptrs());
46 }
47};
48
50{
51public:
53 PDEFields& fields,
54 const bool has_overset,
55 const bool mesh_mapping,
56 const std::string& prefix = "diffusion")
57 : m_pdefields(fields)
58 , m_density(fields.repo.get_field("density"))
59 , m_options(prefix, m_pdefields.field.name() + "_" + prefix)
60 , m_mesh_mapping(mesh_mapping)
61 {
62 amrex::LPInfo isolve = m_options.lpinfo();
63 amrex::LPInfo iapply;
64
65 iapply.setMaxCoarseningLevel(0);
66 const auto& mesh = m_pdefields.repo.mesh();
67
68 const auto& bclo = diffusion::get_diffuse_tensor_bc(
69 m_pdefields.field, amrex::Orientation::low);
70 const auto& bchi = diffusion::get_diffuse_tensor_bc(
71 m_pdefields.field, amrex::Orientation::high);
72
73 if (!has_overset) {
74
75 m_solver_scalar = std::make_unique<amrex::MLABecLaplacian>(
76 mesh.Geom(0, mesh.finestLevel()),
77 mesh.boxArray(0, mesh.finestLevel()),
78 mesh.DistributionMap(0, mesh.finestLevel()), isolve,
79 amrex::Vector<
80 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
81 AMREX_SPACEDIM);
82 m_applier_scalar = std::make_unique<amrex::MLABecLaplacian>(
83 mesh.Geom(0, mesh.finestLevel()),
84 mesh.boxArray(0, mesh.finestLevel()),
85 mesh.DistributionMap(0, mesh.finestLevel()), iapply,
86 amrex::Vector<
87 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
88 AMREX_SPACEDIM);
89 } else {
90 auto imask =
91 fields.repo.get_int_field("mask_cell").vec_const_ptrs();
92 m_solver_scalar = std::make_unique<amrex::MLABecLaplacian>(
93 mesh.Geom(0, mesh.finestLevel()),
94 mesh.boxArray(0, mesh.finestLevel()),
95 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve,
96 amrex::Vector<
97 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
98 AMREX_SPACEDIM);
99 m_applier_scalar = std::make_unique<amrex::MLABecLaplacian>(
100 mesh.Geom(0, mesh.finestLevel()),
101 mesh.boxArray(0, mesh.finestLevel()),
102 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply,
103 amrex::Vector<
104 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
105 AMREX_SPACEDIM);
106 }
107
108 m_solver_scalar->setMaxOrder(m_options.max_order);
110
111 m_solver_scalar->setDomainBC(bclo, bchi);
112 m_applier_scalar->setDomainBC(bclo, bchi);
113 }
114
115 template <typename Scheme>
116 void compute_diff_term(const FieldState fstate)
117 {
118 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
120 : fstate;
121 const auto& repo = m_pdefields.repo;
122 const int nlevels = repo.num_active_levels();
123 const auto& geom = repo.mesh().Geom();
124
125 bool diff_for_RHS(fstate == amr_wind::FieldState::New);
126 auto& divtau = m_pdefields.diff_term.state(tau_state);
127 const auto& density = m_density.state(fstate);
128 const auto& viscosity = m_pdefields.mueff;
129 Field const* mesh_detJ =
130 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
131 : nullptr;
132 std::unique_ptr<ScratchField> rho_times_detJ =
133 m_mesh_mapping ? repo.create_scratch_field(
135 : nullptr;
136
137 const amrex::Real alpha = 0.0;
138 const amrex::Real beta = -1.0;
139 m_applier_scalar->setScalars(alpha, beta);
140
141 for (int lev = 0; lev < nlevels; ++lev) {
142 m_applier_scalar->setLevelBC(lev, &m_pdefields.field(lev));
143
144 // A coeffs
145 if (m_mesh_mapping) {
146 (*rho_times_detJ)(lev).setVal(0.0);
147 amrex::MultiFab::AddProduct(
148 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
149 0, 0, 1, m_density.num_grow()[0]);
150 m_applier_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
151 } else {
152 m_applier_scalar->setACoeffs(lev, density(lev));
153 }
154
155 // B coeffs
157 geom[lev], viscosity(lev));
158 if (m_mesh_mapping) {
160 }
161 m_applier_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
162 }
163
164 amrex::MLMG mlmg(*m_applier_scalar);
165 mlmg.apply(divtau.vec_ptrs(), m_pdefields.field.vec_ptrs());
166
167 if (!diff_for_RHS) {
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();
171
172 amrex::ParallelFor(
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;
179 });
180 }
181 amrex::Gpu::streamSynchronize();
182 }
183 }
184
185 void linsys_solve(const amrex::Real dt)
186 {
187 const FieldState fstate = FieldState::New;
188 auto& repo = m_pdefields.repo;
189 const auto& geom = repo.mesh().Geom();
190 const auto& field = m_pdefields.field;
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);
195 const auto& viscosity = m_pdefields.mueff;
196 Field const* mesh_detJ =
197 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
198 : nullptr;
199 std::unique_ptr<ScratchField> rho_times_detJ =
200 m_mesh_mapping ? repo.create_scratch_field(
202 : nullptr;
203
204 const amrex::Real alpha = 1.0;
205 const amrex::Real beta = dt;
206 m_solver_scalar->setScalars(alpha, beta);
207
208 for (int lev = 0; lev < nlevels; ++lev) {
209 m_solver_scalar->setLevelBC(lev, &m_pdefields.field(lev));
210
211 // A coeffs
212 if (m_mesh_mapping) {
213 (*rho_times_detJ)(lev).setVal(0.0);
214 amrex::MultiFab::AddProduct(
215 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
216 0, 0, 1, m_density.num_grow()[0]);
217 m_solver_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
218 } else {
219 m_solver_scalar->setACoeffs(lev, density(lev));
220 }
221
222 // B coeffs
224 geom[lev], viscosity(lev));
225 if (m_mesh_mapping) {
227 }
228 m_solver_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
229 }
230
231 // Always multiply with rho since there is no diffusion term for density
232 for (int lev = 0; lev < nlevels; ++lev) {
233 auto& rhs = (*rhs_ptr)(lev);
234
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();
238
239 amrex::ParallelFor(
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);
245 });
246 }
247 amrex::Gpu::streamSynchronize();
248
249 amrex::MLMG mlmg(*m_solver_scalar);
250 m_options(mlmg);
251 mlmg.solve(
252 m_pdefields.field.vec_ptrs(), rhs_ptr->vec_const_ptrs(),
254
255 io::print_mlmg_info(field.name() + "_multicomponent_solve", mlmg);
256 }
257
258protected:
262 bool m_mesh_mapping{false};
263
264 std::unique_ptr<amrex::MLABecLaplacian> m_solver_scalar;
265 std::unique_ptr<amrex::MLABecLaplacian> m_applier_scalar;
266};
267
269{
270public:
272 PDEFields& fields,
273 const bool has_overset,
274 const bool mesh_mapping,
275 const std::string& prefix = "diffusion")
276 : m_pdefields(fields)
277 , m_density(fields.repo.get_field("density"))
278 , m_options(prefix, m_pdefields.field.name() + "_" + prefix)
279 , m_mesh_mapping(mesh_mapping)
280 {
281 amrex::LPInfo isolve = m_options.lpinfo();
282 amrex::LPInfo iapply;
283
284 iapply.setMaxCoarseningLevel(0);
285 const auto& mesh = m_pdefields.repo.mesh();
286
287 const auto& bclo = diffusion::get_diffuse_tensor_bc(
288 m_pdefields.field, amrex::Orientation::low);
289 const auto& bchi = diffusion::get_diffuse_tensor_bc(
290 m_pdefields.field, amrex::Orientation::high);
291
292 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
293 if (!has_overset) {
294 m_solver_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
295 mesh.Geom(0, mesh.finestLevel()),
296 mesh.boxArray(0, mesh.finestLevel()),
297 mesh.DistributionMap(0, mesh.finestLevel()), isolve);
298 m_applier_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
299 mesh.Geom(0, mesh.finestLevel()),
300 mesh.boxArray(0, mesh.finestLevel()),
301 mesh.DistributionMap(0, mesh.finestLevel()), iapply);
302 } else {
303 auto imask =
304 fields.repo.get_int_field("mask_cell").vec_const_ptrs();
305 m_solver_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
306 mesh.Geom(0, mesh.finestLevel()),
307 mesh.boxArray(0, mesh.finestLevel()),
308 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve);
309 m_applier_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
310 mesh.Geom(0, mesh.finestLevel()),
311 mesh.boxArray(0, mesh.finestLevel()),
312 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply);
313 }
314
315 m_solver_scalar[i]->setMaxOrder(m_options.max_order);
316 m_applier_scalar[i]->setMaxOrder(m_options.max_order);
317
318 m_solver_scalar[i]->setDomainBC(bclo[i], bchi[i]);
319 m_applier_scalar[i]->setDomainBC(bclo[i], bchi[i]);
320 }
321 }
322
323 template <typename Scheme>
324 void compute_diff_term(const FieldState fstate)
325 {
326 auto tau_state = std::is_same<Scheme, fvm::Godunov>::value
328 : fstate;
329 const auto& repo = m_pdefields.repo;
330 const int nlevels = repo.num_active_levels();
331 const auto& geom = repo.mesh().Geom();
332
333 auto& divtau = m_pdefields.diff_term.state(tau_state);
334 bool diff_for_RHS(fstate == amr_wind::FieldState::New);
335 const auto& density = m_density.state(fstate);
336 const auto& viscosity = m_pdefields.mueff;
337 Field const* mesh_detJ =
338 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
339 : nullptr;
340 std::unique_ptr<ScratchField> rho_times_detJ =
341 m_mesh_mapping ? repo.create_scratch_field(
343 : nullptr;
344
345 const amrex::Real alpha = 0.0;
346 const amrex::Real beta = -1.0;
347
348 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
349 m_applier_scalar[i]->setScalars(alpha, beta);
350
351 for (int lev = 0; lev < nlevels; ++lev) {
352
353 m_applier_scalar[i]->setLevelBC(
354 lev, &m_pdefields.field.subview(i)(lev));
355
356 if (m_mesh_mapping) {
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]);
361 }
362
363 // A coeffs
364 if (m_mesh_mapping) {
365 m_applier_scalar[i]->setACoeffs(
366 lev, (*rho_times_detJ)(lev));
367 } else {
368 m_applier_scalar[i]->setACoeffs(lev, density(lev));
369 }
370
371 // B coeffs
373 geom[lev], viscosity(lev));
374
375 if (m_mesh_mapping) {
377 }
378
379 m_applier_scalar[i]->setBCoeffs(
380 lev, amrex::GetArrOfConstPtrs(b));
381 }
382
383 auto divtau_comp = divtau.subview(i);
384 auto vel_comp = m_pdefields.field.subview(i);
385
386 amrex::MLMG mlmg(*m_applier_scalar[i]);
387 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
388 }
389
390 if (!diff_for_RHS) {
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();
394
395 amrex::ParallelFor(
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;
402 });
403 }
404 amrex::Gpu::streamSynchronize();
405 }
406 }
407
408 void linsys_solve(const amrex::Real dt)
409 {
410 const FieldState fstate = FieldState::New;
411 auto& repo = m_pdefields.repo;
412 const auto& field = m_pdefields.field;
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);
417 const auto& viscosity = m_pdefields.mueff;
418 const auto& geom = repo.mesh().Geom();
419
420 Field const* mesh_detJ =
421 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
422 : nullptr;
423 std::unique_ptr<ScratchField> rho_times_detJ =
424 m_mesh_mapping ? repo.create_scratch_field(
426 : nullptr;
427
428 const amrex::Real alpha = 1.0;
429 const amrex::Real beta = dt;
430
431 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
432 m_solver_scalar[i]->setScalars(alpha, beta);
433
434 for (int lev = 0; lev < nlevels; ++lev) {
435
436 if (m_mesh_mapping) {
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]);
441 }
442
443 m_solver_scalar[i]->setLevelBC(
444 lev, &m_pdefields.field.subview(i)(lev));
445
446 // A coeffs
447 if (m_mesh_mapping) {
448 m_solver_scalar[i]->setACoeffs(lev, (*rho_times_detJ)(lev));
449 } else {
450 m_solver_scalar[i]->setACoeffs(lev, density(lev));
451 }
452
453 // B coeffs
455 geom[lev], viscosity(lev));
456 if (m_mesh_mapping) {
458 }
459
460 m_solver_scalar[i]->setBCoeffs(
461 lev, amrex::GetArrOfConstPtrs(b));
462 }
463 }
464
465 // Always multiply with rho since there is no diffusion term for density
466 for (int lev = 0; lev < nlevels; ++lev) {
467 auto& rhs = (*rhs_ptr)(lev);
468
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();
472
473 amrex::ParallelFor(
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);
479 });
480 }
481 amrex::Gpu::streamSynchronize();
482
483 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
484 auto vel_comp = m_pdefields.field.subview(i);
485 auto rhs_ptr_comp = rhs_ptr->subview(i);
486
487 amrex::MLMG mlmg(*m_solver_scalar[i]);
488 m_options(mlmg);
489
490 mlmg.solve(
491 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
493
495 field.name() + std::to_string(i) + "_solve", mlmg);
496 }
497 }
498
499protected:
503 bool m_mesh_mapping{false};
504
505 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
507 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
509};
510
514template <typename Scheme>
515struct DiffusionOp<ICNS, Scheme>
516{
517 std::unique_ptr<ICNSDiffTensorOp> m_tensor_op;
518 std::unique_ptr<ICNSDiffScalarOp> m_scalar_op;
519 std::unique_ptr<ICNSDiffScalarSegregatedOp> m_scalar_segregated_op;
520
521 static_assert(
522 ICNS::ndim == AMREX_SPACEDIM,
523 "DiffusionOp invoked for scalar PDE type");
524
525 bool use_segregated_op = false;
526
528 PDEFields& fields, const bool has_overset, const bool mesh_mapping)
529 {
530 bool use_tensor_op = true;
531
532 amrex::ParmParse pp(fields.field.name() + "_diffusion");
533 pp.query("use_tensor_operator", use_tensor_op);
534 pp.query("use_segregated_operator", use_segregated_op);
535
536 if (use_tensor_op && use_segregated_op) {
537 amrex::Abort(
538 "Tensor and segregated operators should not be enabled "
539 "simultaneously.");
540 }
541
542 if (use_tensor_op) {
543 m_tensor_op = std::make_unique<ICNSDiffTensorOp>(
544 fields, has_overset, mesh_mapping);
545 } else {
546 if (use_segregated_op) {
547 m_scalar_segregated_op =
548 std::make_unique<ICNSDiffScalarSegregatedOp>(
549 fields, has_overset, mesh_mapping);
550 } else {
551 m_scalar_op = std::make_unique<ICNSDiffScalarOp>(
552 fields, has_overset, mesh_mapping);
553 }
554 }
555 }
556
557 void compute_diff_term(const FieldState fstate)
558 {
559 if (m_tensor_op) {
560 m_tensor_op->compute_diff_term<Scheme>(fstate);
561 } else {
562 if (use_segregated_op) {
563 m_scalar_segregated_op->compute_diff_term<Scheme>(fstate);
564 } else {
565 m_scalar_op->compute_diff_term<Scheme>(fstate);
566 }
567 }
568 }
569
570 void linsys_solve(const amrex::Real dt)
571 {
572 if (m_tensor_op) {
573 m_tensor_op->linsys_solve(dt);
574 } else {
575 if (use_segregated_op) {
576 m_scalar_segregated_op->linsys_solve(dt);
577 } else {
578 m_scalar_op->linsys_solve(dt);
579 }
580 }
581 }
582};
583
584} // namespace amr_wind::pde
585
586#endif /* ICNS_DIFFUSION_H */
Definition Field.H:116
ViewField< Field > subview(const int scomp=0, const int ncomp=1)
Definition Field.H:296
const amrex::IntVect & num_grow() const
Ghost cells.
Definition Field.H:137
const std::string & name() const
Name of this field (including state information)
Definition Field.H:125
int num_comp() const
Number of components for this field.
Definition Field.H:134
amrex::Vector< amrex::MultiFab * > vec_ptrs() noexcept
Return a vector of MultiFab pointers for all levels.
Definition Field.cpp:142
Field & state(const FieldState fstate)
Return field at a different time state.
Definition Field.cpp:114
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
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
const amrex::AmrCore & mesh() const
Return a reference to the underlying AMR mesh instance.
Definition FieldRepo.H:358
amrex::Vector< const amrex::iMultiFab * > vec_const_ptrs() const noexcept
Definition IntField.cpp:46
Definition DiffusionOps.H:30
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
Definition icns_diffusion.H:50
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
Definition icns_diffusion.H:269
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
Definition icns_diffusion.H:16
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)
@ New
Same as FieldState::NP1.
void print_mlmg_info(const std::string &solve_name, const amrex::MLMG &mlmg)
Definition console_io.cpp:168
Definition AdvOp_Godunov.H:16
Definition console_io.cpp:25
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
amrex::LPInfo & lpinfo()
Linear operator options during construction.
Definition MLMGOptions.H:48
int max_order
Definition MLMGOptions.H:51
amrex::Real abs_tol
Absolute tolerance for convergence checks.
Definition MLMGOptions.H:57
amrex::Real rel_tol
Relative tolerance for convergence of MLMG solvers.
Definition MLMGOptions.H:54
DiffusionOp(PDEFields &fields, const bool has_overset, const bool mesh_mapping)
Definition icns_diffusion.H:527
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
Definition PDEOps.H:167
Definition icns.H:34
static constexpr int ndim
Definition icns.H:40
Definition PDEFields.H:27
Field & diff_term
Diffusion term for this PDE.
Definition PDEFields.H:41
Field & mueff
Effective visocity field (e.g., velocity_mueff)
Definition PDEFields.H:36
FieldRepo & repo
Reference to the field repository instance.
Definition PDEFields.H:31
Field & field
Solution variable (e.g., velocity, temperature)
Definition PDEFields.H:34