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