/home/runner/work/kynema-sgf/kynema-sgf/src/equation_systems/icns/icns_diffusion.H Source File

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/equation_systems/icns/icns_diffusion.H Source File
Kynema-SGF 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 kynema_sgf::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 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(
134 1, m_density.num_grow()[0], FieldLoc::CELL)
135 : nullptr;
136
137 const amrex::Real alpha = 0.0_rt;
138 const amrex::Real beta = -1.0_rt;
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_rt);
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
168 void linsys_solve(const amrex::Real dt)
169 {
170 const FieldState fstate = FieldState::New;
171 auto& repo = m_pdefields.repo;
172 const auto& geom = repo.mesh().Geom();
173 const auto& field = m_pdefields.field;
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);
178 const auto& viscosity = m_pdefields.mueff;
179 Field const* mesh_detJ =
180 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
181 : nullptr;
182 std::unique_ptr<ScratchField> rho_times_detJ =
183 m_mesh_mapping ? repo.create_scratch_field(
184 1, m_density.num_grow()[0], FieldLoc::CELL)
185 : nullptr;
186
187 const amrex::Real alpha = 1.0_rt;
188 const amrex::Real beta = dt;
189 m_solver_scalar->setScalars(alpha, beta);
190
191 for (int lev = 0; lev < nlevels; ++lev) {
192 m_solver_scalar->setLevelBC(lev, &m_pdefields.field(lev));
193
194 // A coeffs
195 if (m_mesh_mapping) {
196 (*rho_times_detJ)(lev).setVal(0.0_rt);
197 amrex::MultiFab::AddProduct(
198 (*rho_times_detJ)(lev), density(lev), 0, (*mesh_detJ)(lev),
199 0, 0, 1, m_density.num_grow()[0]);
200 m_solver_scalar->setACoeffs(lev, (*rho_times_detJ)(lev));
201 } else {
202 m_solver_scalar->setACoeffs(lev, density(lev));
203 }
204
205 // B coeffs
207 geom[lev], viscosity(lev));
208 if (m_mesh_mapping) {
210 }
211 m_solver_scalar->setBCoeffs(lev, amrex::GetArrOfConstPtrs(b));
212 }
213
214 // Always multiply with rho since there is no diffusion term for density
215 for (int lev = 0; lev < nlevels; ++lev) {
216 auto& rhs = (*rhs_ptr)(lev);
217
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();
221
222 amrex::ParallelFor(
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);
227 });
228 }
229 amrex::Gpu::streamSynchronize();
230
231 amrex::MLMG mlmg(*m_solver_scalar);
232 m_options(mlmg);
233 mlmg.solve(
234 m_pdefields.field.vec_ptrs(), rhs_ptr->vec_const_ptrs(),
235 m_options.rel_tol, m_options.abs_tol);
236
237 io::print_mlmg_info(field.name() + "_multicomponent_solve", mlmg);
238 }
239
240protected:
244 bool m_mesh_mapping{false};
245
246 std::unique_ptr<amrex::MLABecLaplacian> m_solver_scalar;
247 std::unique_ptr<amrex::MLABecLaplacian> m_applier_scalar;
248};
249
251{
252public:
254 PDEFields& fields,
255 const bool has_overset,
256 const bool mesh_mapping,
257 const std::string& prefix = "diffusion")
258 : m_pdefields(fields)
259 , m_density(fields.repo.get_field("density"))
260 , m_options(prefix, m_pdefields.field.name() + "_" + prefix)
261 , m_mesh_mapping(mesh_mapping)
262 {
263 amrex::LPInfo isolve = m_options.lpinfo();
264 amrex::LPInfo iapply;
265
266 iapply.setMaxCoarseningLevel(0);
267 const auto& mesh = m_pdefields.repo.mesh();
268
269 const auto& bclo = diffusion::get_diffuse_tensor_bc(
270 m_pdefields.field, amrex::Orientation::low);
271 const auto& bchi = diffusion::get_diffuse_tensor_bc(
272 m_pdefields.field, amrex::Orientation::high);
273
274 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
275 if (!has_overset) {
276 m_solver_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
277 mesh.Geom(0, mesh.finestLevel()),
278 mesh.boxArray(0, mesh.finestLevel()),
279 mesh.DistributionMap(0, mesh.finestLevel()), isolve);
280 m_applier_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
281 mesh.Geom(0, mesh.finestLevel()),
282 mesh.boxArray(0, mesh.finestLevel()),
283 mesh.DistributionMap(0, mesh.finestLevel()), iapply);
284 } else {
285 auto imask =
286 fields.repo.get_int_field("mask_cell").vec_const_ptrs();
287 m_solver_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
288 mesh.Geom(0, mesh.finestLevel()),
289 mesh.boxArray(0, mesh.finestLevel()),
290 mesh.DistributionMap(0, mesh.finestLevel()), imask, isolve);
291 m_applier_scalar[i] = std::make_unique<amrex::MLABecLaplacian>(
292 mesh.Geom(0, mesh.finestLevel()),
293 mesh.boxArray(0, mesh.finestLevel()),
294 mesh.DistributionMap(0, mesh.finestLevel()), imask, iapply);
295 }
296
297 m_solver_scalar[i]->setMaxOrder(m_options.max_order);
298 m_applier_scalar[i]->setMaxOrder(m_options.max_order);
299
300 m_solver_scalar[i]->setDomainBC(bclo[i], bchi[i]);
301 m_applier_scalar[i]->setDomainBC(bclo[i], bchi[i]);
302 }
303 }
304
305 template <typename Scheme>
306 void compute_diff_term(const FieldState fstate)
307 {
308 auto tau_state =
309 std::is_same_v<Scheme, fvm::Godunov> ? FieldState::New : fstate;
310 const auto& repo = m_pdefields.repo;
311 const int nlevels = repo.num_active_levels();
312 const auto& geom = repo.mesh().Geom();
313
314 auto& divtau = m_pdefields.diff_term.state(tau_state);
315 const auto& density = m_density.state(fstate);
316 const auto& viscosity = m_pdefields.mueff;
317 Field const* mesh_detJ =
318 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
319 : nullptr;
320 std::unique_ptr<ScratchField> rho_times_detJ =
321 m_mesh_mapping ? repo.create_scratch_field(
322 1, m_density.num_grow()[0], FieldLoc::CELL)
323 : nullptr;
324
325 const amrex::Real alpha = 0.0_rt;
326 const amrex::Real beta = -1.0_rt;
327
328 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
329 m_applier_scalar[i]->setScalars(alpha, beta);
330
331 for (int lev = 0; lev < nlevels; ++lev) {
332
333 m_applier_scalar[i]->setLevelBC(
334 lev, &m_pdefields.field.subview(i)(lev));
335
336 if (m_mesh_mapping) {
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]);
341 }
342
343 // A coeffs
344 if (m_mesh_mapping) {
345 m_applier_scalar[i]->setACoeffs(
346 lev, (*rho_times_detJ)(lev));
347 } else {
348 m_applier_scalar[i]->setACoeffs(lev, density(lev));
349 }
350
351 // B coeffs
353 geom[lev], viscosity(lev));
354
355 if (m_mesh_mapping) {
357 }
358
359 m_applier_scalar[i]->setBCoeffs(
360 lev, amrex::GetArrOfConstPtrs(b));
361 }
362
363 auto divtau_comp = divtau.subview(i);
364 auto vel_comp = m_pdefields.field.subview(i);
365
366 amrex::MLMG mlmg(*m_applier_scalar[i]);
367 mlmg.apply(divtau_comp.vec_ptrs(), vel_comp.vec_ptrs());
368 }
369 }
370
371 void linsys_solve(const amrex::Real dt)
372 {
373 const FieldState fstate = FieldState::New;
374 auto& repo = m_pdefields.repo;
375 const auto& field = m_pdefields.field;
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);
380 const auto& viscosity = m_pdefields.mueff;
381 const auto& geom = repo.mesh().Geom();
382
383 Field const* mesh_detJ =
384 m_mesh_mapping ? &(repo.get_mesh_mapping_det_j(FieldLoc::CELL))
385 : nullptr;
386 std::unique_ptr<ScratchField> rho_times_detJ =
387 m_mesh_mapping ? repo.create_scratch_field(
388 1, m_density.num_grow()[0], FieldLoc::CELL)
389 : nullptr;
390
391 const amrex::Real alpha = 1.0_rt;
392 const amrex::Real beta = dt;
393
394 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
395 m_solver_scalar[i]->setScalars(alpha, beta);
396
397 for (int lev = 0; lev < nlevels; ++lev) {
398
399 if (m_mesh_mapping) {
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]);
404 }
405
406 m_solver_scalar[i]->setLevelBC(
407 lev, &m_pdefields.field.subview(i)(lev));
408
409 // A coeffs
410 if (m_mesh_mapping) {
411 m_solver_scalar[i]->setACoeffs(lev, (*rho_times_detJ)(lev));
412 } else {
413 m_solver_scalar[i]->setACoeffs(lev, density(lev));
414 }
415
416 // B coeffs
418 geom[lev], viscosity(lev));
419 if (m_mesh_mapping) {
421 }
422
423 m_solver_scalar[i]->setBCoeffs(
424 lev, amrex::GetArrOfConstPtrs(b));
425 }
426 }
427
428 // Always multiply with rho since there is no diffusion term for density
429 for (int lev = 0; lev < nlevels; ++lev) {
430 auto& rhs = (*rhs_ptr)(lev);
431
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();
435
436 amrex::ParallelFor(
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);
441 });
442 }
443 amrex::Gpu::streamSynchronize();
444
445 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
446 auto vel_comp = m_pdefields.field.subview(i);
447 auto rhs_ptr_comp = rhs_ptr->subview(i);
448
449 amrex::MLMG mlmg(*m_solver_scalar[i]);
450 m_options(mlmg);
451
452 mlmg.solve(
453 vel_comp.vec_ptrs(), rhs_ptr_comp.vec_const_ptrs(),
454 m_options.rel_tol, m_options.abs_tol);
455
457 field.name() + std::to_string(i) + "_solve", mlmg);
458 }
459 }
460
461protected:
465 bool m_mesh_mapping{false};
466
467 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
469 amrex::Array<std::unique_ptr<amrex::MLABecLaplacian>, AMREX_SPACEDIM>
471};
472
476template <typename Scheme>
477struct DiffusionOp<ICNS, Scheme>
478{
479 std::unique_ptr<ICNSDiffTensorOp> m_tensor_op;
480 std::unique_ptr<ICNSDiffScalarOp> m_scalar_op;
481 std::unique_ptr<ICNSDiffScalarSegregatedOp> m_scalar_segregated_op;
482
483 static_assert(
484 ICNS::ndim == AMREX_SPACEDIM,
485 "DiffusionOp invoked for scalar PDE type");
486
487 bool use_segregated_op = false;
488
490 PDEFields& fields, const bool has_overset, const bool mesh_mapping)
491 {
492 bool use_tensor_op = true;
493
494 amrex::ParmParse pp(fields.field.name() + "_diffusion");
495 pp.query("use_tensor_operator", use_tensor_op);
496 pp.query("use_segregated_operator", use_segregated_op);
497
498 if (use_tensor_op && use_segregated_op) {
499 amrex::Abort(
500 "Tensor and segregated operators should not be enabled "
501 "simultaneously.");
502 }
503
504 if (use_tensor_op) {
505 m_tensor_op = std::make_unique<ICNSDiffTensorOp>(
506 fields, has_overset, mesh_mapping);
507 } else {
508 if (use_segregated_op) {
510 std::make_unique<ICNSDiffScalarSegregatedOp>(
511 fields, has_overset, mesh_mapping);
512 } else {
513 m_scalar_op = std::make_unique<ICNSDiffScalarOp>(
514 fields, has_overset, mesh_mapping);
515 }
516 }
517 }
518
519 void compute_diff_term(const FieldState fstate)
520 {
521 if (m_tensor_op) {
522 m_tensor_op->compute_diff_term<Scheme>(fstate);
523 } else {
524 if (use_segregated_op) {
525 m_scalar_segregated_op->compute_diff_term<Scheme>(fstate);
526 } else {
527 m_scalar_op->compute_diff_term<Scheme>(fstate);
528 }
529 }
530 }
531
532 void linsys_solve(const amrex::Real dt)
533 {
534 if (m_tensor_op) {
535 m_tensor_op->linsys_solve(dt);
536 } else {
537 if (use_segregated_op) {
538 m_scalar_segregated_op->linsys_solve(dt);
539 } else {
540 m_scalar_op->linsys_solve(dt);
541 }
542 }
543 }
544};
545
546} // namespace kynema_sgf::pde
547
548#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
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
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