PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Metric.h
Go to the documentation of this file.
1#ifndef METRIC_H
2#define METRIC_H
3
4// Include necessary headers
5#include <petsc.h>
6#include "variables.h" // Common type definitions
7#include "logging.h" // Logging macros and definitions
8#include <stdlib.h>
9#include "io.h"
10#include "setup.h" // For SetDMDAProcLayout
11
12/**
13 * @brief Public interface for `MetricLogicalToPhysical()`.
14 *
15 * @param user Primary `UserCtx` input for the operation.
16 * @param X Parameter `X` passed to `MetricLogicalToPhysical()`.
17 * @param i Parameter `i` passed to `MetricLogicalToPhysical()`.
18 * @param j Parameter `j` passed to `MetricLogicalToPhysical()`.
19 * @param k Parameter `k` passed to `MetricLogicalToPhysical()`.
20 * @param xi Parameter `xi` passed to `MetricLogicalToPhysical()`.
21 * @param eta Parameter `eta` passed to `MetricLogicalToPhysical()`.
22 * @param zta Parameter `zta` passed to `MetricLogicalToPhysical()`.
23 * @param Xp Parameter `Xp` passed to `MetricLogicalToPhysical()`.
24 * @return PetscErrorCode 0 on success.
25 */
26PetscErrorCode MetricLogicalToPhysical(UserCtx *user, const Cmpnts ***X,
27 PetscInt i,PetscInt j,PetscInt k,
28 PetscReal xi,PetscReal eta,PetscReal zta,
29 Cmpnts *Xp);
30
31/**
32 * @brief Public interface for `MetricGetCellVertices()`.
33 *
34 * @param user Primary `UserCtx` input for the operation.
35 * @param X Parameter `X` passed to `MetricGetCellVertices()`.
36 * @param i Parameter `i` passed to `MetricGetCellVertices()`.
37 * @param j Parameter `j` passed to `MetricGetCellVertices()`.
38 * @param k Parameter `k` passed to `MetricGetCellVertices()`.
39 * @param V Parameter `V` passed to `MetricGetCellVertices()`.
40 * @return PetscErrorCode 0 on success.
41 */
42PetscErrorCode MetricGetCellVertices(UserCtx *user,const Cmpnts ***X,
43 PetscInt i,PetscInt j,PetscInt k,
44 Cmpnts V[8]);
45
46/**
47 * @brief Public interface for `MetricJacobian()`.
48 *
49 * @param user Primary `UserCtx` input for the operation.
50 * @param X Parameter `X` passed to `MetricJacobian()`.
51 * @param i Parameter `i` passed to `MetricJacobian()`.
52 * @param j Parameter `j` passed to `MetricJacobian()`.
53 * @param k Parameter `k` passed to `MetricJacobian()`.
54 * @param xi Parameter `xi` passed to `MetricJacobian()`.
55 * @param eta Parameter `eta` passed to `MetricJacobian()`.
56 * @param zta Parameter `zta` passed to `MetricJacobian()`.
57 * @param J Parameter `J` passed to `MetricJacobian()`.
58 * @param detJ Parameter `detJ` passed to `MetricJacobian()`.
59 * @return PetscErrorCode 0 on success.
60 */
61PetscErrorCode MetricJacobian(UserCtx *user,const Cmpnts ***X,
62 PetscInt i,PetscInt j,PetscInt k,
63 PetscReal xi,PetscReal eta,PetscReal zta,
64 PetscReal J[3][3],PetscReal *detJ);
65
66/**
67 * @brief Public interface for `MetricVelocityContravariant()`.
68 *
69 * @param J Parameter `J` passed to `MetricVelocityContravariant()`.
70 * @param detJ Parameter `detJ` passed to `MetricVelocityContravariant()`.
71 * @param u Parameter `u` passed to `MetricVelocityContravariant()`.
72 * @param uc Parameter `uc` passed to `MetricVelocityContravariant()`.
73 * @return PetscErrorCode 0 on success.
74 */
75PetscErrorCode MetricVelocityContravariant(const PetscReal J[3][3],
76 PetscReal detJ,
77 const PetscReal u[3],PetscReal uc[3]);
78
79/**
80 * @brief Computes the unit normal vectors and areas of the three faces of a computational cell.
81 *
82 * Given the metric vectors (csi, eta, zet), this function calculates the geometric
83 * properties of the cell faces aligned with the i, j, and k directions.
84 *
85 * @param csi Parameter `csi` passed to `CalculateFaceNormalAndArea()`.
86 * @param eta Parameter `eta` passed to `CalculateFaceNormalAndArea()`.
87 * @param zet Parameter `zet` passed to `CalculateFaceNormalAndArea()`.
88 * @param ni Output:
89 * @param nj Output:
90 * @param nk Output:
91 * @param Ai Output:
92 * @param Aj Output:
93 * @param Ak Output:
94 * @return PetscErrorCode 0 on success.
95 */
96PetscErrorCode CalculateFaceNormalAndArea(Cmpnts csi, Cmpnts eta, Cmpnts zet, double ni[3], double nj[3], double nk[3], double *Ai, double *Aj, double *Ak);
97
98/**
99 * @brief Inverts the 3x3 covariant metric tensor to obtain the contravariant metric tensor.
100 *
101 * In curvilinear coordinates, the input matrix `g` contains the dot products of the
102 * covariant basis vectors (e.g., g_ij = e_i . e_j). Its inverse, `G`, is the
103 * contravariant metric tensor, which is essential for transforming vectors and tensors
104 * between coordinate systems.
105 *
106 * @param covariantTensor Input: A 3x3 matrix representing the covariant metric tensor.
107 * @param contravariantTensor Output: A 3x3 matrix where the inverted result is stored.
108 * @return PetscErrorCode 0 on success.
109 */
110PetscErrorCode InvertCovariantMetricTensor(double covariantTensor[3][3], double contravariantTensor[3][3]);
111
112/**
113 * @brief Computes characteristic length scales (dx, dy, dz) for a curvilinear cell.
114 *
115 * For a non-uniform, non-orthogonal cell, there is no single "dx". This function
116 * computes an effective length scale in each Cartesian direction based on the cell
117 * volume and the areas of its faces.
118 *
119 * @param ajc The
120 * @param csi Parameter `csi` passed to `ComputeCellCharacteristicLengthScale()`.
121 * @param eta Parameter `eta` passed to `ComputeCellCharacteristicLengthScale()`.
122 * @param zet Parameter `zet` passed to `ComputeCellCharacteristicLengthScale()`.
123 * @param dx Output:
124 * @param dy Output:
125 * @param dz Output:
126 * @return PetscErrorCode 0 on success.
127 */
128PetscErrorCode ComputeCellCharacteristicLengthScale(PetscReal ajc, Cmpnts csi, Cmpnts eta, Cmpnts zet, double *dx, double *dy, double *dz);
129
130/**
131 * @brief Applies periodic boundary corrections to cell centers (Cent) and grid spacing (GridSpace).
132 *
133 * This function handles the special logic needed when periodic boundaries are present.
134 * For coarse grids (cgrid), it directly copies from ghost regions. For fine grids,
135 * it calculates the boundary values using the GridSpace information.
136 *
137 * @param user The UserCtx containing grid and field data.
138 * @return PetscErrorCode 0 on success.
139 */
141
142/**
143 * @brief Applies periodic boundary corrections to i-face centers (Centx).
144 *
145 * Only X-direction periodicity affects Centx. This function must be called
146 * after Centx has been initially computed but before it's used for metric calculations.
147 *
148 * @param user The UserCtx containing grid and field data.
149 * @return PetscErrorCode 0 on success.
150 */
152
153/**
154 * @brief Applies periodic boundary corrections to j-face centers (Centy).
155 *
156 * Only Y-direction periodicity affects Centy. This function must be called
157 * after Centy has been initially computed but before it's used for metric calculations.
158 *
159 * @param user The UserCtx containing grid and field data.
160 * @return PetscErrorCode 0 on success.
161 */
163
164/**
165 * @brief Applies periodic boundary corrections to k-face centers (Centz).
166 *
167 * Only Z-direction periodicity affects Centz. This function must be called
168 * after Centz has been initially computed but before it's used for metric calculations.
169 *
170 * @param user The UserCtx containing grid and field data.
171 * @return PetscErrorCode 0 on success.
172 */
174
175/**
176 * @brief Computes the primary face metric components (Csi, Eta, Zet), including
177 * boundary extrapolation, and stores them in the corresponding global Vec
178 * members of the UserCtx structure (user->Csi, user->Eta, user->Zet).
179 *
180 * This is a self-contained routine that performs the following steps:
181 * 1. Obtains local ghosted nodal coordinates using DMGetCoordinatesLocal.
182 * 2. Calculates metrics for INTERIOR faces where finite difference stencils are valid.
183 * 3. EXTRAPOLATES metrics for faces on the physical domain boundaries by copying
184 * from the nearest computed interior face.
185 * 4. Assembles the global `user->Csi`, `user->Eta`, `user->Zet` Vecs.
186 * 5. Updates the local ghosted `user->lCsi`, `user->lEta`, `user->lZet` Vecs.
187 *
188 * @param[in,out] user Pointer to the UserCtx structure.
189 *
190 * @return PetscErrorCode 0 on success.
191 *
192 * @note
193 * - This function is a complete "compute and make ready" unit for Csi, Eta, and Zet.
194 * - It's recommended to call `VecZeroEntries` on user->Csi, Eta, Zet before this
195 * if they might contain old data.
196 */
197PetscErrorCode ComputeFaceMetrics(UserCtx *user);
198
199/**
200 * @brief Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells
201 * and stores it in `user->Aj`. This version includes boundary extrapolation.
202 *
203 * Nodal coordinates are obtained internally.
204 * Refer to previous Doxygen comments for details on physical locations and
205 * storage convention (`aj_arr[k_n][j_n][i_n]` for cell `C(i_n-1,j_n-1,k_n-1)`).
206 *
207 * @param[in,out] user Pointer to the UserCtx structure.
208 *
209 * @return PetscErrorCode 0 on success.
210 */
211PetscErrorCode ComputeCellCenteredJacobianInverse(UserCtx *user);
212
213/**
214 * @brief Ensure a **right-handed** metric basis (`Csi`, `Eta`, `Zet`) and a
215 * **positive Jacobian** (`Aj`) over the whole domain.
216 *
217 * The metric-generation kernels are completely algebraic, so they will happily
218 * deliver a *left-handed* basis if the mesh file enumerates nodes in the
219 * opposite ζ-direction.
220 * This routine makes the orientation explicit and—if needed—repairs it
221 * **once per run**:
222 *
223 * | Step | Action |
224 * |------|--------|
225 * | 1 | Compute global `Aj_min`, `Aj_max`. |
226 * | 2 | **Mixed signs** (`Aj_min < 0 && Aj_max > 0`) &rarr; abort: the mesh is topologically inconsistent. |
227 * | 3 | **All negative** (`Aj_max < 0`) &rarr; flip <br>`Csi`, `Eta`, `Zet`, `Aj` & update local ghosts. |
228 * | 4 | Store `user->orientation = ±1` so BC / IC routines can apply sign-aware logic if they care about inlet direction. |
229 *
230 * @param[in,out] user Fully initialised #UserCtx that already contains
231 * `Csi`, `Eta`, `Zet`, `Aj`, their **local** ghosts, and
232 * valid distributed DMs.
233 *
234 * @return `0` on success or a PETSc error code on failure.
235 *
236 * @note Call **immediately after** `ComputeCellCenteredJacobianInverse()` and
237 * before any routine that differentiates or applies BCs.
238 *
239 * @note Author metadata intentionally omitted in API docs.
240 */
241PetscErrorCode CheckAndFixGridOrientation(UserCtx *user);
242
243/**
244 * @brief Computes the physical location of cell centers and the spacing between them.
245 *
246 * This function calculates two key geometric properties from the nodal coordinates:
247 * 1. `Cent`: A vector field storing the (x,y,z) coordinates of the center of each grid cell.
248 * 2. `GridSpace`: A vector field storing the physical distance between adjacent
249 * cell centers in the i, j, and k computational directions.
250 *
251 * It is a direct adaptation of the corresponding logic from the legacy `FormMetrics`.
252 *
253 * @param user The UserCtx for a specific grid level. The function populates `user->Cent` and `user->GridSpace`.
254 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
255 */
256PetscErrorCode ComputeCellCentersAndSpacing(UserCtx *user);
257
258/**
259 * @brief Computes metrics centered on constant-i faces (i-faces).
260 *
261 * This function calculates the metric terms (ICsi, IEta, IZet) and the inverse
262 * Jacobian (IAj) located at the center of each i-face. The stencils use
263 * i-face-centered coordinates (`Centx`) which must be computed first.
264 * The logic is a direct adaptation of the legacy FormMetrics function.
265 *
266 * @param user The UserCtx for a specific grid level. Populates user->ICsi, etc.
267 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
268 */
269PetscErrorCode ComputeIFaceMetrics(UserCtx *user);
270
271/**
272 * @brief Computes metrics centered on constant-j faces (j-faces).
273 *
274 * This function calculates the metric terms (`JCsi`, `JEta`, `JZet`) and the
275 * inverse Jacobian (`JAj`) located at the geometric center of each constant-j
276 * face. This is a critical step for staggered-grid finite difference schemes.
277 *
278 * The process is a direct and faithful refactoring of the corresponding logic
279 * from the legacy `FormMetrics` function:
280 * 1. It first calculates the physical (x,y,z) coordinates of the center of
281 * each i-face and stores them in the `user->Centy` vector.
282 * 2. It then uses a boundary-aware, second-order finite difference stencil on
283 * the `Centy` field to compute the derivatives (e.g., d(x)/d(csi)).
284 * - Central differences are used in the grid interior.
285 * - One-sided differences are used at the physical domain boundaries.
286 * 3. Finally, these derivatives are used to compute the final metric terms and
287 * the inverse Jacobian, which are stored in their respective `Vec` objects.
288 *
289 * @param user The UserCtx for a specific grid level. This function populates
290 * the `user->JCsi`, `user->JEta`, `user->JZet`, and `user->JAj` vectors.
291 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
292 */
293PetscErrorCode ComputeJFaceMetrics(UserCtx *user);
294
295/**
296 * @brief Computes metrics centered on constant-k faces (k-faces).
297 *
298 * This function calculates the metric terms (`KCsi`, `KEta`, `KZet`) and the
299 * inverse Jacobian (`KAj`) located at the geometric center of each constant-j
300 * face. This is a critical step for staggered-grid finite difference schemes.
301 *
302 * The process is a direct and faithful refactoring of the corresponding logic
303 * from the legacy `FormMetrics` function:
304 * 1. It first calculates the physical (x,y,z) coordinates of the center of
305 * each i-face and stores them in the `user->Centz` vector.
306 * 2. It then uses a boundary-aware, second-order finite difference stencil on
307 * the `Centz` field to compute the derivatives (e.g., d(x)/d(csi)).
308 * - Central differences are used in the grid interior.
309 * - One-sided differences are used at the physical domain boundaries.
310 * 3. Finally, these derivatives are used to compute the final metric terms and
311 * the inverse Jacobian, which are stored in their respective `Vec` objects.
312 *
313 * @param user The UserCtx for a specific grid level. This function populates
314 * the `user->KCsi`, `user->KEta`, `user->KZet`, and `user->KAj` vectors.
315 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
316 */
317PetscErrorCode ComputeKFaceMetrics(UserCtx *user);
318
319/**
320 * @brief Performs a diagnostic check on the divergence of the face area metric vectors.
321 *
322 * For a closed cell, the sum of the face area vectors should be zero (Gauss's
323 * divergence theorem). This function computes a measure of this divergence and
324 * reports the maximum value over the domain. A small value indicates a
325 * well-formed grid. This is a direct adaptation of the legacy function.
326 *
327 * @param user The UserCtx for a specific grid level (typically the finest).
328 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
329 */
330PetscErrorCode ComputeMetricsDivergence(UserCtx *user);
331
332/**
333 * @brief Computes the max-min values of the grid metrics.
334 *
335 * This function serves as a diagnostic tool to assess the quality of the grid
336 * metrics. It calculates the bounds of the face metrics (Csi, Eta, Zet).
337 *
338 * @param user The UserCtx, containing all necessary grid data.
339 * @return PetscErrorCode
340 */
341PetscErrorCode ComputeMetricNorms(UserCtx *user);
342
343/**
344 * @brief Orchestrates the calculation of all grid metrics.
345 *
346 * This function iterates through every UserCtx in the multigrid and multi-block
347 * hierarchy. For each context, it calls a series of modern, modular helper
348 * functions to compute the face metrics (Csi, Eta, Zet), the cell-centered
349 * inverse Jacobian (Aj), and to validate the grid's orientation.
350 *
351 * @param simCtx The master SimCtx, containing the configured UserCtx hierarchy.
352 * @return PetscErrorCode
353 */
354PetscErrorCode CalculateAllGridMetrics(SimCtx *simCtx);
355#endif /* METRIC_H */
PetscErrorCode CalculateAllGridMetrics(SimCtx *simCtx)
Orchestrates the calculation of all grid metrics.
Definition Metric.c:2168
PetscErrorCode MetricJacobian(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, PetscReal xi, PetscReal eta, PetscReal zta, PetscReal J[3][3], PetscReal *detJ)
Public interface for MetricJacobian().
Definition Metric.c:106
PetscErrorCode InvertCovariantMetricTensor(double covariantTensor[3][3], double contravariantTensor[3][3])
Inverts the 3x3 covariant metric tensor to obtain the contravariant metric tensor.
Definition Metric.c:203
PetscErrorCode ComputeMetricNorms(UserCtx *user)
Computes the max-min values of the grid metrics.
Definition Metric.c:2036
PetscErrorCode ApplyPeriodicCorrectionsToCellCentersAndSpacing(UserCtx *user)
Applies periodic boundary corrections to cell centers (Cent) and grid spacing (GridSpace).
Definition Metric.c:392
PetscErrorCode CheckAndFixGridOrientation(UserCtx *user)
Ensure a right-handed metric basis (Csi, Eta, Zet) and a positive Jacobian (Aj) over the whole domain...
Definition Metric.c:314
PetscErrorCode ComputeCellCentersAndSpacing(UserCtx *user)
Computes the physical location of cell centers and the spacing between them.
Definition Metric.c:1186
PetscErrorCode MetricGetCellVertices(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, Cmpnts V[8])
Public interface for MetricGetCellVertices().
Definition Metric.c:26
PetscErrorCode ComputeJFaceMetrics(UserCtx *user)
Computes metrics centered on constant-j faces (j-faces).
Definition Metric.c:1511
PetscErrorCode ComputeMetricsDivergence(UserCtx *user)
Performs a diagnostic check on the divergence of the face area metric vectors.
Definition Metric.c:1944
PetscErrorCode ComputeFaceMetrics(UserCtx *user)
Computes the primary face metric components (Csi, Eta, Zet), including boundary extrapolation,...
Definition Metric.c:830
PetscErrorCode ComputeCellCharacteristicLengthScale(PetscReal ajc, Cmpnts csi, Cmpnts eta, Cmpnts zet, double *dx, double *dy, double *dz)
Computes characteristic length scales (dx, dy, dz) for a curvilinear cell.
Definition Metric.c:283
PetscErrorCode CalculateFaceNormalAndArea(Cmpnts csi, Cmpnts eta, Cmpnts zet, double ni[3], double nj[3], double nk[3], double *Ai, double *Aj, double *Ak)
Computes the unit normal vectors and areas of the three faces of a computational cell.
Definition Metric.c:237
PetscErrorCode ApplyPeriodicCorrectionsToIFaceCenter(UserCtx *user)
Applies periodic boundary corrections to i-face centers (Centx).
Definition Metric.c:591
PetscErrorCode ComputeCellCenteredJacobianInverse(UserCtx *user)
Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells and stores it in u...
Definition Metric.c:1042
PetscErrorCode ComputeKFaceMetrics(UserCtx *user)
Computes metrics centered on constant-k faces (k-faces).
Definition Metric.c:1719
PetscErrorCode ApplyPeriodicCorrectionsToJFaceCenter(UserCtx *user)
Applies periodic boundary corrections to j-face centers (Centy).
Definition Metric.c:671
PetscErrorCode MetricLogicalToPhysical(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, PetscReal xi, PetscReal eta, PetscReal zta, Cmpnts *Xp)
Public interface for MetricLogicalToPhysical().
Definition Metric.c:77
PetscErrorCode ComputeIFaceMetrics(UserCtx *user)
Computes metrics centered on constant-i faces (i-faces).
Definition Metric.c:1289
PetscErrorCode ApplyPeriodicCorrectionsToKFaceCenter(UserCtx *user)
Applies periodic boundary corrections to k-face centers (Centz).
Definition Metric.c:751
PetscErrorCode MetricVelocityContravariant(const PetscReal J[3][3], PetscReal detJ, const PetscReal u[3], PetscReal uc[3])
Public interface for MetricVelocityContravariant().
Definition Metric.c:166
Public interface for data input/output routines.
Logging utilities and macros for PETSc-based applications.
Main header file for a complex fluid dynamics solver.
A 3D point or vector with PetscScalar components.
Definition variables.h:100
The master context for the entire simulation.
Definition variables.h:643
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811