|
PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
|
Go to the source code of this file.
Macros | |
| #define | __FUNCT__ "MetricGetCellVertices" |
| #define | __FUNCT__ "TrilinearBlend" |
| #define | __FUNCT__ "MetricLogicalToPhysical" |
| #define | __FUNCT__ "MetricJacobian" |
| #define | __FUNCT__ "MetricVelocityContravariant" |
| #define | __FUNCT__ "InvertCovariantMetricTensor" |
| #define | __FUNCT__ "CalculateFaceNormalAndArea" |
| #define | __FUNCT__ "ComputeCellCharacteristicLengthScale" |
| #define | __FUNCT__ "CheckAndFixGridOrientation" |
| #define | __FUNCT__ "ApplyPeriodicCorrectionsToCellCentersAndSpacing" |
| #define | __FUNCT__ "ApplyPeriodicCorrectionsToIFaceCenter" |
| #define | __FUNCT__ "ApplyPeriodicCorrectionsToJFaceCenter" |
| #define | __FUNCT__ "ApplyPeriodicCorrectionsToKFaceCenter" |
| #define | __FUNCT__ "ComputeFaceMetrics" |
| #define | __FUNCT__ "ComputeCellCenteredJacobianInverse" |
| #define | __FUNCT__ "ComputeCellCentersAndSpacing" |
| #define | __FUNCT__ "ComputeIFaceMetrics" |
| #define | __FUNCT__ "ComputeJFaceMetrics" |
| #define | __FUNCT__ "ComputeJFaceMetrics" |
| #define | __FUNCT__ "ComputeMetricsDivergence" |
| #define | __FUNCT__ "ComputeMetricNorms" |
| #define | __FUNCT__ "ComputeAllGridMetrics" |
Functions | |
| 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). | |
| static void | TrilinearBlend (const Cmpnts V[8], PetscReal xi, PetscReal eta, PetscReal zta, Cmpnts *Xp) |
| Map logical coordinates to physical space using trilinear basis. | |
| 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). | |
| 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). | |
| 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). | |
| PetscErrorCode | InvertCovariantMetricTensor (double covariantTensor[3][3], double contravariantTensor[3][3]) |
| Inverts the 3x3 covariant metric tensor to obtain the contravariant metric tensor. | |
| 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. | |
| 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. | |
| PetscErrorCode | CheckAndFixGridOrientation (UserCtx *user) |
Ensure a right-handed metric basis (Csi, Eta, Zet) and a positive Jacobian (Aj) over the whole domain. | |
| PetscErrorCode | ApplyPeriodicCorrectionsToCellCentersAndSpacing (UserCtx *user) |
| Applies periodic boundary corrections to cell centers (Cent) and grid spacing (GridSpace). | |
| PetscErrorCode | ApplyPeriodicCorrectionsToIFaceCenter (UserCtx *user) |
| Applies periodic boundary corrections to i-face centers (Centx). | |
| PetscErrorCode | ApplyPeriodicCorrectionsToJFaceCenter (UserCtx *user) |
| Applies periodic boundary corrections to j-face centers (Centy). | |
| PetscErrorCode | ApplyPeriodicCorrectionsToKFaceCenter (UserCtx *user) |
| Applies periodic boundary corrections to k-face centers (Centz). | |
| PetscErrorCode | ComputeFaceMetrics (UserCtx *user) |
| Computes the primary face metric components (Csi, Eta, Zet), including boundary extrapolation, and stores them in the corresponding global Vec members of the UserCtx structure (user->Csi, user->Eta, user->Zet). | |
| PetscErrorCode | ComputeCellCenteredJacobianInverse (UserCtx *user) |
Calculates the cell-centered inverse Jacobian determinant (1/J), including boundary extrapolation, stores it in user->Aj, assembles user->Aj, and updates user->lAj. | |
| PetscErrorCode | ComputeCellCentersAndSpacing (UserCtx *user) |
| Computes the physical location of cell centers and the spacing between them. | |
| PetscErrorCode | ComputeIFaceMetrics (UserCtx *user) |
| Computes metrics centered on constant-i faces (i-faces). | |
| PetscErrorCode | ComputeJFaceMetrics (UserCtx *user) |
| Computes metrics centered on constant-j faces (j-faces). | |
| PetscErrorCode | ComputeKFaceMetrics (UserCtx *user) |
| Computes metrics centered on constant-k faces (k-faces). | |
| static PetscInt | Gidx (PetscInt i, PetscInt j, PetscInt k, UserCtx *user) |
| PetscErrorCode | ComputeMetricsDivergence (UserCtx *user) |
| Computes the divergence of the grid metrics and identifies the maximum value. | |
| PetscErrorCode | ComputeMetricNorms (UserCtx *user) |
| Computes the max-min values of the grid metrics. | |
| PetscErrorCode | CalculateAllGridMetrics (SimCtx *simCtx) |
| Orchestrates the calculation of all grid metrics. | |
| #define __FUNCT__ "ApplyPeriodicCorrectionsToCellCentersAndSpacing" |
| 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).
Vertices are returned in the standard trilinear ordering: bit 0 → x-corner, bit 1 → y-corner, bit 2 → z-corner. (000 = origin of the cell, 111 = far corner.)
Definition at line 27 of file Metric.c.
|
inlinestatic |
Map logical coordinates to physical space using trilinear basis.
| [in] | V | Array of the eight vertex coordinates (MetricGetCellVertices). |
| [in] | xi,eta,zta | Logical coordinates in [0,1]. |
| [out] | Xp | Physical coordinate. |
Definition at line 54 of file Metric.c.
| 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 at line 77 of file Metric.c.
| 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).
J = [ x_ξ x_η x_ζ ]
[ y_ξ y_η y_ζ ]
[ z_ξ z_η z_ζ ]
This is handy for converting physical velocities (u,v,w) into contravariant components and for volume weighting.
Definition at line 110 of file Metric.c.
| 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 at line 167 of file Metric.c.
| PetscErrorCode InvertCovariantMetricTensor | ( | double | covariantTensor[3][3], |
| double | contravariantTensor[3][3] | ||
| ) |
Inverts the 3x3 covariant metric tensor to obtain the contravariant metric tensor.
In curvilinear coordinates, the input matrix g contains the dot products of the covariant basis vectors (e.g., g_ij = e_i . e_j). Its inverse, G, is the contravariant metric tensor, which is essential for transforming vectors and tensors between coordinate systems.
| covariantTensor | Input: A 3x3 matrix representing the covariant metric tensor. |
| contravariantTensor | Output: A 3x3 matrix where the inverted result is stored. |
Definition at line 212 of file Metric.c.
| 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.
Given the metric vectors (csi, eta, zet), this function calculates the geometric properties of the cell faces aligned with the i, j, and k directions.
| csi,eta,zet | The metric vectors at the cell center. |
| ni | Output: A 3-element array for the unit normal vector of the i-face. |
| nj | Output: A 3-element array for the unit normal vector of the j-face. |
| nk | Output: A 3-element array for the unit normal vector of the k-face. |
| Ai | Output: Pointer to store the area of the i-face. |
| Aj | Output: Pointer to store the area of the j-face. |
| Ak | Output: Pointer to store the area of the k-face. |
Definition at line 257 of file Metric.c.
| 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.
For a non-uniform, non-orthogonal cell, there is no single "dx". This function computes an effective length scale in each Cartesian direction based on the cell volume and the areas of its faces.
| ajc | The Jacobian of the grid transformation (1 / cell volume). |
| csi,eta,zet | The metric vectors at the cell center. |
| dx | Output: Pointer to store the characteristic length in the x-direction. |
| dy | Output: Pointer to store the characteristic length in the y-direction. |
| dz | Output: Pointer to store the characteristic length in the z-direction. |
Definition at line 313 of file Metric.c.
| PetscErrorCode CheckAndFixGridOrientation | ( | UserCtx * | user | ) |
Ensure a right-handed metric basis (Csi, Eta, Zet) and a positive Jacobian (Aj) over the whole domain.
The metric-generation kernels are completely algebraic, so they will happily deliver a left-handed basis if the mesh file enumerates nodes in the opposite ζ-direction.
This routine makes the orientation explicit and—if needed—repairs it once per run:
| Step | Action |
|---|---|
| 1 | Compute global Aj_min, Aj_max. |
| 2 | Mixed signs (Aj_min < 0 && Aj_max > 0) → abort: the mesh is topologically inconsistent. |
| 3 | All negative (Aj_max < 0) → flip Csi, Eta, Zet, Aj & update local ghosts. |
| 4 | Store user->orientation = ±1 so BC / IC routines can apply sign-aware logic if they care about inlet direction. |
| [in,out] | user | Fully initialised UserCtx that already contains Csi, Eta, Zet, Aj, their local ghosts, and valid distributed DMs. |
0 on success or a PETSc error code on failure.ComputeCellCenteredJacobianInverse() and before any routine that differentiates or applies BCs.Definition at line 368 of file Metric.c.
| PetscErrorCode ApplyPeriodicCorrectionsToCellCentersAndSpacing | ( | UserCtx * | user | ) |
Applies periodic boundary corrections to cell centers (Cent) and grid spacing (GridSpace).
This function handles the special logic needed when periodic boundaries are present. For coarse grids (cgrid), it directly copies from ghost regions. For fine grids, it calculates the boundary values using the GridSpace information.
| user | The UserCtx containing grid and field data. |
Definition at line 452 of file Metric.c.
| PetscErrorCode ApplyPeriodicCorrectionsToIFaceCenter | ( | UserCtx * | user | ) |
Applies periodic boundary corrections to i-face centers (Centx).
Only X-direction periodicity affects Centx. This function must be called after Centx has been initially computed but before it's used for metric calculations.
| user | The UserCtx containing grid and field data. |
Definition at line 656 of file Metric.c.
| PetscErrorCode ApplyPeriodicCorrectionsToJFaceCenter | ( | UserCtx * | user | ) |
Applies periodic boundary corrections to j-face centers (Centy).
Only Y-direction periodicity affects Centy. This function must be called after Centy has been initially computed but before it's used for metric calculations.
| user | The UserCtx containing grid and field data. |
Definition at line 742 of file Metric.c.
| PetscErrorCode ApplyPeriodicCorrectionsToKFaceCenter | ( | UserCtx * | user | ) |
Applies periodic boundary corrections to k-face centers (Centz).
Applies periodic boundary corrections to i-face centers (Centx).
Only Z-direction periodicity affects Centz. This function must be called after Centz has been initially computed but before it's used for metric calculations.
| user | The UserCtx containing grid and field data. |
Definition at line 828 of file Metric.c.
| PetscErrorCode ComputeFaceMetrics | ( | UserCtx * | user | ) |
Computes the primary face metric components (Csi, Eta, Zet), including boundary extrapolation, and stores them in the corresponding global Vec members of the UserCtx structure (user->Csi, user->Eta, user->Zet).
This is a self-contained routine that performs the following steps:
user->Csi, user->Eta, user->Zet Vecs.user->lCsi, user->lEta, user->lZet Vecs.This function calculates the face area vectors, which are fundamental to the finite volume method on a curvilinear grid. The process is performed in two main stages to ensure robustness and correctness across the entire domain, including physical boundaries.
Stage 1: Interior Face Calculation The core of the function uses centered finite-difference stencils on the nodal coordinates (coor) to compute the metric terms. For example, csi at a node (k,j,i) is calculated using a 2x2 stencil of nodes in the J-K plane.
Crucially, these stencils are only valid for interior nodes/faces where all required neighboring nodes exist. The loops are intentionally constructed (e.g., i_loop_start = 1) to skip the calculation directly on the domain's physical boundaries (e.g., at i=0, j=0, etc.), as the stencil would require out-of-bounds data.
Stage 2: Boundary Face Extrapolation After the interior metrics are computed, the values for the faces lying on the physical boundaries of the domain are populated. This is done via zero-order extrapolation, which simply means copying the values from the nearest valid interior face layer.
For example, on a rank owning the global i=0 boundary, zet_arr[k][j][0] is explicitly set to be equal to zet_arr[k][j][1]. Similarly, on a rank owning the i=mx-1 boundary, zet_arr[k][j][mx-1] is set to zet_arr[k][j][mx-2].
This two-stage approach ensures that every node in the physical domain, including the boundaries, has a valid metric value associated with it.
The function concludes by assembling the global PETSc vectors and updating the local ghosted versions, making the computed metrics ready for use by the solver.
| [in,out] | user | Pointer to the UserCtx structure. |
VecZeroEntries on user->Csi, Eta, Zet before this if they might contain old data.Definition at line 956 of file Metric.c.
| PetscErrorCode ComputeCellCenteredJacobianInverse | ( | UserCtx * | user | ) |
Calculates the cell-centered inverse Jacobian determinant (1/J), including boundary extrapolation, stores it in user->Aj, assembles user->Aj, and updates user->lAj.
Calculates the cell-centered inverse Jacobian determinant (1/J) for INTERIOR cells and stores it in user->Aj.
| [in,out] | user | Pointer to the UserCtx structure. |
Definition at line 1171 of file Metric.c.
| PetscErrorCode ComputeCellCentersAndSpacing | ( | UserCtx * | user | ) |
Computes the physical location of cell centers and the spacing between them.
This function calculates two key geometric properties from the nodal coordinates:
Cent: A vector field storing the (x,y,z) coordinates of the center of each grid cell.GridSpace: A vector field storing the physical distance between adjacent cell centers in the i, j, and k computational directions.It is a direct adaptation of the corresponding logic from the legacy FormMetrics.
| user | The UserCtx for a specific grid level. The function populates user->Cent and user->GridSpace. |
Definition at line 1326 of file Metric.c.
| PetscErrorCode ComputeIFaceMetrics | ( | UserCtx * | user | ) |
Computes metrics centered on constant-i faces (i-faces).
This function calculates the metric terms (ICsi, IEta, IZet) and the inverse Jacobian (IAj) located at the geometric center of each constant-i face. This is a critical step for staggered-grid finite difference schemes.
The process is a direct and faithful refactoring of the corresponding logic from the legacy FormMetrics function:
user->Centx vector.Centx field to compute the derivatives (e.g., d(x)/d(csi)).Vec objects.| user | The UserCtx for a specific grid level. This function populates the user->ICsi, user->IEta, user->IZet, and user->IAj vectors. |
Definition at line 1450 of file Metric.c.
| PetscErrorCode ComputeJFaceMetrics | ( | UserCtx * | user | ) |
Computes metrics centered on constant-j faces (j-faces).
This function calculates the metric terms (JCsi, JEta, JZet) and the inverse Jacobian (JAj) located at the geometric center of each constant-j face. This is a critical step for staggered-grid finite difference schemes.
The process is a direct and faithful refactoring of the corresponding logic from the legacy FormMetrics function:
user->Centy vector.Centy field to compute the derivatives (e.g., d(x)/d(csi)).Vec objects.| user | The UserCtx for a specific grid level. This function populates the user->JCsi, user->JEta, user->JZet, and user->JAj vectors. |
Definition at line 1692 of file Metric.c.
| PetscErrorCode ComputeKFaceMetrics | ( | UserCtx * | user | ) |
Computes metrics centered on constant-k faces (k-faces).
This function calculates the metric terms (KCsi, KEta, KZet) and the inverse Jacobian (KAj) located at the geometric center of each constant-k face. This is a critical step for staggered-grid finite difference schemes.
The process is a direct and faithful refactoring of the corresponding logic from the legacy FormMetrics function:
user->Centz vector.Centz field to compute the derivatives (e.g., d(x)/d(csi)).Vec objects.| user | The UserCtx for a specific grid level. This function populates the user->KCsi, user->KEta, user->KZet, and user->KAj vectors. |
Definition at line 1920 of file Metric.c.
|
static |
Definition at line 2122 of file Metric.c.
| PetscErrorCode ComputeMetricsDivergence | ( | UserCtx * | user | ) |
Computes the divergence of the grid metrics and identifies the maximum value.
Performs a diagnostic check on the divergence of the face area metric vectors.
This function serves as a diagnostic tool to assess the quality of the grid metrics. It calculates the divergence of the face metrics (Csi, Eta, Zet) and scales it by the inverse of the cell Jacobian. The maximum divergence value is then located, and its grid coordinates are printed to the console, helping to pinpoint areas of potential grid quality issues.
| user | The UserCtx, containing all necessary grid data. |
Definition at line 2152 of file Metric.c.
| PetscErrorCode ComputeMetricNorms | ( | UserCtx * | user | ) |
Computes the max-min values of the grid metrics.
This function serves as a diagnostic tool to assess the quality of the grid metrics. It calculates the bounds of the face metrics (Csi, Eta, Zet).
| user | The UserCtx, containing all necessary grid data. |
Definition at line 2249 of file Metric.c.
| PetscErrorCode CalculateAllGridMetrics | ( | SimCtx * | simCtx | ) |
Orchestrates the calculation of all grid metrics.
This function iterates through every UserCtx in the multigrid and multi-block hierarchy. For each context, it calls a series of modern, modular helper functions to compute the face metrics (Csi, Eta, Zet), the cell-centered inverse Jacobian (Aj), and to validate the grid's orientation.
Definition at line 2388 of file Metric.c.