/home/runner/work/amr-wind/amr-wind/amr-wind/core/FieldFillPatchOps.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/core/FieldFillPatchOps.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
FieldFillPatchOps.H
Go to the documentation of this file.
1#ifndef FIELDFILLPATCHOPS_H
2#define FIELDFILLPATCHOPS_H
3
9
10#include "AMReX_AmrCore.H"
11#include "AMReX_MultiFab.H"
12#include "AMReX_REAL.H"
13#include "AMReX_PhysBCFunct.H"
14#include "AMReX_FillPatchUtil.H"
15
31namespace amr_wind {
32
39{
40public:
42
43 virtual ~FieldFillPatchOpsBase() = default;
44
47 virtual void fillpatch(
48 const int lev,
49 const amrex::Real time,
50 amrex::MultiFab& mfab,
51 const amrex::IntVect& nghost,
52 const FieldState fstate = FieldState::New) = 0;
53
57 const int lev,
58 const amrex::Real time,
59 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& mfabs,
60 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& ffabs,
61 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& cfabs,
62 const amrex::IntVect& nghost,
63 const amrex::Vector<amrex::BCRec>& fillpatch_bcrec,
64 const amrex::Vector<amrex::BCRec>& physbc_bcrec,
65 const FieldState fstate = FieldState::New) = 0;
66
69 const int lev,
70 const amrex::Real time,
71 amrex::MultiFab& mfab,
72 const amrex::IntVect& nghost,
73 const FieldState fstate = FieldState::New) = 0;
74
76 virtual void fillphysbc(
77 const int lev,
78 const amrex::Real time,
79 amrex::MultiFab& mfab,
80 const amrex::IntVect& nghost,
81 const FieldState fstate = FieldState::New) = 0;
82
83 virtual void set_inflow(
84 const int lev,
85 const amrex::Real time,
86 amrex::MultiFab& mfab,
87 const amrex::IntVect& nghost,
88 const FieldState fstate = FieldState::New) = 0;
89
91 const int lev,
92 const amrex::Real time,
93 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM> mfabs) = 0;
94};
95
101{
102public:
103 FieldFillConstScalar(Field& /*unused*/, amrex::Real fill_val)
104 : m_fill_val(fill_val)
105 {}
106
108 int /*lev*/,
109 amrex::Real /*time*/,
110 amrex::MultiFab& mfab,
111 const amrex::IntVect& /*nghost*/,
112 const FieldState /*fstate*/) override
113 {
114 mfab.setVal(m_fill_val);
115 }
116
118 int /*lev*/,
119 amrex::Real /*time*/,
120 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& mfabs,
121 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& /*ffabs*/,
122 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& /*cfabs*/,
123 const amrex::IntVect& /*nghost*/,
124 const amrex::Vector<amrex::BCRec>& /*bcrec*/,
125 const amrex::Vector<amrex::BCRec>& /*bcrec*/,
126 const FieldState /*fstate*/) override
127 {
128 for (const auto& mfab : mfabs) {
129 mfab->setVal(m_fill_val);
130 }
131 }
132
134 int /*lev*/,
135 amrex::Real /*time*/,
136 amrex::MultiFab& mfab,
137 const amrex::IntVect& /*nghost*/,
138 const FieldState /*fstate*/) override
139 {
140 mfab.setVal(m_fill_val);
141 }
142
144 int /*lev*/,
145 amrex::Real /*time*/,
146 amrex::MultiFab& mfab,
147 const amrex::IntVect& /*nghost*/,
148 const FieldState /*fstate*/) override
149 {
150 mfab.setVal(m_fill_val);
151 }
152
154 int /*lev*/,
155 amrex::Real /*time*/,
156 amrex::MultiFab& /*mfab*/,
157 const amrex::IntVect& /*nghost*/,
158 const FieldState /*fstate*/) override
159 {
160 amrex::Abort("FieldFillConstScalar::set_inflow is not implemented");
161 }
162
163private:
164 amrex::Real m_fill_val;
165};
166
172template <typename BCOpCreator>
174{
175public:
186 Field& field,
187 const amrex::AmrCore& mesh,
188 const SimTime& time,
191 : m_time(time)
192 , m_mesh(mesh)
193 , m_field(field)
194 , m_op(field)
195 , m_mapper(field_impl::get_interpolation_operator(itype))
196 , m_face_mapper(field_impl::get_interpolation_operator(face_itype))
197 {
199 }
200
202 Field& field,
203 const amrex::AmrCore& mesh,
204 const SimTime& time,
205 const BCOpCreator& bc_op,
208 : m_time(time)
209 , m_mesh(mesh)
210 , m_field(field)
211 , m_op(bc_op)
212 , m_mapper(field_impl::get_interpolation_operator(itype))
213 , m_face_mapper(field_impl::get_interpolation_operator(face_itype))
214 {
216 }
217
224 amrex::Vector<amrex::MultiFab*> get_mfab_vec(int lev)
225 {
226 const int nstates = amrex::min(m_field.num_time_states(), 2);
227 amrex::Vector<amrex::MultiFab*> ret;
228
229 // The states in the FieldInfo data are ordered from newest to oldest,
230 // so swap the order
231 for (int i = nstates - 1; i >= 0; --i) {
232 const auto fstate = static_cast<FieldState>(i);
233 ret.push_back(&m_field.state(fstate)(lev));
234 }
235 return ret;
236 }
237
238 // Version that does no interpolation in time
239
241 int lev,
242 amrex::Real time,
243 amrex::MultiFab& mfab,
244 const amrex::IntVect& nghost,
245 const FieldState fstate = FieldState::New) override
246 {
247 auto& fld = m_field.state(fstate);
248 if (lev == 0) {
249 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> physbc(
250 m_mesh.Geom(lev), m_field.bcrec(), bc_functor());
251
252 amrex::FillPatchSingleLevel(
253 mfab, nghost, time, {&fld(lev)}, {time}, 0, 0,
254 m_field.num_comp(), m_mesh.Geom(lev), physbc, 0);
255 } else {
256 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> cphysbc(
257 m_mesh.Geom(lev - 1), m_field.bcrec(), bc_functor());
258
259 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> fphysbc(
260 m_mesh.Geom(lev), m_field.bcrec(), bc_functor());
261
262 amrex::FillPatchTwoLevels(
263 mfab, nghost, time, {&fld(lev - 1)}, {time}, {&fld(lev)},
264 {time}, 0, 0, m_field.num_comp(), m_mesh.Geom(lev - 1),
265 m_mesh.Geom(lev), cphysbc, 0, fphysbc, 0,
266 m_mesh.refRatio(lev - 1), m_mapper, m_field.bcrec(), 0);
267 }
268 }
269
271 int lev,
272 amrex::Real time,
273 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& mfabs,
274 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& ffabs,
275 amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>& cfabs,
276 const amrex::IntVect& nghost,
277 const amrex::Vector<amrex::BCRec>& fillpatch_bcrec,
278 const amrex::Vector<amrex::BCRec>& physbc_bcrec,
279 const FieldState /*fstate = FieldState::New*/) override
280 {
281
282 if (lev == 0) {
283 for (int i = 0; i < static_cast<int>(mfabs.size()); i++) {
284 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> physbc(
285 m_mesh.Geom(lev), physbc_bcrec, bc_functor_face(i));
286 amrex::FillPatchSingleLevel(
287 *mfabs[i], nghost, time, {ffabs[i]}, {time}, 0, 0, 1,
288 m_mesh.Geom(lev), physbc, i);
289 }
290 } else {
291 AMREX_D_TERM(
292 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> cphysbcx(
293 m_mesh.Geom(lev - 1), physbc_bcrec, bc_functor_face(0));
294 , amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> cphysbcy(
295 m_mesh.Geom(lev - 1), physbc_bcrec, bc_functor_face(1));
296 , amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> cphysbcz(
297 m_mesh.Geom(lev - 1), physbc_bcrec, bc_functor_face(2)););
298
299 AMREX_D_TERM(
300 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> fphysbcx(
301 m_mesh.Geom(lev), physbc_bcrec, bc_functor_face(0));
302 , amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> fphysbcy(
303 m_mesh.Geom(lev), physbc_bcrec, bc_functor_face(1));
304 , amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> fphysbcz(
305 m_mesh.Geom(lev), physbc_bcrec, bc_functor_face(2)););
306
307 amrex::Array<
308 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>>,
309 AMREX_SPACEDIM>
310 cphysbc_arr = {AMREX_D_DECL(cphysbcx, cphysbcy, cphysbcz)};
311
312 amrex::Array<
313 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>>,
314 AMREX_SPACEDIM>
315 fphysbc_arr = {AMREX_D_DECL(fphysbcx, fphysbcy, fphysbcz)};
316
317 amrex::Array<int, AMREX_SPACEDIM> idx = {AMREX_D_DECL(0, 1, 2)};
318 const amrex::Array<amrex::Vector<amrex::BCRec>, AMREX_SPACEDIM>
319 bcrec_arr = {AMREX_D_DECL(
320 fillpatch_bcrec, fillpatch_bcrec, fillpatch_bcrec)};
321
322 amrex::FillPatchTwoLevels(
323 mfabs, nghost, time, {cfabs}, {time}, {ffabs}, {time}, 0, 0, 1,
324 m_mesh.Geom(lev - 1), m_mesh.Geom(lev), cphysbc_arr, idx,
325 fphysbc_arr, idx, m_mesh.refRatio(lev - 1), m_face_mapper,
326 bcrec_arr, idx);
327 }
328 }
329
331 int lev,
332 amrex::Real time,
333 amrex::MultiFab& mfab,
334 const amrex::IntVect& nghost,
335 const FieldState fstate = FieldState::New) override
336 {
337 const auto& fld = m_field.state(fstate);
338 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> cphysbc(
339 m_mesh.Geom(lev - 1), m_field.bcrec(), bc_functor());
340
341 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> fphysbc(
342 m_mesh.Geom(lev), m_field.bcrec(), bc_functor());
343
344 amrex::InterpFromCoarseLevel(
345 mfab, nghost, time, fld(lev - 1), 0, 0, m_field.num_comp(),
346 m_mesh.Geom(lev - 1), m_mesh.Geom(lev), cphysbc, 0, fphysbc, 0,
347 m_mesh.refRatio(lev - 1), m_mapper, m_field.bcrec(), 0);
348 }
349
351 int lev,
352 amrex::Real time,
353 amrex::MultiFab& mfab,
354 const amrex::IntVect& nghost,
355 const FieldState /*fstate*/) override
356 {
357 amrex::PhysBCFunct<amrex::GpuBndryFuncFab<Functor>> physbc(
358 m_mesh.Geom(lev), m_field.bcrec(), bc_functor());
359 physbc.FillBoundary(mfab, 0, m_field.num_comp(), nghost, time, 0);
360 }
361
363 int lev,
364 amrex::Real time,
365 amrex::MultiFab& mfab,
366 const amrex::IntVect& nghost,
367 const FieldState /*fstate*/) override
368 {
369 const int ng = nghost[0];
370 const auto& bctype = m_field.bc_type();
371 const auto& geom = m_mesh.Geom(lev);
372 const auto& gdata = geom.data();
373 const auto& domain = geom.growPeriodicDomain(ng);
374 const auto& bcfunc = bc_functor();
375 const auto& ncomp = m_field.num_comp();
376
377 for (amrex::OrientationIter oit; oit != nullptr; ++oit) {
378 const auto ori = oit();
379 if ((bctype[ori] != BC::mass_inflow) &&
380 (bctype[ori] != BC::mass_inflow_outflow)) {
381 continue;
382 }
383
384 const int idir = ori.coordDir();
385 const auto& dbx =
386 ori.isLow() ? amrex::adjCellLo(domain, idir, nghost[idir])
387 : amrex::adjCellHi(domain, idir, nghost[idir]);
388
389#ifdef AMREX_USE_OMP
390#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
391#endif
392 for (amrex::MFIter mfi(mfab); mfi.isValid(); ++mfi) {
393 const auto& gbx = amrex::grow(mfi.validbox(), nghost);
394 const auto& bx = gbx & dbx;
395 if (!bx.ok()) {
396 continue;
397 }
398
399 const auto& marr = mfab[mfi].array();
400 amrex::ParallelFor(
401 bx, [=] AMREX_GPU_DEVICE(
402 const int i, const int j, const int k) noexcept {
403 for (int n = 0; n < ncomp; ++n) {
404 bcfunc.set_inflow(
405 {i, j, k}, marr, gdata, time, ori, n, 0, 0);
406 }
407 });
408 }
409 }
410 }
411
413 const int lev,
414 const amrex::Real time,
415 const amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM> mfabs) override
416 {
417 const auto& bctype = m_field.bc_type();
418 const auto& geom = m_mesh.Geom(lev);
419 const auto& gdata = geom.data();
420 const auto& bcfunc = bc_functor();
421
422 for (amrex::OrientationIter oit; oit != nullptr; ++oit) {
423 const auto ori = oit();
424 if ((bctype[ori] != BC::mass_inflow) &&
425 (bctype[ori] != BC::mass_inflow_outflow)) {
426 continue;
427 }
428
429 const int idir = ori.coordDir();
430 const auto& domain_box = geom.Domain();
431 for (int fdir = 0; fdir < AMREX_SPACEDIM; ++fdir) {
432
433 // Only face-normal velocities populated here
434 if (idir != fdir) {
435 continue;
436 }
437 const auto& dbx = ori.isLow() ? amrex::bdryLo(domain_box, idir)
438 : amrex::bdryHi(domain_box, idir);
439
440 auto& mfab = *mfabs[fdir];
441#ifdef AMREX_USE_OMP
442#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
443#endif
444 for (amrex::MFIter mfi(mfab); mfi.isValid(); ++mfi) {
445 const auto& vbx = mfi.validbox();
446 const auto& bx = vbx & dbx;
447 if (!bx.ok()) {
448 continue;
449 }
450
451 const auto& marr = mfab[mfi].array();
452 amrex::FArrayBox tmp_fab(
453 bx, mfab.nComp(), amrex::The_Async_Arena());
454 amrex::Array4<amrex::Real> const& tmp_marr =
455 tmp_fab.array();
456
457 amrex::ParallelFor(
458 bx,
459 [=] AMREX_GPU_DEVICE(
460 const int i, const int j, const int k) noexcept {
461 bcfunc.set_inflow(
462 {i, j, k}, tmp_marr, gdata, time, ori, 0, 0,
463 fdir);
464
465 const bool is_outflow =
466 (ori.isLow() && tmp_marr(i, j, k) < 0) ||
467 (ori.isHigh() && tmp_marr(i, j, k) > 0);
468 if ((bctype[ori] == BC::mass_inflow_outflow) &&
469 is_outflow) {
470 marr(i, j, k) = marr(i, j, k);
471 } else {
472 marr(i, j, k) = tmp_marr(i, j, k);
473 }
474 });
475 }
476 }
477 }
478 }
479
480protected:
481 Functor bc_functor() { return m_op(); }
482
483 Functor bc_functor_face(const int face_dir) { return m_op(face_dir); }
484
486 {
487 if (std::any_of(
488 m_mesh.refRatio().begin(), m_mesh.refRatio().end(),
489 [](amrex::IntVect iv) { return iv != amrex::IntVect(2); })) {
490 amrex::Print()
491 << "WARNING: Changing the face interpolator to FaceLinear "
492 "because a refinement ratio is different than 2"
493 << std::endl;
496 }
497 }
498
500 const amrex::AmrCore& m_mesh;
502
504
506 amrex::Interpolater* m_mapper;
507
509 amrex::Interpolater* m_face_mapper;
510};
511
512} // namespace amr_wind
513
514#endif /* FIELDFILLPATCHOPS_H */
Definition FieldFillPatchOps.H:101
void fillpatch_from_coarse(int, amrex::Real, amrex::MultiFab &mfab, const amrex::IntVect &, const FieldState) override
Implementation that handles filling patches from a coarse to fine level.
Definition FieldFillPatchOps.H:133
void set_inflow(int, amrex::Real, amrex::MultiFab &, const amrex::IntVect &, const FieldState) override
Definition FieldFillPatchOps.H:153
amrex::Real m_fill_val
Definition FieldFillPatchOps.H:164
void fillphysbc(int, amrex::Real, amrex::MultiFab &mfab, const amrex::IntVect &, const FieldState) override
Implementation that handles filling physical boundary conditions.
Definition FieldFillPatchOps.H:143
FieldFillConstScalar(Field &, amrex::Real fill_val)
Definition FieldFillPatchOps.H:103
void fillpatch_sibling_fields(int, amrex::Real, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &mfabs, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &, const amrex::IntVect &, const amrex::Vector< amrex::BCRec > &, const amrex::Vector< amrex::BCRec > &, const FieldState) override
Definition FieldFillPatchOps.H:117
void fillpatch(int, amrex::Real, amrex::MultiFab &mfab, const amrex::IntVect &, const FieldState) override
Definition FieldFillPatchOps.H:107
Definition FieldFillPatchOps.H:39
virtual void fillpatch_sibling_fields(const int lev, const amrex::Real time, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &mfabs, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &ffabs, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &cfabs, const amrex::IntVect &nghost, const amrex::Vector< amrex::BCRec > &fillpatch_bcrec, const amrex::Vector< amrex::BCRec > &physbc_bcrec, const FieldState fstate=FieldState::New)=0
virtual ~FieldFillPatchOpsBase()=default
virtual void set_inflow_sibling_fields(const int lev, const amrex::Real time, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > mfabs)=0
virtual void fillphysbc(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState fstate=FieldState::New)=0
Implementation that handles filling physical boundary conditions.
virtual void set_inflow(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState fstate=FieldState::New)=0
virtual void fillpatch_from_coarse(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState fstate=FieldState::New)=0
Implementation that handles filling patches from a coarse to fine level.
virtual void fillpatch(const int lev, const amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState fstate=FieldState::New)=0
Definition FieldFillPatchOps.H:174
amrex::Vector< amrex::MultiFab * > get_mfab_vec(int lev)
Definition FieldFillPatchOps.H:224
FieldFillPatchOps(Field &field, const amrex::AmrCore &mesh, const SimTime &time, FieldInterpolator itype=FieldInterpolator::CellConsLinear, FieldInterpolator face_itype=FieldInterpolator::FaceDivFree)
Definition FieldFillPatchOps.H:185
const SimTime & m_time
Definition FieldFillPatchOps.H:499
const BCOpCreator m_op
Definition FieldFillPatchOps.H:503
void check_face_mapper()
Definition FieldFillPatchOps.H:485
Field & m_field
Definition FieldFillPatchOps.H:501
typename BCOpCreator::FunctorType Functor
Definition FieldFillPatchOps.H:176
const amrex::AmrCore & m_mesh
Definition FieldFillPatchOps.H:500
amrex::Interpolater * m_mapper
Function that handles interpolation from coarse to fine level.
Definition FieldFillPatchOps.H:506
Functor bc_functor_face(const int face_dir)
Definition FieldFillPatchOps.H:483
void set_inflow(int lev, amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState) override
Definition FieldFillPatchOps.H:362
Functor bc_functor()
Definition FieldFillPatchOps.H:481
void set_inflow_sibling_fields(const int lev, const amrex::Real time, const amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > mfabs) override
Definition FieldFillPatchOps.H:412
void fillpatch(int lev, amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState fstate=FieldState::New) override
Definition FieldFillPatchOps.H:240
void fillpatch_from_coarse(int lev, amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState fstate=FieldState::New) override
Implementation that handles filling patches from a coarse to fine level.
Definition FieldFillPatchOps.H:330
FieldFillPatchOps(Field &field, const amrex::AmrCore &mesh, const SimTime &time, const BCOpCreator &bc_op, FieldInterpolator itype=FieldInterpolator::CellConsLinear, FieldInterpolator face_itype=FieldInterpolator::FaceDivFree)
Definition FieldFillPatchOps.H:201
amrex::Interpolater * m_face_mapper
Function that handles interpolation from coarse to fine level for faces.
Definition FieldFillPatchOps.H:509
void fillphysbc(int lev, amrex::Real time, amrex::MultiFab &mfab, const amrex::IntVect &nghost, const FieldState) override
Implementation that handles filling physical boundary conditions.
Definition FieldFillPatchOps.H:350
void fillpatch_sibling_fields(int lev, amrex::Real time, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &mfabs, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &ffabs, amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > &cfabs, const amrex::IntVect &nghost, const amrex::Vector< amrex::BCRec > &fillpatch_bcrec, const amrex::Vector< amrex::BCRec > &physbc_bcrec, const FieldState) override
Definition FieldFillPatchOps.H:270
Definition Field.H:116
int num_time_states() const
Number of exact time states available for this field.
Definition Field.H:140
const amrex::GpuArray< BC, AMREX_SPACEDIM *2 > & bc_type() const
Definition Field.H:171
int num_comp() const
Number of components for this field.
Definition Field.H:134
Field & state(const FieldState fstate)
Return field at a different time state.
Definition Field.cpp:114
amrex::Vector< amrex::BCRec > & bcrec() const
Return reference to host view of BCRec array.
Definition Field.H:188
Definition SimTime.H:30
amrex::Interpolater * get_interpolation_operator(const FieldInterpolator itype)
Definition FieldUtils.H:85
FieldState
Definition FieldDescTypes.H:14
FieldInterpolator
Definition FieldDescTypes.H:37
@ New
Same as FieldState::NP1.
@ FaceLinear
Linear face interpolation.
@ FaceDivFree
Divergence free face interpolation.
@ CellConsLinear
Linear interpolation.
@ mass_inflow_outflow
@ mass_inflow
Definition BCInterface.cpp:7
Definition FieldBCOps.H:228
DirichletOp< InflowOpType, WallOpType > FunctorType
Definition FieldBCOps.H:231