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