/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 =
43 std::is_same_v<Scheme, fvm::Godunov> ? FieldState::New : fstate;
44 auto& divtau = this->m_pdefields.diff_term.state(tau_state);
45
46 amrex::MLMG mlmg(*this->m_applier);
47 mlmg.apply(divtau.vec_ptrs(), this->m_pdefields.field.vec_ptrs());
48 }
49};
50
52{
53public:
55 PDEFields& fields,
56 const bool has_overset,
57 const bool mesh_mapping,
58 const std::string& prefix = "diffusion")
59 : m_pdefields(fields)
60 , m_density(fields.repo.get_field("density"))
61 , m_options(prefix, m_pdefields.field.name() + "_" + prefix)
62 , m_mesh_mapping(mesh_mapping)
63 {
64 amrex::LPInfo isolve = m_options.lpinfo();
65 amrex::LPInfo iapply;
66
67 iapply.setMaxCoarseningLevel(0);
68 const auto& mesh = m_pdefields.repo.mesh();
69
70 const auto& bclo = diffusion::get_diffuse_tensor_bc(
71 m_pdefields.field, amrex::Orientation::low);
72 const auto& bchi = diffusion::get_diffuse_tensor_bc(
73 m_pdefields.field, amrex::Orientation::high);
74
75 if (!has_overset) {
76
77 m_solver_scalar = std::make_unique<amrex::MLABecLaplacian>(
78 mesh.Geom(0, mesh.finestLevel()),
79 mesh.boxArray(0, mesh.finestLevel()),
80 mesh.DistributionMap(0, mesh.finestLevel()), isolve,
81 amrex::Vector<
82 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
83 AMREX_SPACEDIM);
84 m_applier_scalar = std::make_unique<amrex::MLABecLaplacian>(
85 mesh.Geom(0, mesh.finestLevel()),
86 mesh.boxArray(0, mesh.finestLevel()),
87 mesh.DistributionMap(0, mesh.finestLevel()), iapply,
88 amrex::Vector<
89 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
90 AMREX_SPACEDIM);
91 } else {
92 auto imask =
93 fields.repo.get_int_field("mask_cell").vec_const_ptrs();
94 m_solver_scalar = std::make_unique<amrex::MLABecLaplacian>(
95 mesh.Geom(0, mesh.finestLevel()),
96 mesh.boxArray(0, mesh.finestLevel()),
97 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve,
98 amrex::Vector<
99 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
100 AMREX_SPACEDIM);
101 m_applier_scalar = std::make_unique<amrex::MLABecLaplacian>(
102 mesh.Geom(0, mesh.finestLevel()),
103 mesh.boxArray(0, mesh.finestLevel()),
104 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply,
105 amrex::Vector<
106 amrex::FabFactory<amrex::MLABecLaplacian::FAB> const*>(),
107 AMREX_SPACEDIM);
108 }
109
110 m_solver_scalar->setMaxOrder(m_options.max_order);
111 m_applier_scalar->setMaxOrder(m_options.max_order);
112
113 m_solver_scalar->setDomainBC(bclo, bchi);
114 m_applier_scalar->setDomainBC(bclo, bchi);
115 }
116
117 template <typename Scheme>
118 void compute_diff_term(const FieldState fstate)
119 {
120 auto tau_state =
121 std::is_same_v<Scheme, fvm::Godunov> ? FieldState::New : fstate;
122 const auto& repo = m_pdefields.repo;
123 const int nlevels = repo.num_active_levels();
124 const auto& geom = repo.mesh().Geom();
125
126 bool diff_for_RHS(fstate == amr_wind::FieldState::New);
127 auto& divtau = m_pdefields.diff_term.state(tau_state);
128 const auto& density = m_density.state(fstate);
129 const auto& viscosity = m_pdefields.mueff;
130 Field const* mesh_detJ =
131 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
132 : nullptr;
133 std::unique_ptr<ScratchField> rho_times_detJ =
134 m_mesh_mapping ? repo.create_scratch_field(
135 1, m_density.num_grow()[0], FieldLoc::CELL)
136 : nullptr;
137
138 const amrex::Real alpha = 0.0_rt;
139 const amrex::Real beta = -1.0_rt;
140 m_applier_scalar->setScalars(alpha, beta);
141
142 for (int lev = 0; lev < nlevels; ++lev) {
143 m_applier_scalar->setLevelBC(lev, &m_pdefields.field(lev));
144
145 // A coeffs
146 if (m_mesh_mapping) {
147 (*rho_times_detJ)(lev).setVal(0.0_rt);
148 amrex::MultiFab::AddProduct(
149 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
150 0, 0, 1, m_density.num_grow()[0]);
151 m_applier_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
152 } else {
153 m_applier_scalar->setACoeffs(lev, density(lev));
154 }
155
156 // B coeffs
158 geom[lev], viscosity(lev));
159 if (m_mesh_mapping) {
161 }
162 m_applier_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
163 }
164
165 amrex::MLMG mlmg(*m_applier_scalar);
166 mlmg.apply(divtau.vec_ptrs(), m_pdefields.field.vec_ptrs());
167
168 if (!diff_for_RHS) {
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();
172
173 amrex::ParallelFor(
174 divtau(lev),
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;
180 });
181 }
182 amrex::Gpu::streamSynchronize();
183 }
184 }
185
186 void linsys_solve(const amrex::Real dt)
187 {
188 const FieldState fstate = FieldState::New;
189 auto& repo = m_pdefields.repo;
190 const auto& geom = repo.mesh().Geom();
191 const auto& field = m_pdefields.field;
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);
196 const auto& viscosity = m_pdefields.mueff;
197 Field const* mesh_detJ =
198 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
199 : nullptr;
200 std::unique_ptr<ScratchField> rho_times_detJ =
201 m_mesh_mapping ? repo.create_scratch_field(
202 1, m_density.num_grow()[0], FieldLoc::CELL)
203 : nullptr;
204
205 const amrex::Real alpha = 1.0_rt;
206 const amrex::Real beta = dt;
207 m_solver_scalar->setScalars(alpha, beta);
208
209 for (int lev = 0; lev < nlevels; ++lev) {
210 m_solver_scalar->setLevelBC(lev, &m_pdefields.field(lev));
211
212 // A coeffs
213 if (m_mesh_mapping) {
214 (*rho_times_detJ)(lev).setVal(0.0_rt);
215 amrex::MultiFab::AddProduct(
216 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
217 0, 0, 1, m_density.num_grow()[0]);
218 m_solver_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
219 } else {
220 m_solver_scalar->setACoeffs(lev, density(lev));
221 }
222
223 // B coeffs
225 geom[lev], viscosity(lev));
226 if (m_mesh_mapping) {
228 }
229 m_solver_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
230 }
231
232 // Always multiply with rho since there is no diffusion term for density
233 for (int lev = 0; lev < nlevels; ++lev) {
234 auto& rhs = (*rhs_ptr)(lev);
235
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();
239
240 amrex::ParallelFor(
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);
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(),
253 m_options.rel_tol, m_options.abs_tol);
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 =
327 std::is_same_v<Scheme, fvm::Godunov> ? FieldState::New : fstate;
328 const auto& repo = m_pdefields.repo;
329 const int nlevels = repo.num_active_levels();
330 const auto& geom = repo.mesh().Geom();
331
332 auto& divtau = m_pdefields.diff_term.state(tau_state);
333 bool diff_for_RHS(fstate == amr_wind::FieldState::New);
334 const auto& density = m_density.state(fstate);
335 const auto& viscosity = m_pdefields.mueff;
336 Field const* mesh_detJ =
337 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
338 : nullptr;
339 std::unique_ptr<ScratchField> rho_times_detJ =
340 m_mesh_mapping ? repo.create_scratch_field(
341 1, m_density.num_grow()[0], FieldLoc::CELL)
342 : nullptr;
343
344 const amrex::Real alpha = 0.0_rt;
345 const amrex::Real beta = -1.0_rt;
346
347 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
348 m_applier_scalar[i]->setScalars(alpha, beta);
349
350 for (int lev = 0; lev < nlevels; ++lev) {
351
352 m_applier_scalar[i]->setLevelBC(
353 lev, &m_pdefields.field.subview(i)(lev));
354
355 if (m_mesh_mapping) {
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]);
360 }
361
362 // A coeffs
363 if (m_mesh_mapping) {
364 m_applier_scalar[i]->setACoeffs(
365 lev, (*rho_times_detJ)(lev));
366 } else {
367 m_applier_scalar[i]->setACoeffs(lev, density(lev));
368 }
369
370 // B coeffs
372 geom[lev], viscosity(lev));
373
374 if (m_mesh_mapping) {
376 }
377
378 m_applier_scalar[i]->setBCoeffs(
379 lev, amrex::GetArrOfConstPtrs(b));
380 }
381
382 auto divtau_comp = divtau.subview(i);
383 auto vel_comp = m_pdefields.field.subview(i);
384
385 amrex::MLMG mlmg(*m_applier_scalar[i]);
386 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
387 }
388
389 if (!diff_for_RHS) {
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();
393
394 amrex::ParallelFor(
395 divtau(lev),
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;
401 });
402 }
403 amrex::Gpu::streamSynchronize();
404 }
405 }
406
407 void linsys_solve(const amrex::Real dt)
408 {
409 const FieldState fstate = FieldState::New;
410 auto& repo = m_pdefields.repo;
411 const auto& field = m_pdefields.field;
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);
416 const auto& viscosity = m_pdefields.mueff;
417 const auto& geom = repo.mesh().Geom();
418
419 Field const* mesh_detJ =
420 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
421 : nullptr;
422 std::unique_ptr<ScratchField> rho_times_detJ =
423 m_mesh_mapping ? repo.create_scratch_field(
424 1, m_density.num_grow()[0], FieldLoc::CELL)
425 : nullptr;
426
427 const amrex::Real alpha = 1.0_rt;
428 const amrex::Real beta = dt;
429
430 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
431 m_solver_scalar[i]->setScalars(alpha, beta);
432
433 for (int lev = 0; lev < nlevels; ++lev) {
434
435 if (m_mesh_mapping) {
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]);
440 }
441
442 m_solver_scalar[i]->setLevelBC(
443 lev, &m_pdefields.field.subview(i)(lev));
444
445 // A coeffs
446 if (m_mesh_mapping) {
447 m_solver_scalar[i]->setACoeffs(lev, (*rho_times_detJ)(lev));
448 } else {
449 m_solver_scalar[i]->setACoeffs(lev, density(lev));
450 }
451
452 // B coeffs
454 geom[lev], viscosity(lev));
455 if (m_mesh_mapping) {
457 }
458
459 m_solver_scalar[i]->setBCoeffs(
460 lev, amrex::GetArrOfConstPtrs(b));
461 }
462 }
463
464 // Always multiply with rho since there is no diffusion term for density
465 for (int lev = 0; lev < nlevels; ++lev) {
466 auto& rhs = (*rhs_ptr)(lev);
467
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();
471
472 amrex::ParallelFor(
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);
477 });
478 }
479 amrex::Gpu::streamSynchronize();
480
481 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
482 auto vel_comp = m_pdefields.field.subview(i);
483 auto rhs_ptr_comp = rhs_ptr->subview(i);
484
485 amrex::MLMG mlmg(*m_solver_scalar[i]);
486 m_options(mlmg);
487
488 mlmg.solve(
489 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
490 m_options.rel_tol, m_options.abs_tol);
491
493 field.name() + std::to_string(i) + "_solve", mlmg);
494 }
495 }
496
497protected:
501 bool m_mesh_mapping{false};
502
503 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
505 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
507};
508
512template <typename Scheme>
513struct DiffusionOp<ICNS, Scheme>
514{
515 std::unique_ptr<ICNSDiffTensorOp> m_tensor_op;
516 std::unique_ptr<ICNSDiffScalarOp> m_scalar_op;
517 std::unique_ptr<ICNSDiffScalarSegregatedOp> m_scalar_segregated_op;
518
519 static_assert(
520 ICNS::ndim == AMREX_SPACEDIM,
521 "DiffusionOp invoked for scalar PDE type");
522
523 bool use_segregated_op = false;
524
526 PDEFields& fields, const bool has_overset, const bool mesh_mapping)
527 {
528 bool use_tensor_op = true;
529
530 amrex::ParmParse pp(fields.field.name() + "_diffusion");
531 pp.query("use_tensor_operator", use_tensor_op);
532 pp.query("use_segregated_operator", use_segregated_op);
533
534 if (use_tensor_op && use_segregated_op) {
535 amrex::Abort(
536 "Tensor and segregated operators should not be enabled "
537 "simultaneously.");
538 }
539
540 if (use_tensor_op) {
541 m_tensor_op = std::make_unique<ICNSDiffTensorOp>(
542 fields, has_overset, mesh_mapping);
543 } else {
544 if (use_segregated_op) {
546 std::make_unique<ICNSDiffScalarSegregatedOp>(
547 fields, has_overset, mesh_mapping);
548 } else {
549 m_scalar_op = std::make_unique<ICNSDiffScalarOp>(
550 fields, has_overset, mesh_mapping);
551 }
552 }
553 }
554
555 void compute_diff_term(const FieldState fstate)
556 {
557 if (m_tensor_op) {
558 m_tensor_op->compute_diff_term<Scheme>(fstate);
559 } else {
560 if (use_segregated_op) {
561 m_scalar_segregated_op->compute_diff_term<Scheme>(fstate);
562 } else {
563 m_scalar_op->compute_diff_term<Scheme>(fstate);
564 }
565 }
566 }
567
568 void linsys_solve(const amrex::Real dt)
569 {
570 if (m_tensor_op) {
571 m_tensor_op->linsys_solve(dt);
572 } else {
573 if (use_segregated_op) {
574 m_scalar_segregated_op->linsys_solve(dt);
575 } else {
576 m_scalar_op->linsys_solve(dt);
577 }
578 }
579 }
580};
581
582} // namespace amr_wind::pde
583
584#endif /* ICNS_DIFFUSION_H */
Definition Field.H:112
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
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