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