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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/vof/vof_momentum_flux.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
vof_momentum_flux.H
Go to the documentation of this file.
1#ifndef VOF_MOMENTUM_FLUX_H
2#define VOF_MOMENTUM_FLUX_H
3
5
6namespace amr_wind::multiphase {
7
8static void hybrid_fluxes(
9 const FieldRepo& repo,
10 const int ncomp,
11 const amrex::Gpu::DeviceVector<int>& iconserv,
12 ScratchField& flux_x,
13 ScratchField& flux_y,
14 ScratchField& flux_z,
15 const Field& dof_field,
16 const Field& dof_nph,
17 const Field& src_term,
18 const Field& rho_o,
19 const Field& rho_nph,
20 const Field& u_mac,
21 const Field& v_mac,
22 const Field& w_mac,
23 amrex::Vector<amrex::BCRec> const& velbc,
24 amrex::BCRec const* velbc_d,
25 amrex::Vector<amrex::BCRec> const& rhobc,
26 amrex::BCRec const* rhobc_d,
27 const amrex::Real dt,
28 godunov::scheme mflux_scheme,
29 bool allow_inflow_on_outflow,
30 bool use_forces_in_trans)
31{
32 // Get geometry
33 const auto& geom = repo.mesh().Geom();
34 // Get advected alpha fields
35 // At this point in the solve, they have been converted to advected density
36 const auto& advrho_x = repo.get_field("advalpha_x");
37 const auto& advrho_y = repo.get_field("advalpha_y");
38 const auto& advrho_z = repo.get_field("advalpha_z");
39 // Get VOF
40 const auto& vof = repo.get_field("vof").state(amr_wind::FieldState::Old);
41
42 // Create scratch arrays for local flux storage
43 auto ftmp_x =
45 auto ftmp_y =
47 auto ftmp_z =
49
50 auto face_x =
52 auto face_y =
54 auto face_z =
56
57 // Create iconserv = 0 array for density, avoid multiplication by face vel
58 amrex::Gpu::DeviceVector<int> idnsty;
59 idnsty.resize(1, -1);
60
61 for (int lev = 0; lev < repo.num_active_levels(); ++lev) {
62 // Set up temporary arrays for fluxes that will replace other quantities
63 // form multifab for transport variable and source term
64 amrex::MultiFab q(
65 dof_field(lev).boxArray(), dof_field(lev).DistributionMap(), ncomp,
67 amrex::MultiFab::Copy(
68 q, dof_field(lev), 0, 0, ncomp, fvm::Godunov::nghost_state);
69 amrex::MultiFab fq(
70 src_term(lev).boxArray(), src_term(lev).DistributionMap(), ncomp,
72 amrex::MultiFab::Copy(
73 fq, src_term(lev), 0, 0, ncomp, fvm::Godunov::nghost_src);
74
75 for (int idim = 0; idim < dof_field.num_comp(); ++idim) {
76 amrex::MultiFab::Multiply(
77 q, rho_o(lev), 0, idim, 1, fvm::Godunov::nghost_state);
78 amrex::MultiFab::Multiply(
79 fq, rho_o(lev), 0, idim, 1, fvm::Godunov::nghost_src);
80 }
81
82 amrex::MultiFab frho(
83 src_term(lev).boxArray(), src_term(lev).DistributionMap(), 1,
85 frho.setVal(0.0);
86
87 // Set up temporary arrays for NPH at boundaries
88 amrex::MultiFab q_nph(
89 dof_field(lev).boxArray(), dof_field(lev).DistributionMap(), ncomp,
91 amrex::MultiFab::Copy(
92 q_nph, dof_nph(lev), 0, 0, ncomp, fvm::Godunov::nghost_state);
93 // Density at NPH is fully calculated a priori
94 // Multiply to get momentum at NPH
95 for (int idim = 0; idim < dof_field.num_comp(); ++idim) {
96 amrex::MultiFab::Multiply(q_nph, rho_nph(lev), 0, idim, 1, 1);
97 }
98
99 std::string advection_type = "Godunov";
100 int limiter_type = PPM::UPWIND; // default
101 if (mflux_scheme == godunov::scheme::MINMOD) {
102 limiter_type = PPM::MINMOD;
103 } else if (mflux_scheme != godunov::scheme::UPWIND) {
104 amrex::Abort("Unknown limiter type in vof_momentum_flux.H");
105 }
106
107 const bool known_edge_state = false;
108 const bool godunov_use_ppm = true;
109 const bool is_velocity = false;
110 const bool fluxes_are_area_weighted = false;
111
112 amrex::MFItInfo mfi_info;
113 if (amrex::Gpu::notInLaunchRegion()) {
114 mfi_info.EnableTiling(amrex::IntVect(1024, 1024, 1024))
115 .SetDynamic(true);
116 }
117#ifdef AMREX_USE_OMP
118#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
119#endif
120 for (amrex::MFIter mfi(dof_field(lev), mfi_info); mfi.isValid();
121 ++mfi) {
122 const auto& bx = mfi.tilebox();
123
124 // Arrays for easy reference (w for working array)
125 auto F_x = flux_x(lev).array(mfi);
126 auto F_y = flux_y(lev).array(mfi);
127 auto F_z = flux_z(lev).array(mfi);
128 auto Fw_x = (*ftmp_x)(lev).array(mfi);
129 auto Fw_y = (*ftmp_y)(lev).array(mfi);
130 auto Fw_z = (*ftmp_z)(lev).array(mfi);
131
132 // Temporary divu
133 amrex::FArrayBox divufab(amrex::grow(bx, 1), 1);
134 divufab.setVal<amrex::RunOn::Device>(0.0);
135 const auto& divu = divufab.array();
136
137 // Compute momentum flux quantities (multiplied by face velocity)
138 HydroUtils::ComputeFluxesOnBoxFromState(
139 bx, ncomp, mfi, q.const_array(mfi), q_nph.const_array(mfi),
140 Fw_x, Fw_y, Fw_z, (*face_x)(lev).array(mfi),
141 (*face_y)(lev).array(mfi), (*face_z)(lev).array(mfi),
142 known_edge_state, u_mac(lev).const_array(mfi),
143 v_mac(lev).const_array(mfi), w_mac(lev).const_array(mfi), divu,
144 fq.const_array(mfi), geom[lev], dt, velbc, velbc_d,
145 iconserv.data(), godunov_use_ppm, use_forces_in_trans,
146 is_velocity, fluxes_are_area_weighted, advection_type,
147 limiter_type, allow_inflow_on_outflow);
148
149 const auto& bxg1 = amrex::grow(bx, 1);
150 const auto& xbx = amrex::surroundingNodes(bx, 0);
151 const auto& ybx = amrex::surroundingNodes(bx, 1);
152 const auto& zbx = amrex::surroundingNodes(bx, 2);
153
154 // Where interface is present, replace current momentum flux
155 // quantities with those from mflux_scheme
156 auto volfrac = vof(lev).const_array(mfi);
157 amrex::ParallelFor(
158 bxg1, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
159 bool vf_i = multiphase::interface_band(i, j, k, volfrac);
160 // Fluxes in each direction, check with flag too
161 if (xbx.contains(i, j, k)) {
162 bool vf_nb =
163 multiphase::interface_band(i - 1, j, k, volfrac);
164 if (vf_i || vf_nb) {
165 for (int n = 0; n < ncomp; ++n) {
166 F_x(i, j, k, n) = Fw_x(i, j, k, n);
167 }
168 }
169 }
170 if (ybx.contains(i, j, k)) {
171 bool vf_nb =
172 multiphase::interface_band(i, j - 1, k, volfrac);
173 if (vf_i || vf_nb) {
174 for (int n = 0; n < ncomp; ++n) {
175 F_y(i, j, k, n) = Fw_y(i, j, k, n);
176 }
177 }
178 }
179 if (zbx.contains(i, j, k)) {
180 bool vf_nb =
181 multiphase::interface_band(i, j, k - 1, volfrac);
182 if (vf_i || vf_nb) {
183 for (int n = 0; n < ncomp; ++n) {
184 F_z(i, j, k, n) = Fw_z(i, j, k, n);
185 }
186 }
187 }
188 });
189
190 // Compute density flux quantities (not multiplied by face velocity)
191 HydroUtils::ComputeFluxesOnBoxFromState(
192 bx, 1, mfi, rho_o(lev).const_array(mfi), Fw_x, Fw_y, Fw_z,
193 (*face_x)(lev).array(mfi), (*face_y)(lev).array(mfi),
194 (*face_z)(lev).array(mfi), known_edge_state,
195 u_mac(lev).const_array(mfi), v_mac(lev).const_array(mfi),
196 w_mac(lev).const_array(mfi), divu, frho.const_array(mfi),
197 geom[lev], dt, rhobc, rhobc_d, idnsty.data(), godunov_use_ppm,
198 use_forces_in_trans, is_velocity, fluxes_are_area_weighted,
199 advection_type, limiter_type, allow_inflow_on_outflow);
200
201 // At this point in the solve, these are advected
202 auto ar_x = advrho_x(lev).const_array(mfi);
203 auto ar_y = advrho_y(lev).const_array(mfi);
204 auto ar_z = advrho_z(lev).const_array(mfi);
205
206 // When interface is present, divide by interpolated density to get
207 // Favre-averaged flux value. Multiply all fluxes with advected
208 // face density to create mass-consistent fluxes
209 amrex::ParallelFor(
210 bxg1, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
211 bool vf_i = multiphase::interface_band(i, j, k, volfrac);
212 // Fluxes in each direction, check with flag too
213 if (xbx.contains(i, j, k)) {
214 bool vf_nb =
215 multiphase::interface_band(i - 1, j, k, volfrac);
216 for (int n = 0; n < ncomp; ++n) {
217 if (vf_i || vf_nb) {
218 F_x(i, j, k, n) /= Fw_x(i, j, k, 0);
219 }
220 F_x(i, j, k, n) *= ar_x(i, j, k);
221 }
222 }
223 if (ybx.contains(i, j, k)) {
224 bool vf_nb =
225 multiphase::interface_band(i, j - 1, k, volfrac);
226 for (int n = 0; n < ncomp; ++n) {
227 if (vf_i || vf_nb) {
228 F_y(i, j, k, n) /= Fw_y(i, j, k, 0);
229 }
230 F_y(i, j, k, n) *= ar_y(i, j, k);
231 }
232 }
233 if (zbx.contains(i, j, k)) {
234 bool vf_nb =
235 multiphase::interface_band(i, j, k - 1, volfrac);
236 for (int n = 0; n < ncomp; ++n) {
237 if (vf_i || vf_nb) {
238 F_z(i, j, k, n) /= Fw_z(i, j, k, 0);
239 }
240 F_z(i, j, k, n) *= ar_z(i, j, k);
241 }
242 }
243 });
244
245 amrex::Gpu::streamSynchronize();
246 }
247 }
248}
249} // namespace amr_wind::multiphase
250
251#endif
Definition Field.H:116
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
Definition FieldRepo.H:86
std::unique_ptr< ScratchField > create_scratch_field(const std::string &name, const int ncomp=1, const int nghost=0, const FieldLoc floc=FieldLoc::CELL) const
Definition FieldRepo.cpp:302
Field & get_field(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:149
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
const amrex::AmrCore & mesh() const
Return a reference to the underlying AMR mesh instance.
Definition FieldRepo.H:358
Definition ScratchField.H:30
@ ZFACE
Face-centered in z-direction.
@ XFACE
Face-centered in x-direction (e.g., face normal velocity)
@ YFACE
Face-centered in y-direction.
@ Old
Same as FieldState::N.
Definition height_functions.H:8
static void hybrid_fluxes(const FieldRepo &repo, const int ncomp, const amrex::Gpu::DeviceVector< int > &iconserv, ScratchField &flux_x, ScratchField &flux_y, ScratchField &flux_z, const Field &dof_field, const Field &dof_nph, const Field &src_term, const Field &rho_o, const Field &rho_nph, const Field &u_mac, const Field &v_mac, const Field &w_mac, amrex::Vector< amrex::BCRec > const &velbc, amrex::BCRec const *velbc_d, amrex::Vector< amrex::BCRec > const &rhobc, amrex::BCRec const *rhobc_d, const amrex::Real dt, godunov::scheme mflux_scheme, bool allow_inflow_on_outflow, bool use_forces_in_trans)
Definition vof_momentum_flux.H:8
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool interface_band(const int i, const int j, const int k, amrex::Array4< amrex::Real const > const &volfrac, const int n_band=1, const amrex::Real tiny=constants::TIGHT_TOL) noexcept
Definition volume_fractions.H:492
scheme
Definition Godunov.H:11
static constexpr int nghost_state
Number of ghost in the state variable.
Definition SchemeTraits.H:19
static constexpr int nghost_src
Number of ghost cells in the source term variable.
Definition SchemeTraits.H:21