DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
EulerP_Evaluator_impl.cpp
Go to the documentation of this file.
1/** @file EulerP_Evaluator_impl.cpp
2 * @brief Host-backend explicit specializations of Evaluator_impl kernel methods.
3 *
4 * Provides the DeviceBackend::Host specialization for each Evaluator_impl static method.
5 * Each method contains an OpenMP-parallelized loop that calls the corresponding
6 * device-callable kernel function from EulerP_Evaluator_impl_common.hxx.
7 *
8 * For the CUDA backend specializations, see the separate CUDA compilation unit.
9 */
10#include "DNDS/Defines.hpp"
12#include "EulerP/EulerP.hpp"
15
16namespace DNDS::EulerP
17{
18 static constexpr DeviceBackend B = DeviceBackend::Host;
19
20 /**
21 * @brief Host specialization: Green-Gauss gradient reconstruction.
22 *
23 * First loop (serial): generates boundary ghost values for all boundary faces.
24 * Second loop (OpenMP parallel): computes the Green-Gauss cell gradient for all owned cells.
25 */
26 template <>
28 {
29 using namespace Geom;
30
32
39
40 int nVarsScalar = uScalar.size();
41 int nVars = nVarsFlow + nVarsScalar;
42
43 auto &mesh = self_view.fv.mesh;
44 auto &fv = self_view.fv;
45 auto &bcHandler = self_view.bcHandler;
46 auto &phy = self_view.physics;
47
48 DNDS_check_throw(faceBCBuffer.father.data());
49 DNDS_check_throw(faceBCBuffer.father.Size() >= mesh.NumFace());
50
51 /*********************** */
52 // bc handling
53 for (index iBnd = 0; iBnd < mesh.NumBnd(); iBnd++)
54 {
55 RecGradient_GGRec_Kernel_BndVal(self_view, arg.portable, iBnd, mesh.NumBnd(), nVars, nVarsScalar);
56 }
57
58 /*********************** */
59 // rec
60#if defined(DNDS_DIST_MT_USE_OMP)
61# pragma omp parallel for schedule(runtime)
62#endif
63 for (index iCell = 0; iCell < mesh.NumCell(); iCell++)
64 {
65 RecGradient_GGRec_Kernel_GG(self_view, arg.portable, iCell, mesh.NumCell(), nVars, nVarsScalar);
66 }
67 }
68
69 /**
70 * @brief Host specialization: Barth-Jespersen gradient limiter.
71 *
72 * First loop (OpenMP parallel): applies the flow-variable limiter to all owned cells.
73 * Second loop (OpenMP parallel, conditional): applies the scalar-variable limiter if
74 * there are transported scalars.
75 */
76 template <>
78 {
79 using namespace Geom;
80
86
87 int nVarsScalar = uScalar.size();
88 int nVars = nVarsFlow + nVarsScalar;
89
90 auto &mesh = self_view.fv.mesh;
91 auto &fv = self_view.fv;
92 auto &bcHandler = self_view.bcHandler;
93 auto &phy = self_view.physics;
94
95 /*********************** */
96 // limit
97
98#if defined(DNDS_DIST_MT_USE_OMP)
99# pragma omp parallel for schedule(runtime)
100#endif
101 for (index iCell = 0; iCell < mesh.NumCell(); iCell++)
102 {
103 RecGradient_BarthLimiter_Kernel_FlowPart(self_view, arg.portable, iCell, mesh.NumCell(), nVars, nVarsScalar);
104 }
105
106 // for scalars
107
108 if (nVarsScalar)
109#if defined(DNDS_DIST_MT_USE_OMP)
110# pragma omp parallel for schedule(runtime)
111#endif
112 for (index iCell = 0; iCell < mesh.NumCell(); iCell++)
113 {
114 RecGradient_BarthLimiter_Kernel_ScalarPart(self_view, arg.portable, iCell, mesh.NumCell(), nVars, nVarsScalar);
115 }
116 }
117
118 /**
119 * @brief Host specialization: conservative-to-primitive conversion with viscosity.
120 *
121 * OpenMP-parallel loop over all points in the u array, calling Cons2PrimMu_Kernel
122 * for each cell.
123 */
124 template <>
126 {
127 using namespace Geom;
128
132
133 int nVarsScalar = uScalar.size();
134 int nVars = nVarsFlow + nVarsScalar;
135
136 auto &mesh = self_view.fv.mesh;
137 auto &fv = self_view.fv;
138 auto &bcHandler = self_view.bcHandler;
139 auto &phy = self_view.physics;
140
141#if defined(DNDS_DIST_MT_USE_OMP)
142# pragma omp parallel for schedule(runtime)
143#endif
144 for (index iPt = 0; iPt < u.Size(); iPt++)
145 {
146 Cons2PrimMu_Kernel(self_view, arg.portable, iPt, u.Size(), nVars, nVarsScalar);
147 }
148 }
149
150 /**
151 * @brief Host specialization: conservative-to-primitive conversion (no gradients/viscosity).
152 *
153 * OpenMP-parallel loop over all points in the u array, calling Cons2Prim_Kernel
154 * for each cell.
155 */
156 template <>
158 {
159 using namespace Geom;
160
164
165 int nVarsScalar = uScalar.size();
166 int nVars = nVarsFlow + nVarsScalar;
167
168 auto &mesh = self_view.fv.mesh;
169 auto &fv = self_view.fv;
170 auto &bcHandler = self_view.bcHandler;
171 auto &phy = self_view.physics;
172
173#if defined(DNDS_DIST_MT_USE_OMP)
174# pragma omp parallel for schedule(runtime)
175#endif
176 for (index iPt = 0; iPt < u.Size(); iPt++)
177 {
178 Cons2Prim_Kernel(self_view, arg.portable, iPt, u.Size(), nVars, nVarsScalar);
179 }
180 }
181
182 /**
183 * @brief Host specialization: per-face eigenvalue estimation.
184 *
185 * OpenMP-parallel loop over all processor-local faces, calling
186 * EstEigenDt_GetFaceLam_Kernel for each face.
187 */
188 template <>
190 {
191 using namespace Geom;
192
195
202
203 auto &mesh = self_view.fv.mesh;
204 auto &fv = self_view.fv;
205 auto &bcHandler = self_view.bcHandler;
206 auto &phy = self_view.physics;
207
208 int nVarsScalar = 0; // no uScalar here
209 int nVars = nVarsFlow + nVarsScalar;
210
211#if defined(DNDS_DIST_MT_USE_OMP)
212# pragma omp parallel for schedule(runtime)
213#endif
214 for (index iFace = 0; iFace < mesh.NumFaceProc(); iFace++)
215 {
216 EstEigenDt_GetFaceLam_Kernel(self_view, arg.portable, iFace, mesh.NumFaceProc(), nVars, nVarsScalar);
217 }
218 }
219
220 /**
221 * @brief Host specialization: face eigenvalue accumulation to cell time steps.
222 *
223 * OpenMP-parallel loop over all owned cells, calling
224 * EstEigenDt_FaceLam2CellDt_Kernel for each cell.
225 */
226 template <>
228 {
229
232
233 auto &mesh = self_view.fv.mesh;
234 auto &fv = self_view.fv;
235 auto &bcHandler = self_view.bcHandler;
236 auto &phy = self_view.physics;
237
238 int nVarsScalar = 0; // no uScalar here
239 int nVars = nVarsFlow + nVarsScalar;
240#if defined(DNDS_DIST_MT_USE_OMP)
241# pragma omp parallel for schedule(runtime)
242#endif
243 for (index iCell = 0; iCell < mesh.NumCell(); iCell++)
244 {
245 EstEigenDt_FaceLam2CellDt_Kernel(self_view, arg.portable, iCell, mesh.NumCell(), nVars, nVarsScalar);
246 }
247 }
248
249 /**
250 * @brief Host specialization: 2nd-order face value reconstruction.
251 *
252 * OpenMP-parallel loop over all processor-local faces, calling
253 * RecFace2nd_Kernel for each face.
254 */
255 template <>
257 {
258 using namespace Geom;
261
263
264 auto &mesh = self_view.fv.mesh;
265 auto &fv = self_view.fv;
266 auto &bcHandler = self_view.bcHandler;
267 auto &phy = self_view.physics;
268
269 int nVarsScalar = uScalar.size();
270 int nVars = nVarsScalar + nVarsFlow;
271
272#if defined(DNDS_DIST_MT_USE_OMP)
273# pragma omp parallel for schedule(runtime)
274#endif
275 for (index iFace = 0; iFace < mesh.NumFaceProc(); iFace++)
276 {
277 RecFace2nd_Kernel(self_view, arg.portable, iFace, mesh.NumFaceProc(), nVars, nVarsScalar);
278 }
279 }
280
281 /**
282 * @brief Host specialization: 2nd-order Roe flux evaluation and face-to-cell RHS scatter.
283 *
284 * First loop (OpenMP parallel): computes per-face numerical flux via Flux2nd_Kernel_FluxFace.
285 * Second loop (OpenMP parallel): scatters face fluxes to cell RHS via Flux2nd_Kernel_Face2Cell.
286 */
287 template <>
289 {
290 using namespace Geom;
293
295
296 auto &mesh = self_view.fv.mesh;
297 auto &fv = self_view.fv;
298 auto &bcHandler = self_view.bcHandler;
299 auto &phy = self_view.physics;
300
301 int nVarsScalar = uScalarFL.size();
302 int nVars = nVarsScalar + nVarsFlow;
303
304#if defined(DNDS_DIST_MT_USE_OMP)
305# pragma omp parallel for schedule(runtime)
306#endif
307 for (index iFace = 0; iFace < mesh.NumFaceProc(); iFace++)
308 {
309 Flux2nd_Kernel_FluxFace(self_view, arg.portable, iFace, mesh.NumFaceProc(), nVars, nVarsScalar);
310 }
311
312#if defined(DNDS_DIST_MT_USE_OMP)
313# pragma omp parallel for schedule(runtime)
314#endif
315 for (index iCell = 0; iCell < mesh.NumCell(); iCell++)
316 {
317 Flux2nd_Kernel_Face2Cell(self_view, arg.portable, iCell, mesh.NumCell(), nVars, nVarsScalar);
318 }
319 }
320}
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
Device memory abstraction layer with backend-specific storage and factory creation.
#define DNDS_check_throw(expr)
Runtime check active in both debug and release builds. Throws std::runtime_error if expr evaluates to...
Definition Errors.hpp:89
Core type definitions and utilities for the EulerP alternative Navier-Stokes evaluator module.
Backend-specific implementation layer for EulerP Evaluator kernel dispatch.
#define DNDS_EULERP_IMPL_ARG_GET_REF_PORTABLE(member)
#define DNDS_EULERP_IMPL_ARG_GET_REF(member)
Shared device-callable kernel implementations for the EulerP Evaluator.
Namespace for the EulerP alternative evaluator module with GPU support.
Definition EulerP.hpp:29
DNDS_DEVICE void EstEigenDt_FaceLam2CellDt_Kernel(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::EstEigenDt_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
Per-cell kernel converting face eigenvalues to a local CFL time step.
DNDS_DEVICE void Flux2nd_Kernel_FluxFace(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Flux2nd_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
2nd-order Roe inviscid flux computation kernel (per-face).
DNDS_DEVICE void RecGradient_BarthLimiter_Kernel_FlowPart(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
Barth-Jespersen gradient limiter kernel for flow variables (per-cell).
DNDS_DEVICE void RecFace2nd_Kernel(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecFace2nd_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
2nd-order face value reconstruction kernel (per-face).
DNDS_DEVICE void RecGradient_GGRec_Kernel_GG(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
Green-Gauss gradient reconstruction kernel (per-cell).
DNDS_DEVICE void Cons2Prim_Kernel(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Cons2Prim_Arg::Portable &arg, index iPt, index iPtEnd, int nVars, int nVarsScalar)
Conservative-to-primitive conversion without gradient transformation or viscosity (per-cell).
DNDS_DEVICE void RecGradient_BarthLimiter_Kernel_ScalarPart(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
Barth-Jespersen gradient limiter kernel for transported scalar variables (per-cell).
DNDS_DEVICE void Flux2nd_Kernel_Face2Cell(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Flux2nd_Arg::Portable &arg, index iCell, index iCellEnd, int nVars, int nVarsScalar)
Scatters face fluxes to cell RHS residual (per-cell).
DNDS_DEVICE void RecGradient_GGRec_Kernel_BndVal(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::RecGradient_Arg::Portable &arg, index iBnd, index iBndEnd, int nVars, int nVarsScalar)
Generates boundary ghost values for Green-Gauss gradient reconstruction.
DNDS_DEVICE void EstEigenDt_GetFaceLam_Kernel(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::EstEigenDt_Arg::Portable &arg, index iFace, index iFaceEnd, int nVars, int nVarsScalar)
Per-face eigenvalue estimation kernel for time-step computation.
DNDS_DEVICE void Cons2PrimMu_Kernel(EvaluatorDeviceView< B > &self_view, typename Evaluator_impl< B >::Cons2PrimMu_Arg::Portable &arg, index iPt, index iPtEnd, int nVars, int nVarsScalar)
Conservative-to-primitive conversion with gradient transformation and viscosity (per-cell).
DeviceBackend
Enumerates the backends a DeviceStorage / Array can live on.
@ Host
Plain CPU memory.
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:107
Device-side argument struct for conservative-to-primitive + viscosity kernel.
struct DNDS::EulerP::Evaluator_impl::Cons2PrimMu_Arg::Portable portable
Device-side argument struct for conservative-to-primitive conversion (no gradients/viscosity).
struct DNDS::EulerP::Evaluator_impl::Cons2Prim_Arg::Portable portable
Device-side argument struct for eigenvalue estimation and time-step computation.
struct DNDS::EulerP::Evaluator_impl::EstEigenDt_Arg::Portable portable
Device-side argument struct for 2nd-order flux evaluation and RHS accumulation.
struct DNDS::EulerP::Evaluator_impl::Flux2nd_Arg::Portable portable
Device-side argument struct for 2nd-order face value reconstruction.
struct DNDS::EulerP::Evaluator_impl::RecFace2nd_Arg::Portable portable
Device-side argument struct for gradient reconstruction kernels.
struct DNDS::EulerP::Evaluator_impl::RecGradient_Arg::Portable portable
static void Flux2nd(Flux2nd_Arg &arg)
Evaluates 2nd-order Roe flux per face and scatters to cell RHS.
static void Cons2PrimMu(Cons2PrimMu_Arg &arg)
Executes conservative-to-primitive conversion with viscosity computation.
static void EstEigenDt_GetFaceLam(EstEigenDt_Arg &arg)
First pass: computes per-face eigenvalue estimates from cell states.
static void Cons2Prim(Cons2Prim_Arg &arg)
Executes conservative-to-primitive conversion (no gradients/viscosity).
static void EstEigenDt_FaceLam2CellDt(EstEigenDt_Arg &arg)
Second pass: accumulates face eigenvalues to cells and computes local dt.
static void RecFace2nd(RecFace2nd_Arg &arg)
Executes 2nd-order face value reconstruction from cell-centered data.
static void RecGradient_BarthLimiter(RecGradient_Arg &arg)
Barth-Jespersen gradient limiter applied to reconstructed gradients.
static void RecGradient_GGRec(RecGradient_Arg &arg)
Green-Gauss gradient reconstruction: boundary ghost values + cell gradient computation.