/home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/vof/vof_advection.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/vof/vof_advection.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
vof_advection.H
Go to the documentation of this file.
1
2#ifndef VOF_ADVECTION_H
3#define VOF_ADVECTION_H
4
9
10namespace amr_wind::pde {
11
15template <>
16struct AdvectionOp<VOF, fvm::Godunov>
17{
19 CFDSim& sim,
20 PDEFields& fields_in,
21 bool /*unused*/,
22 bool /*unused*/,
23 bool /*unused*/,
24 bool /*unused*/)
25 : m_time(sim.time())
26 , fields(fields_in)
27 , u_mac(fields_in.repo.get_field("u_mac"))
28 , v_mac(fields_in.repo.get_field("v_mac"))
29 , w_mac(fields_in.repo.get_field("w_mac"))
30 {
31 amrex::ParmParse pp_multiphase("VOF");
32 pp_multiphase.query("remove_debris", m_rm_debris);
33 pp_multiphase.query("replace_masked", m_replace_mask);
34
35 // Setup density factor arrays for multiplying velocity flux
37 {"advalpha_x", "advalpha_y", "advalpha_z"}, 1,
39 }
40
42 const FieldState /*unused*/,
43 const amrex::Real /*unused*/,
44 const amrex::Real /*unused*/)
45 {}
46
47 void operator()(const FieldState /*unused*/, const amrex::Real dt)
48 {
49 static_assert(
50 VOF::ndim == 1, "Invalid number of components for scalar");
51
52 auto& repo = fields.repo;
53 const auto& geom = repo.mesh().Geom();
54 const int nlevels = repo.num_active_levels();
55
56 auto& aa_x = repo.get_field("advalpha_x");
57 auto& aa_y = repo.get_field("advalpha_y");
58 auto& aa_z = repo.get_field("advalpha_z");
59
60 // Old and new states
61 auto& dof_old = fields.field.state(amr_wind::FieldState::Old);
62 auto& dof_new = fields.field;
63 // Working state of vof is nph, to keep others untouched during step
64 auto& dof_field = fields.field.state(amr_wind::FieldState::NPH);
65
66 // Initialize as old state values
67 for (int lev = 0; lev < repo.num_active_levels(); ++lev) {
68 amrex::MultiFab::Copy(
69 dof_field(lev), dof_old(lev), 0, 0, dof_field.num_comp(),
70 dof_field.num_grow());
71 }
72
73 //
74 // Advect volume using Implicit Eulerian Sweeping method with PLIC
75 // reconstruction
76 //
77
78 auto flux_x =
79 repo.create_scratch_field(1, 0, amr_wind::FieldLoc::XFACE);
80 auto flux_y =
81 repo.create_scratch_field(1, 0, amr_wind::FieldLoc::YFACE);
82 auto flux_z =
83 repo.create_scratch_field(1, 0, amr_wind::FieldLoc::ZFACE);
84
85 // Scratch field for fluxC
86 auto fluxC = repo.create_scratch_field(1, 0, amr_wind::FieldLoc::CELL);
87
88 // Define the sweep time
89 isweep += 1;
90 if (isweep > 3) {
91 isweep = 1;
92 }
93
94 amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>> fluxes(
95 nlevels);
96 amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>> advas(
97 nlevels);
98 for (int lev = 0; lev < nlevels; ++lev) {
99 fluxes[lev][0] = &(*flux_x)(lev);
100 fluxes[lev][1] = &(*flux_y)(lev);
101 fluxes[lev][2] = &(*flux_z)(lev);
102 advas[lev][0] = &aa_x(lev);
103 advas[lev][1] = &aa_y(lev);
104 advas[lev][2] = &aa_z(lev);
105 }
106
107 // Split advection step 1, with cmask calculation
109 isweep, 0, nlevels, dof_field, fluxes, (*fluxC), advas, u_mac,
110 v_mac, w_mac, dof_field.bc_type(), geom, m_time.new_time(), dt,
111 m_rm_debris);
112 // (copy old boundaries to working state)
113 // Split advection step 2
115 isweep, 1, nlevels, dof_field, fluxes, (*fluxC), advas, u_mac,
116 v_mac, w_mac, dof_field.bc_type(), geom, m_time.new_time(), dt,
117 m_rm_debris);
118 // (copy old boundaries to working state)
119 // Split advection step 3
121 isweep, 2, nlevels, dof_field, fluxes, (*fluxC), advas, u_mac,
122 v_mac, w_mac, dof_field.bc_type(), geom, m_time.new_time(), dt,
123 m_rm_debris);
124
125 // Replace masked cells using overset
126 if (repo.int_field_exists("iblank_cell") && m_replace_mask) {
127 auto& f_iblank = repo.get_int_field("iblank_cell");
129 nlevels, f_iblank, dof_field, dof_new);
130 }
131
132 // Copy working version of vof to new state
133 for (int lev = 0; lev < repo.num_active_levels(); ++lev) {
134 amrex::MultiFab::Copy(
135 dof_new(lev), dof_field(lev), 0, 0, dof_field.num_comp(),
136 dof_field.num_grow());
137 }
138 }
139
145 int isweep = 0;
146 bool m_rm_debris{true};
147 bool m_replace_mask{true};
148 // Lagrangian transport is deprecated, only Eulerian is supported
149};
150
151} // namespace amr_wind::pde
152#endif
Definition CFDSim.H:54
Definition Field.H:116
amrex::Vector< Field * > declare_face_normal_field(const amrex::Vector< std::string > &names, const int ncomp=1, const int ngrow=0, const int nstates=1)
Definition FieldRepo.H:225
Definition SimTime.H:30
FieldState
Definition FieldDescTypes.H:14
@ ZFACE
Face-centered in z-direction.
@ XFACE
Face-centered in x-direction (e.g., face normal velocity)
@ CELL
Cell-centered (default)
@ YFACE
Face-centered in y-direction.
@ NPH
State at (n + 1/2) (intermediate) timestep.
@ Old
Same as FieldState::N.
void split_advection_step(int isweep, int iorder, int nlevels, Field &dof_field, amrex::Vector< amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > > const &fluxes, ScratchField &fluxC, amrex::Vector< amrex::Array< amrex::MultiFab *, AMREX_SPACEDIM > > const &advas, Field const &u_mac, Field const &v_mac, Field const &w_mac, amrex::GpuArray< BC, AMREX_SPACEDIM *2 > BCs, amrex::Vector< amrex::Geometry > geom, amrex::Real time, amrex::Real dt, bool rm_debris)
Definition SplitAdvection.cpp:11
static void replace_masked_vof(const int nlevels, amr_wind::IntField &f_iblank, amr_wind::Field &f_vof, amr_wind::Field &f_vof_new)
Definition vof_hybridsolver_ops.H:10
Definition AdvOp_Godunov.H:16
static constexpr int nghost_mac
Number of ghost cells in the MAC face variables.
Definition SchemeTraits.H:23
Field & u_mac
Definition vof_advection.H:142
AdvectionOp(CFDSim &sim, PDEFields &fields_in, bool, bool, bool, bool)
Definition vof_advection.H:18
PDEFields & fields
Definition vof_advection.H:141
Field & w_mac
Definition vof_advection.H:144
Field & v_mac
Definition vof_advection.H:143
void preadvect(const FieldState, const amrex::Real, const amrex::Real)
Definition vof_advection.H:41
const SimTime & m_time
Definition vof_advection.H:140
void operator()(const FieldState, const amrex::Real dt)
Definition vof_advection.H:47
Definition PDEOps.H:163
Definition PDEFields.H:27
FieldRepo & repo
Reference to the field repository instance.
Definition PDEFields.H:31
Definition vof.H:17
static constexpr int ndim
Definition vof.H:24