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
12PetscErrorCode MetricLogicalToPhysical(UserCtx *user, const Cmpnts ***X,
13 PetscInt i,PetscInt j,PetscInt k,
14 PetscReal xi,PetscReal eta,PetscReal zta,
15 Cmpnts *Xp);
16
17PetscErrorCode MetricGetCellVertices(UserCtx *user,const Cmpnts ***X,
18 PetscInt i,PetscInt j,PetscInt k,
19 Cmpnts V[8]);
20
21PetscErrorCode MetricJacobian(UserCtx *user,const Cmpnts ***X,
22 PetscInt i,PetscInt j,PetscInt k,
23 PetscReal xi,PetscReal eta,PetscReal zta,
24 PetscReal J[3][3],PetscReal *detJ);
25
26PetscErrorCode MetricVelocityContravariant(const PetscReal J[3][3],
27 PetscReal detJ,
28 const PetscReal u[3],PetscReal uc[3]);
29
30/**
31 * @brief Computes the unit normal vectors and areas of the three faces of a computational cell.
32 *
33 * Given the metric vectors (csi, eta, zet), this function calculates the geometric
34 * properties of the cell faces aligned with the i, j, and k directions.
35 *
36 * @param csi, eta, zet The metric vectors at the cell center.
37 * @param ni Output: A 3-element array for the unit normal vector of the i-face.
38 * @param nj Output: A 3-element array for the unit normal vector of the j-face.
39 * @param nk Output: A 3-element array for the unit normal vector of the k-face.
40 * @param Ai Output: Pointer to store the area of the i-face.
41 * @param Aj Output: Pointer to store the area of the j-face.
42 * @param Ak Output: Pointer to store the area of the k-face.
43 * @return PetscErrorCode 0 on success.
44 */
45PetscErrorCode CalculateFaceNormalAndArea(Cmpnts csi, Cmpnts eta, Cmpnts zet, double ni[3], double nj[3], double nk[3], double *Ai, double *Aj, double *Ak);
46
47/**
48 * @brief Inverts the 3x3 covariant metric tensor to obtain the contravariant metric tensor.
49 *
50 * In curvilinear coordinates, the input matrix `g` contains the dot products of the
51 * covariant basis vectors (e.g., g_ij = e_i . e_j). Its inverse, `G`, is the
52 * contravariant metric tensor, which is essential for transforming vectors and tensors
53 * between coordinate systems.
54 *
55 * @param covariantTensor Input: A 3x3 matrix representing the covariant metric tensor.
56 * @param contravariantTensor Output: A 3x3 matrix where the inverted result is stored.
57 * @return PetscErrorCode 0 on success.
58 */
59PetscErrorCode InvertCovariantMetricTensor(double covariantTensor[3][3], double contravariantTensor[3][3]);
60
61/**
62 * @brief Computes characteristic length scales (dx, dy, dz) for a curvilinear cell.
63 *
64 * For a non-uniform, non-orthogonal cell, there is no single "dx". This function
65 * computes an effective length scale in each Cartesian direction based on the cell
66 * volume and the areas of its faces.
67 *
68 * @param ajc The Jacobian of the grid transformation (1 / cell volume).
69 * @param csi, eta, zet The metric vectors at the cell center.
70 * @param dx Output: Pointer to store the characteristic length in the x-direction.
71 * @param dy Output: Pointer to store the characteristic length in the y-direction.
72 * @param dz Output: Pointer to store the characteristic length in the z-direction.
73 * @return PetscErrorCode 0 on success.
74 */
75PetscErrorCode ComputeCellCharacteristicLengthScale(PetscReal ajc, Cmpnts csi, Cmpnts eta, Cmpnts zet, double *dx, double *dy, double *dz);
76
77
78/**
79 * @brief Computes the primary face metric components (Csi, Eta, Zet), including
80 * boundary extrapolation, and stores them in the corresponding global Vec
81 * members of the UserCtx structure (user->Csi, user->Eta, user->Zet).
82 *
83 * This is a self-contained routine that performs the following steps:
84 * 1. Obtains local ghosted nodal coordinates using DMGetCoordinatesLocal.
85 * 2. Calculates metrics for INTERIOR faces where finite difference stencils are valid.
86 * 3. EXTRAPOLATES metrics for faces on the physical domain boundaries by copying
87 * from the nearest computed interior face.
88 * 4. Assembles the global `user->Csi`, `user->Eta`, `user->Zet` Vecs.
89 * 5. Updates the local ghosted `user->lCsi`, `user->lEta`, `user->lZet` Vecs.
90 *
91 * @param[in,out] user Pointer to the UserCtx structure.
92 *
93 * @return PetscErrorCode 0 on success.
94 *
95 * @note
96 * - This function is a complete "compute and make ready" unit for Csi, Eta, and Zet.
97 * - It's recommended to call `VecZeroEntries` on user->Csi, Eta, Zet before this
98 * if they might contain old data.
99 */
100PetscErrorCode ComputeFaceMetrics(UserCtx *user);
101
102/**
103 * @brief Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells
104 * and stores it in `user->Aj`. This version includes boundary extrapolation.
105 *
106 * Nodal coordinates are obtained internally.
107 * Refer to previous Doxygen comments for details on physical locations and
108 * storage convention (`aj_arr[k_n][j_n][i_n]` for cell `C(i_n-1,j_n-1,k_n-1)`).
109 *
110 * @param[in,out] user Pointer to the UserCtx structure.
111 *
112 * @return PetscErrorCode 0 on success.
113 */
114PetscErrorCode ComputeCellCenteredJacobianInverse(UserCtx *user);
115
116/**
117 * @brief Ensure a **right-handed** metric basis (`Csi`, `Eta`, `Zet`) and a
118 * **positive Jacobian** (`Aj`) over the whole domain.
119 *
120 * The metric-generation kernels are completely algebraic, so they will happily
121 * deliver a *left-handed* basis if the mesh file enumerates nodes in the
122 * opposite ζ-direction.
123 * This routine makes the orientation explicit and—if needed—repairs it
124 * **once per run**:
125 *
126 * | Step | Action |
127 * |------|--------|
128 * | 1 | Compute global `Aj_min`, `Aj_max`. |
129 * | 2 | **Mixed signs** (`Aj_min < 0 && Aj_max > 0`) &rarr; abort: the mesh is topologically inconsistent. |
130 * | 3 | **All negative** (`Aj_max < 0`) &rarr; flip <br>`Csi`, `Eta`, `Zet`, `Aj` & update local ghosts. |
131 * | 4 | Store `user->orientation = ±1` so BC / IC routines can apply sign-aware logic if they care about inlet direction. |
132 *
133 * @param[in,out] user Fully initialised #UserCtx that already contains
134 * `Csi`, `Eta`, `Zet`, `Aj`, their **local** ghosts, and
135 * valid distributed DMs.
136 *
137 * @return `0` on success or a PETSc error code on failure.
138 *
139 * @note Call **immediately after** `ComputeCellCenteredJacobianInverse()` and
140 * before any routine that differentiates or applies BCs.
141 *
142 * @author <your name>
143 */
144PetscErrorCode CheckAndFixGridOrientation(UserCtx *user);
145
146/**
147 * @brief Computes the physical location of cell centers and the spacing between them.
148 *
149 * This function calculates two key geometric properties from the nodal coordinates:
150 * 1. `Cent`: A vector field storing the (x,y,z) coordinates of the center of each grid cell.
151 * 2. `GridSpace`: A vector field storing the physical distance between adjacent
152 * cell centers in the i, j, and k computational directions.
153 *
154 * It is a direct adaptation of the corresponding logic from the legacy `FormMetrics`.
155 *
156 * @param user The UserCtx for a specific grid level. The function populates `user->Cent` and `user->GridSpace`.
157 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
158 */
159PetscErrorCode ComputeCellCentersAndSpacing(UserCtx *user);
160
161/**
162 * @brief Computes metrics centered on constant-i faces (i-faces).
163 *
164 * This function calculates the metric terms (ICsi, IEta, IZet) and the inverse
165 * Jacobian (IAj) located at the center of each i-face. The stencils use
166 * i-face-centered coordinates (`Centx`) which must be computed first.
167 * The logic is a direct adaptation of the legacy FormMetrics function.
168 *
169 * @param user The UserCtx for a specific grid level. Populates user->ICsi, etc.
170 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
171 */
172PetscErrorCode ComputeIFaceMetrics(UserCtx *user);
173
174/**
175 * @brief Computes metrics centered on constant-j faces (j-faces).
176 *
177 * This function calculates the metric terms (`JCsi`, `JEta`, `JZet`) and the
178 * inverse Jacobian (`JAj`) located at the geometric center of each constant-j
179 * face. This is a critical step for staggered-grid finite difference schemes.
180 *
181 * The process is a direct and faithful refactoring of the corresponding logic
182 * from the legacy `FormMetrics` function:
183 * 1. It first calculates the physical (x,y,z) coordinates of the center of
184 * each i-face and stores them in the `user->Centy` vector.
185 * 2. It then uses a boundary-aware, second-order finite difference stencil on
186 * the `Centy` field to compute the derivatives (e.g., d(x)/d(csi)).
187 * - Central differences are used in the grid interior.
188 * - One-sided differences are used at the physical domain boundaries.
189 * 3. Finally, these derivatives are used to compute the final metric terms and
190 * the inverse Jacobian, which are stored in their respective `Vec` objects.
191 *
192 * @param user The UserCtx for a specific grid level. This function populates
193 * the `user->JCsi`, `user->JEta`, `user->JZet`, and `user->JAj` vectors.
194 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
195 */
196PetscErrorCode ComputeJFaceMetrics(UserCtx *user);
197
198/**
199 * @brief Computes metrics centered on constant-k faces (k-faces).
200 *
201 * This function calculates the metric terms (`KCsi`, `KEta`, `KZet`) and the
202 * inverse Jacobian (`KAj`) located at the geometric center of each constant-j
203 * face. This is a critical step for staggered-grid finite difference schemes.
204 *
205 * The process is a direct and faithful refactoring of the corresponding logic
206 * from the legacy `FormMetrics` function:
207 * 1. It first calculates the physical (x,y,z) coordinates of the center of
208 * each i-face and stores them in the `user->Centz` vector.
209 * 2. It then uses a boundary-aware, second-order finite difference stencil on
210 * the `Centz` field to compute the derivatives (e.g., d(x)/d(csi)).
211 * - Central differences are used in the grid interior.
212 * - One-sided differences are used at the physical domain boundaries.
213 * 3. Finally, these derivatives are used to compute the final metric terms and
214 * the inverse Jacobian, which are stored in their respective `Vec` objects.
215 *
216 * @param user The UserCtx for a specific grid level. This function populates
217 * the `user->KCsi`, `user->KEta`, `user->KZet`, and `user->KAj` vectors.
218 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
219 */
220PetscErrorCode ComputeKFaceMetrics(UserCtx *user);
221
222/**
223 * @brief Performs a diagnostic check on the divergence of the face area metric vectors.
224 *
225 * For a closed cell, the sum of the face area vectors should be zero (Gauss's
226 * divergence theorem). This function computes a measure of this divergence and
227 * reports the maximum value over the domain. A small value indicates a
228 * well-formed grid. This is a direct adaptation of the legacy function.
229 *
230 * @param user The UserCtx for a specific grid level (typically the finest).
231 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
232 */
233PetscErrorCode ComputeMetricsDivergence(UserCtx *user);
234
235/**
236 * @brief Computes the max-min values of the grid metrics.
237 *
238 * This function serves as a diagnostic tool to assess the quality of the grid
239 * metrics. It calculates the bounds of the face metrics (Csi, Eta, Zet).
240 *
241 * @param user The UserCtx, containing all necessary grid data.
242 * @return PetscErrorCode
243 */
244PetscErrorCode ComputeMetricNorms(UserCtx *user);
245
246/**
247 * @brief Orchestrates the calculation of all grid metrics.
248 *
249 * This function iterates through every UserCtx in the multigrid and multi-block
250 * hierarchy. For each context, it calls a series of modern, modular helper
251 * functions to compute the face metrics (Csi, Eta, Zet), the cell-centered
252 * inverse Jacobian (Aj), and to validate the grid's orientation.
253 *
254 * @param simCtx The master SimCtx, containing the configured UserCtx hierarchy.
255 * @return PetscErrorCode
256 */
257PetscErrorCode CalculateAllGridMetrics(SimCtx *simCtx);
258#endif /* METRIC_H */
PetscErrorCode CalculateAllGridMetrics(SimCtx *simCtx)
Orchestrates the calculation of all grid metrics.
Definition Metric.c:1809
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)
Compute Jacobian matrix and its determinant at (xi,eta,zta).
Definition Metric.c:110
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:212
PetscErrorCode ComputeMetricNorms(UserCtx *user)
Computes the max-min values of the grid metrics.
Definition Metric.c:1670
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:368
PetscErrorCode ComputeCellCentersAndSpacing(UserCtx *user)
Computes the physical location of cell centers and the spacing between them.
Definition Metric.c:829
PetscErrorCode MetricGetCellVertices(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, Cmpnts V[8])
Extract the eight vertex coordinates of the hexahedral cell (i,j,k).
Definition Metric.c:27
PetscErrorCode ComputeJFaceMetrics(UserCtx *user)
Computes metrics centered on constant-j faces (j-faces).
Definition Metric.c:1151
PetscErrorCode ComputeMetricsDivergence(UserCtx *user)
Performs a diagnostic check on the divergence of the face area metric vectors.
Definition Metric.c:1573
PetscErrorCode ComputeFaceMetrics(UserCtx *user)
Computes the primary face metric components (Csi, Eta, Zet), including boundary extrapolation,...
Definition Metric.c:464
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:313
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:257
PetscErrorCode ComputeCellCenteredJacobianInverse(UserCtx *user)
Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells and stores it in u...
Definition Metric.c:676
PetscErrorCode ComputeKFaceMetrics(UserCtx *user)
Computes metrics centered on constant-k faces (k-faces).
Definition Metric.c:1359
PetscErrorCode MetricLogicalToPhysical(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, PetscReal xi, PetscReal eta, PetscReal zta, Cmpnts *Xp)
Public wrapper: map (cell index, ξ,η,ζ) to (x,y,z).
Definition Metric.c:77
PetscErrorCode ComputeIFaceMetrics(UserCtx *user)
Computes metrics centered on constant-i faces (i-faces).
Definition Metric.c:929
PetscErrorCode MetricVelocityContravariant(const PetscReal J[3][3], PetscReal detJ, const PetscReal u[3], PetscReal uc[3])
Convert physical velocity (u,v,w) to contravariant components (u^xi, u^eta, u^zta).
Definition Metric.c:167
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:538
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661