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

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/equation_systems/AdvOp_Godunov.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
AdvOp_Godunov.H
Go to the documentation of this file.
1#ifndef ADVOP_GODUNOV_H
2#define ADVOP_GODUNOV_H
3
4#include <type_traits>
5
10
11#include "AMReX_Gpu.H"
12#include "AMReX_ParmParse.H"
13#include "AMReX_MultiFabUtil.H"
14#include "hydro_utils.H"
15
17#include "AMReX_REAL.H"
18
19using namespace amrex::literals;
20
21namespace amr_wind::pde {
22
26template <typename PDE>
28 PDE,
29 fvm::Godunov,
30 typename std::enable_if_t<std::is_base_of_v<ScalarTransport, PDE>>>
31{
33 CFDSim& /* sim */,
34 PDEFields& fields_in,
35 bool /* has_overset */,
36 bool /* variable density */,
37 bool /* mesh mapping */,
38 bool /* is_anelastic */)
39 : fields(fields_in)
40 , density(fields_in.repo.get_field("density"))
41 , u_mac(fields_in.repo.get_field("u_mac"))
42 , v_mac(fields_in.repo.get_field("v_mac"))
43 , w_mac(fields_in.repo.get_field("w_mac"))
44 {
45 amrex::ParmParse pp("incflo");
46 pp.query("godunov_type", godunov_type);
47 pp.query("godunov_use_forces_in_trans", godunov_use_forces_in_trans);
48 if (pp.contains("use_ppm") || pp.contains("use_limiter")) {
49 amrex::Abort(
50 "Godunov: use_ppm and use_limiter are deprecated. Please "
51 "update input file");
52 }
53
54 if (amrex::toLower(godunov_type) == "plm") {
56 } else if (amrex::toLower(godunov_type) == "ppm") {
58 } else if (amrex::toLower(godunov_type) == "ppm_nolim") {
60 amrex::Print() << "WARNING: Using advection type ppm_nolim is not "
61 "recommended. Prefer using weno_z."
62 << std::endl;
63 } else if (amrex::toLower(godunov_type) == "bds") {
65 advection_type = "BDS";
66 } else if (
67 amrex::toLower(godunov_type) == "weno" ||
68 amrex::toLower(godunov_type) == "weno_js") {
70 } else if (amrex::toLower(godunov_type) == "weno_z") {
72 } else {
73 amrex::Abort(
74 "Invalid godunov_type specified. For godunov_type select "
75 "between plm, ppm, ppm_nolim, bds, weno_js, and weno_z. If no "
76 "godunov_type is specified, the default weno_z is used.");
77 }
78
79 // Flux calculation used in multiphase portions of domain
80 pp.query("mflux_type", mflux_type);
81 if (amrex::toLower(mflux_type) == "minmod") {
83 } else if (amrex::toLower(mflux_type) == "upwind") {
85 } else {
86 amrex::Abort("Invalid argument entered for mflux_type.");
87 }
88
89 amrex::ParmParse pp_eq{fields_in.field.base_name()};
90 pp_eq.query(
91 "allow_inflow_at_pressure_outflow", m_allow_inflow_on_outflow);
92
93 // TODO: Need iconserv flag to be adjusted???
94 iconserv.resize(PDE::ndim, 1);
95 }
96
98 const FieldState /*unused*/,
99 const amrex::Real /*unused*/,
100 const amrex::Real /*unused*/)
101 {}
102
103 void operator()(const FieldState fstate, const amrex::Real dt)
104 {
105 static_assert(
106 PDE::ndim == 1, "Invalid number of components for scalar");
107 auto& repo = fields.repo;
108 const auto& geom = repo.mesh().Geom();
109
110 const auto& src_term = fields.src_term;
111 // cppcheck-suppress constVariableReference
112 auto& conv_term = fields.conv_term;
113 const auto& dof_field = fields.field.state(fstate);
114 const auto& dof_nph = fields.field.state(amr_wind::FieldState::NPH);
115
116 auto flux_x =
117 repo.create_scratch_field(PDE::ndim, 0, amr_wind::FieldLoc::XFACE);
118 auto flux_y =
119 repo.create_scratch_field(PDE::ndim, 0, amr_wind::FieldLoc::YFACE);
120 auto flux_z =
121 repo.create_scratch_field(PDE::ndim, 0, amr_wind::FieldLoc::ZFACE);
122 auto face_x =
123 repo.create_scratch_field(PDE::ndim, 0, amr_wind::FieldLoc::XFACE);
124 auto face_y =
125 repo.create_scratch_field(PDE::ndim, 0, amr_wind::FieldLoc::YFACE);
126 auto face_z =
127 repo.create_scratch_field(PDE::ndim, 0, amr_wind::FieldLoc::ZFACE);
128
129 // only needed if multiplying by rho below
130 const auto& den = density.state(fstate);
131 const auto& den_nph = density.state(amr_wind::FieldState::NPH);
132
133 const bool mphase_vof = repo.field_exists("vof");
134
135 for (int lev = 0; lev < repo.num_active_levels(); ++lev) {
136 amrex::MFItInfo mfi_info;
137 if (amrex::Gpu::notInLaunchRegion()) {
138 mfi_info.EnableTiling(amrex::IntVect(1024, 1024, 1024))
139 .SetDynamic(true);
140 }
141#ifdef AMREX_USE_OMP
142#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
143#endif
144 for (amrex::MFIter mfi(dof_field(lev), mfi_info); mfi.isValid();
145 ++mfi) {
146 const auto& bx = mfi.tilebox();
147 auto rho_arr = den(lev).array(mfi);
148 auto rho_nph_arr = den_nph(lev).array(mfi);
149 auto tra_arr = dof_field(lev).array(mfi);
150 auto tra_nph_arr = dof_nph(lev).array(mfi);
151 amrex::FArrayBox rhotracfab;
152 amrex::Array4<amrex::Real> rhotrac;
153 amrex::FArrayBox rhotrac_nph_fab;
154 amrex::Array4<amrex::Real> rhotrac_nph;
155 amrex::FArrayBox srctrac_over_rho_fab;
156 amrex::Array4<amrex::Real> srctrac_over_rho;
157 const auto& src_arr = src_term(lev).const_array(mfi);
158
159 if (PDE::multiply_rho) {
160 auto rhotrac_box =
161 amrex::grow(bx, fvm::Godunov::nghost_state);
162 auto srctrac_box =
163 amrex::grow(bx, fvm::Godunov::nghost_src);
164 rhotracfab.resize(
165 rhotrac_box, PDE::ndim, amrex::The_Async_Arena());
166 rhotrac = rhotracfab.array();
167 rhotrac_nph_fab.resize(
168 rhotrac_box, PDE::ndim, amrex::The_Async_Arena());
169 rhotrac_nph = rhotrac_nph_fab.array();
170 if (mphase_vof) {
171 srctrac_over_rho_fab.resize(
172 rhotrac_box, PDE::ndim, amrex::The_Async_Arena());
173 srctrac_over_rho = srctrac_over_rho_fab.array();
174 }
175
176 amrex::ParallelFor(
177 rhotrac_box, PDE::ndim,
178 [=] AMREX_GPU_DEVICE(
179 int i, int j, int k, int n) noexcept {
180 rhotrac(i, j, k, n) =
181 rho_arr(i, j, k) * tra_arr(i, j, k, n);
182 rhotrac_nph(i, j, k, n) =
183 rho_nph_arr(i, j, k) * tra_nph_arr(i, j, k, n);
184 if (mphase_vof && srctrac_box.contains(i, j, k)) {
185 // Need to remove the effect of the multiply_rho
186 // routine -- the setup is different for
187 // mass-consistent multiphase advection
188 srctrac_over_rho(i, j, k, n) =
189 src_arr(i, j, k, n) / rho_nph_arr(i, j, k);
190 }
191 });
192 }
193
200 amrex::FArrayBox tmpfab(
201 amrex::grow(bx, 1), 1, amrex::The_Async_Arena());
202 tmpfab.setVal<amrex::RunOn::Device>(0.0_rt);
203 const auto& divu = tmpfab.array();
204 const bool is_velocity = false;
205 const bool known_edge_state = false;
206 const bool godunov_use_ppm =
209 int limiter_type;
211 limiter_type = PPM::NoLimiter;
213 limiter_type = PPM::WENOZ;
215 limiter_type = PPM::WENO_JS;
216 } else {
217 limiter_type = PPM::default_limiter;
218 }
219
220 // if state is NPH, then n and n+1 are known, and only
221 // spatial extrapolation is performed
222 const amrex::Real dt_extrap =
223 (fstate == FieldState::NPH) ? 0.0_rt : dt;
224
225 HydroUtils::ComputeFluxesOnBoxFromState(
226 bx, PDE::ndim, mfi,
227 ((PDE::multiply_rho && !mphase_vof) ? rhotrac
228 : tra_arr),
229 ((PDE::multiply_rho && !mphase_vof) ? rhotrac_nph
230 : tra_nph_arr),
231 (*flux_x)(lev).array(mfi), (*flux_y)(lev).array(mfi),
232 (*flux_z)(lev).array(mfi), (*face_x)(lev).array(mfi),
233 (*face_y)(lev).array(mfi), (*face_z)(lev).array(mfi),
234 known_edge_state, u_mac(lev).const_array(mfi),
235 v_mac(lev).const_array(mfi),
236 w_mac(lev).const_array(mfi), divu,
237 ((PDE::multiply_rho && mphase_vof)
238 ? srctrac_over_rho
239 : src_term(lev).const_array(mfi)),
240 geom[lev], dt_extrap, dof_field.bcrec(),
241 dof_field.bcrec_device().data(), iconserv.data(),
242 godunov_use_ppm, godunov_use_forces_in_trans,
244 limiter_type, m_allow_inflow_on_outflow);
245 } else {
246 amrex::Abort("Invalid godunov scheme");
247 }
248 }
249 }
250
251 // Multiphase flux operations
252 if (mphase_vof && PDE::multiply_rho) {
253 // Loop levels
255 repo, PDE::ndim, iconserv, (*flux_x), (*flux_y), (*flux_z),
256 dof_field, dof_nph, src_term, den, den_nph, u_mac, v_mac, w_mac,
257 dof_field.bcrec(), dof_field.bcrec_device().data(), den.bcrec(),
258 den.bcrec_device().data(), dt, mflux_scheme,
260 }
261
262 amrex::Vector<amrex::Array<amrex::MultiFab*, AMREX_SPACEDIM>> fluxes(
263 repo.num_active_levels());
264 for (int lev = 0; lev < repo.num_active_levels(); ++lev) {
265 fluxes[lev][0] = &(*flux_x)(lev);
266 fluxes[lev][1] = &(*flux_y)(lev);
267 fluxes[lev][2] = &(*flux_z)(lev);
268 }
269
270 // In order to enforce conservation across coarse-fine boundaries we
271 // must be sure to average down the fluxes before we use them
272 for (int lev = repo.num_active_levels() - 1; lev > 0; --lev) {
273 amrex::IntVect rr =
274 geom[lev].Domain().size() / geom[lev - 1].Domain().size();
275 amrex::average_down_faces(
276 GetArrOfConstPtrs(fluxes[lev]), fluxes[lev - 1], rr,
277 geom[lev - 1]);
278 }
279
280 for (int lev = 0; lev < repo.num_active_levels(); ++lev) {
281#ifdef AMREX_USE_OMP
282#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
283#endif
284 for (amrex::MFIter mfi(dof_field(lev), amrex::TilingIfNotGPU());
285 mfi.isValid(); ++mfi) {
286 const auto& bx = mfi.tilebox();
287
288 HydroUtils::ComputeDivergence(
289 bx, conv_term(lev).array(mfi), (*flux_x)(lev).array(mfi),
290 (*flux_y)(lev).array(mfi), (*flux_z)(lev).array(mfi),
291 PDE::ndim, geom[lev], amrex::Real(-1.0_rt),
293 }
294 }
295 }
296
302 amrex::Gpu::DeviceVector<int> iconserv;
303
306 std::string godunov_type{"weno_z"};
307 std::string mflux_type{"upwind"};
308 const bool fluxes_are_area_weighted{false};
311 std::string advection_type{"Godunov"};
312};
313
314} // namespace amr_wind::pde
315
316#endif /* ADVOP_GODUNOV_H */
Definition CFDSim.H:54
Definition Field.H:116
const std::string & base_name() const
Base name (without state information)
Definition Field.H:128
FieldState
Definition FieldDescTypes.H:14
@ ZFACE
Face-centered in z-direction.
Definition FieldDescTypes.H:32
@ XFACE
Face-centered in x-direction (e.g., face normal velocity)
Definition FieldDescTypes.H:30
@ YFACE
Face-centered in y-direction.
Definition FieldDescTypes.H:31
@ NPH
State at (n + 1/2) (intermediate) timestep.
Definition FieldDescTypes.H:18
Definition SchemeTraits.H:6
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, bool pre_multiplied_src_term=false)
Definition vof_momentum_flux.H:11
Definition AdvOp_Godunov.H:21
scheme
Definition Godunov.H:11
@ BDS
Definition Godunov.H:11
@ PPM
Definition Godunov.H:11
@ UPWIND
Definition Godunov.H:11
@ PLM
Definition Godunov.H:11
@ WENO_JS
Definition Godunov.H:11
@ MINMOD
Definition Godunov.H:11
@ PPM_NOLIM
Definition Godunov.H:11
@ WENOZ
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
void operator()(const FieldState fstate, const amrex::Real dt)
Definition AdvOp_Godunov.H:103
AdvectionOp(CFDSim &, PDEFields &fields_in, bool, bool, bool, bool)
Definition AdvOp_Godunov.H:32
void preadvect(const FieldState, const amrex::Real, const amrex::Real)
Definition AdvOp_Godunov.H:97
Definition PDEFields.H:27
Field & field
Solution variable (e.g., velocity, temperature)
Definition PDEFields.H:34