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 primary face metric components (Csi, Eta, Zet), including
32 * boundary extrapolation, and stores them in the corresponding global Vec
33 * members of the UserCtx structure (user->Csi, user->Eta, user->Zet).
34 *
35 * This is a self-contained routine that performs the following steps:
36 * 1. Obtains local ghosted nodal coordinates using DMGetCoordinatesLocal.
37 * 2. Calculates metrics for INTERIOR faces where finite difference stencils are valid.
38 * 3. EXTRAPOLATES metrics for faces on the physical domain boundaries by copying
39 * from the nearest computed interior face.
40 * 4. Assembles the global `user->Csi`, `user->Eta`, `user->Zet` Vecs.
41 * 5. Updates the local ghosted `user->lCsi`, `user->lEta`, `user->lZet` Vecs.
42 *
43 * @param[in,out] user Pointer to the UserCtx structure.
44 *
45 * @return PetscErrorCode 0 on success.
46 *
47 * @note
48 * - This function is a complete "compute and make ready" unit for Csi, Eta, and Zet.
49 * - It's recommended to call `VecZeroEntries` on user->Csi, Eta, Zet before this
50 * if they might contain old data.
51 */
52PetscErrorCode ComputeFaceMetrics(UserCtx *user);
53
54/**
55 * @brief Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells
56 * and stores it in `user->Aj`. This version includes boundary extrapolation.
57 *
58 * Nodal coordinates are obtained internally.
59 * Refer to previous Doxygen comments for details on physical locations and
60 * storage convention (`aj_arr[k_n][j_n][i_n]` for cell `C(i_n-1,j_n-1,k_n-1)`).
61 *
62 * @param[in,out] user Pointer to the UserCtx structure.
63 *
64 * @return PetscErrorCode 0 on success.
65 */
67
68/**
69 * @brief Ensure a **right-handed** metric basis (`Csi`, `Eta`, `Zet`) and a
70 * **positive Jacobian** (`Aj`) over the whole domain.
71 *
72 * The metric-generation kernels are completely algebraic, so they will happily
73 * deliver a *left-handed* basis if the mesh file enumerates nodes in the
74 * opposite ζ-direction.
75 * This routine makes the orientation explicit and—if needed—repairs it
76 * **once per run**:
77 *
78 * | Step | Action |
79 * |------|--------|
80 * | 1 | Compute global `Aj_min`, `Aj_max`. |
81 * | 2 | **Mixed signs** (`Aj_min < 0 && Aj_max > 0`) &rarr; abort: the mesh is topologically inconsistent. |
82 * | 3 | **All negative** (`Aj_max < 0`) &rarr; flip <br>`Csi`, `Eta`, `Zet`, `Aj` & update local ghosts. |
83 * | 4 | Store `user->orientation = ±1` so BC / IC routines can apply sign-aware logic if they care about inlet direction. |
84 *
85 * @param[in,out] user Fully initialised #UserCtx that already contains
86 * `Csi`, `Eta`, `Zet`, `Aj`, their **local** ghosts, and
87 * valid distributed DMs.
88 *
89 * @return `0` on success or a PETSc error code on failure.
90 *
91 * @note Call **immediately after** `ComputeCellCenteredJacobianInverse()` and
92 * before any routine that differentiates or applies BCs.
93 *
94 * @author <your name>
95 */
96PetscErrorCode CheckAndFixGridOrientation(UserCtx *user);
97
98/**
99 * @brief Computes the physical location of cell centers and the spacing between them.
100 *
101 * This function calculates two key geometric properties from the nodal coordinates:
102 * 1. `Cent`: A vector field storing the (x,y,z) coordinates of the center of each grid cell.
103 * 2. `GridSpace`: A vector field storing the physical distance between adjacent
104 * cell centers in the i, j, and k computational directions.
105 *
106 * It is a direct adaptation of the corresponding logic from the legacy `FormMetrics`.
107 *
108 * @param user The UserCtx for a specific grid level. The function populates `user->Cent` and `user->GridSpace`.
109 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
110 */
111PetscErrorCode ComputeCellCentersAndSpacing(UserCtx *user);
112
113/**
114 * @brief Computes metrics centered on constant-i faces (i-faces).
115 *
116 * This function calculates the metric terms (ICsi, IEta, IZet) and the inverse
117 * Jacobian (IAj) located at the center of each i-face. The stencils use
118 * i-face-centered coordinates (`Centx`) which must be computed first.
119 * The logic is a direct adaptation of the legacy FormMetrics function.
120 *
121 * @param user The UserCtx for a specific grid level. Populates user->ICsi, etc.
122 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
123 */
124PetscErrorCode ComputeIFaceMetrics(UserCtx *user);
125
126/**
127 * @brief Computes metrics centered on constant-j faces (j-faces).
128 *
129 * This function calculates the metric terms (`JCsi`, `JEta`, `JZet`) and the
130 * inverse Jacobian (`JAj`) located at the geometric center of each constant-j
131 * face. This is a critical step for staggered-grid finite difference schemes.
132 *
133 * The process is a direct and faithful refactoring of the corresponding logic
134 * from the legacy `FormMetrics` function:
135 * 1. It first calculates the physical (x,y,z) coordinates of the center of
136 * each i-face and stores them in the `user->Centy` vector.
137 * 2. It then uses a boundary-aware, second-order finite difference stencil on
138 * the `Centy` field to compute the derivatives (e.g., d(x)/d(csi)).
139 * - Central differences are used in the grid interior.
140 * - One-sided differences are used at the physical domain boundaries.
141 * 3. Finally, these derivatives are used to compute the final metric terms and
142 * the inverse Jacobian, which are stored in their respective `Vec` objects.
143 *
144 * @param user The UserCtx for a specific grid level. This function populates
145 * the `user->JCsi`, `user->JEta`, `user->JZet`, and `user->JAj` vectors.
146 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
147 */
148PetscErrorCode ComputeJFaceMetrics(UserCtx *user);
149
150/**
151 * @brief Computes metrics centered on constant-k faces (k-faces).
152 *
153 * This function calculates the metric terms (`KCsi`, `KEta`, `KZet`) and the
154 * inverse Jacobian (`KAj`) located at the geometric center of each constant-j
155 * face. This is a critical step for staggered-grid finite difference schemes.
156 *
157 * The process is a direct and faithful refactoring of the corresponding logic
158 * from the legacy `FormMetrics` function:
159 * 1. It first calculates the physical (x,y,z) coordinates of the center of
160 * each i-face and stores them in the `user->Centz` vector.
161 * 2. It then uses a boundary-aware, second-order finite difference stencil on
162 * the `Centz` field to compute the derivatives (e.g., d(x)/d(csi)).
163 * - Central differences are used in the grid interior.
164 * - One-sided differences are used at the physical domain boundaries.
165 * 3. Finally, these derivatives are used to compute the final metric terms and
166 * the inverse Jacobian, which are stored in their respective `Vec` objects.
167 *
168 * @param user The UserCtx for a specific grid level. This function populates
169 * the `user->KCsi`, `user->KEta`, `user->KZet`, and `user->KAj` vectors.
170 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
171 */
172PetscErrorCode ComputeKFaceMetrics(UserCtx *user);
173
174/**
175 * @brief Performs a diagnostic check on the divergence of the face area metric vectors.
176 *
177 * For a closed cell, the sum of the face area vectors should be zero (Gauss's
178 * divergence theorem). This function computes a measure of this divergence and
179 * reports the maximum value over the domain. A small value indicates a
180 * well-formed grid. This is a direct adaptation of the legacy function.
181 *
182 * @param user The UserCtx for a specific grid level (typically the finest).
183 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
184 */
185PetscErrorCode ComputeMetricsDivergence(UserCtx *user);
186
187/**
188 * @brief Computes the max-min values of the grid metrics.
189 *
190 * This function serves as a diagnostic tool to assess the quality of the grid
191 * metrics. It calculates the bounds of the face metrics (Csi, Eta, Zet).
192 *
193 * @param user The UserCtx, containing all necessary grid data.
194 * @return PetscErrorCode
195 */
196PetscErrorCode ComputeMetricNorms(UserCtx *user);
197
198/**
199 * @brief Orchestrates the calculation of all grid metrics.
200 *
201 * This function iterates through every UserCtx in the multigrid and multi-block
202 * hierarchy. For each context, it calls a series of modern, modular helper
203 * functions to compute the face metrics (Csi, Eta, Zet), the cell-centered
204 * inverse Jacobian (Aj), and to validate the grid's orientation.
205 *
206 * @param simCtx The master SimCtx, containing the configured UserCtx hierarchy.
207 * @return PetscErrorCode
208 */
209PetscErrorCode CalculateAllGridMetrics(SimCtx *simCtx);
210#endif /* METRIC_H */
PetscErrorCode CalculateAllGridMetrics(SimCtx *simCtx)
Orchestrates the calculation of all grid metrics.
Definition Metric.c:1639
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:93
PetscErrorCode ComputeMetricNorms(UserCtx *user)
Computes the max-min values of the grid metrics.
Definition Metric.c:1506
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:195
PetscErrorCode ComputeCellCentersAndSpacing(UserCtx *user)
Computes the physical location of cell centers and the spacing between them.
Definition Metric.c:644
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:24
PetscErrorCode ComputeJFaceMetrics(UserCtx *user)
Computes metrics centered on constant-j faces (j-faces).
Definition Metric.c:942
PetscErrorCode ComputeMetricsDivergence(UserCtx *user)
Performs a diagnostic check on the divergence of the face area metric vectors.
Definition Metric.c:1416
PetscErrorCode ComputeFaceMetrics(UserCtx *user)
Computes the primary face metric components (Csi, Eta, Zet), including boundary extrapolation,...
Definition Metric.c:286
PetscErrorCode ComputeCellCenteredJacobianInverse(UserCtx *user)
Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells and stores it in u...
Definition Metric.c:491
PetscErrorCode ComputeKFaceMetrics(UserCtx *user)
Computes metrics centered on constant-k faces (k-faces).
Definition Metric.c:1145
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:68
PetscErrorCode ComputeIFaceMetrics(UserCtx *user)
Computes metrics centered on constant-i faces (i-faces).
Definition Metric.c:739
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:141
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:99
The master context for the entire simulation.
Definition variables.h:513
User-defined context containing data specific to a single computational grid level.
Definition variables.h:630