|
PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
|
#include <petscpf.h>#include <petscdmswarm.h>#include <stdlib.h>#include <time.h>#include <math.h>#include <petsctime.h>#include <petscdmcomposite.h>#include <petscerror.h>#include <petscsys.h>#include "variables.h"#include "ParticleSwarm.h"#include "walkingsearch.h"#include "grid.h"#include "logging.h"#include "io.h"#include "setup.h"Go to the source code of this file.
Macros | |
| #define | NUM_WEIGHTS 8 |
| #define | InterpolateFieldFromCornerToCenter(field, centfield, user) |
| Generic macro to call the appropriate interpolation function based on the field type. | |
| #define | InterpolateFieldFromCenterToCorner(blockSize, centfield_ptr, corner_ptr, user_ctx) |
| Macro to dispatch to the correct scalar or vector center-to-corner function based on a runtime block size variable. | |
| #define | InterpolateCornerToFaceCenter(corner_arr, faceX_arr, faceY_arr, faceZ_arr, user_ctx) |
| A type-generic macro that interpolates a field from corner nodes to all face centers. | |
| #define | PieceWiseLinearInterpolation(fieldName, fieldPtr, i, j, k, outPtr) |
Macro that calls either the scalar or vector piecewise interpolation function based on the type of the fieldPtr parameter (3D array). | |
| #define | TrilinearInterpolation(fieldName, fieldPtr, i, j, k, a1, a2, a3, outPtr) |
Macro that calls either the scalar or vector trilinear interpolation function based on the type of the fieldPtr parameter (3D array). | |
Functions | |
| PetscErrorCode | PieceWiseLinearInterpolation_Scalar (const char *fieldName, PetscReal ***fieldScal, PetscInt iCell, PetscInt jCell, PetscInt kCell, PetscReal *val) |
| Returns the scalar value from the input cell index without blending. | |
| PetscErrorCode | PieceWiseLinearInterpolation_Vector (const char *fieldName, Cmpnts ***fieldVec, PetscInt iCell, PetscInt jCell, PetscInt kCell, Cmpnts *vec) |
| Returns the vector value from the input cell index without blending. | |
| PetscErrorCode | TrilinearInterpolation_Scalar (const char *fieldName, PetscReal ***fieldScal, PetscInt i, PetscInt j, PetscInt k, PetscReal a1, PetscReal a2, PetscReal a3, PetscReal *val) |
| Computes the trilinear interpolated scalar at a given point. | |
| PetscErrorCode | TrilinearInterpolation_Vector (const char *fieldName, Cmpnts ***fieldVec, PetscInt i, PetscInt j, PetscInt k, PetscReal a1, PetscReal a2, PetscReal a3, Cmpnts *vec) |
| Computes the trilinear interpolated vector (e.g., velocity) at a given point. | |
| PetscErrorCode | InterpolateEulerFieldToSwarm (UserCtx *user, Vec fieldLocal_cellCentered, const char *fieldName, const char *swarmOutFieldName) |
| Interpolates a cell-centered field (scalar or vector) onto DMSwarm particles, using a robust, PETSc-idiomatic two-stage process. | |
| PetscErrorCode | InterpolateAllFieldsToSwarm (UserCtx *user) |
| Interpolates all relevant fields from the DMDA to the DMSwarm. | |
| PetscErrorCode | InterpolateParticleVelocities (UserCtx *user) |
| Interpolates particle velocities using trilinear interpolation. | |
| PetscErrorCode | InterpolateFieldFromCornerToCenter_Scalar (PetscReal ***field_arr, PetscReal ***centfield_arr, UserCtx *user) |
| Safely interpolate a scalar field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx. | |
| PetscErrorCode | InterpolateFieldFromCornerToCenter_Vector (Cmpnts ***field_arr, Cmpnts ***centfield_arr, UserCtx *user) |
| Safely interpolate a vector field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx. | |
| PetscErrorCode | TestCornerToCenterInterpolation (UserCtx *user) |
| Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector. | |
| PetscErrorCode | InterpolateFieldFromCenterToCorner_Vector (Cmpnts ***centfield_arr, Cmpnts ***corner_arr, UserCtx *user) |
| Interpolates a vector field from cell centers to corner nodes. | |
| PetscErrorCode | InterpolateFieldFromCenterToCorner_Scalar (PetscReal ***centfield_arr, PetscReal ***corner_arr, UserCtx *user) |
| Interpolates a scalar field from cell centers to corner nodes. | |
| PetscErrorCode | GetScatterTargetInfo (UserCtx *user, const char *particleFieldName, DM *targetDM, PetscInt *expected_dof) |
| Determines the target Eulerian DM and expected DOF for scattering a given particle field. | |
| PetscErrorCode | GetPersistentLocalVector (UserCtx *user, const char *fieldName, Vec *localVec) |
| Retrieves the persistent local vector (e.g., lPsi, lUcat) for a given field name. | |
| PetscErrorCode | AccumulateParticleField (DM swarm, const char *particleFieldName, DM gridSumDM, Vec gridSumVec) |
| Accumulates a particle field (scalar or vector) into a target grid sum vector. | |
| PetscErrorCode | NormalizeGridVectorByCount (DM countDM, Vec countVec, DM dataDM, Vec sumVec, Vec avgVec) |
| Normalizes a grid vector of sums by a grid vector of counts to produce an average. | |
| PetscErrorCode | ScatterParticleFieldToEulerField (UserCtx *user, const char *particleFieldName, Vec eulerFieldAverageVec) |
| Scatters a particle field (scalar or vector) to the corresponding Eulerian field average. | |
| PetscErrorCode | ScatterAllParticleFieldsToEulerFields (UserCtx *user) |
| Scatters a predefined set of particle fields to their corresponding Eulerian fields. | |
| PetscErrorCode | InterpolateCornerToFaceCenter_Scalar (PetscReal ***corner_arr, PetscReal ***faceX_arr, PetscReal ***faceY_arr, PetscReal ***faceZ_arr, UserCtx *user) |
| Interpolates a scalar field from corner nodes to all face centers. | |
| PetscErrorCode | InterpolateCornerToFaceCenter_Vector (Cmpnts ***corner_arr, Cmpnts ***faceX_arr, Cmpnts ***faceY_arr, Cmpnts ***faceZ_arr, UserCtx *user) |
| Interpolates a vector field from corner nodes to all face centers. | |
| #define NUM_WEIGHTS 8 |
Definition at line 24 of file interpolation.h.
| #define InterpolateFieldFromCornerToCenter | ( | field, | |
| centfield, | |||
| user | |||
| ) |
Generic macro to call the appropriate interpolation function based on the field type.
This macro will select either the scalar or vector interpolation function based on the type of the 'field' pointer. It also performs a compile-time check that 'centfield' is of the same type as 'field'. If the types are not the same, a compile-time error is produced.
Usage: InterpolateFieldFromCornerToCenter(field, centfield, user);
Definition at line 36 of file interpolation.h.
| #define InterpolateFieldFromCenterToCorner | ( | blockSize, | |
| centfield_ptr, | |||
| corner_ptr, | |||
| user_ctx | |||
| ) |
Macro to dispatch to the correct scalar or vector center-to-corner function based on a runtime block size variable.
This macro uses a ternary operator to inspect the runtime value of 'blockSize' and select the appropriate implementation. It expects the input and output pointers to be void* and handles casting them to the correct, strongly-typed pointers.
| blockSize | The runtime block size (1 for scalar, 3 for vector). |
| centfield_ptr | The input void* pointer to the 3D cell-centered data array. |
| corner_ptr | The output void* pointer to the 3D corner data array. |
| user_ctx | The UserCtx structure. |
Definition at line 56 of file interpolation.h.
| #define InterpolateCornerToFaceCenter | ( | corner_arr, | |
| faceX_arr, | |||
| faceY_arr, | |||
| faceZ_arr, | |||
| user_ctx | |||
| ) |
A type-generic macro that interpolates a field from corner nodes to all face centers.
This macro uses C11's _Generic feature to dispatch to the appropriate underlying function based on the type of the input corner_arr.
corner_arr is of type PetscReal***, it calls InterpolateCornerToFaceCenter_Scalar.corner_arr is of type Cmpnts***, it calls InterpolateCornerToFaceCenter_Vector.| corner_arr | Ghosted node-centered array (global indexing). The type of this argument determines which function is called. |
| faceX_arr | Local array for X-faces. |
| faceY_arr | Local array for Y-faces. |
| faceZ_arr | Local array for Z-faces. |
| user_ctx | User context containing DMDA 'fda'. |
Definition at line 79 of file interpolation.h.
| #define PieceWiseLinearInterpolation | ( | fieldName, | |
| fieldPtr, | |||
| i, | |||
| j, | |||
| k, | |||
| outPtr | |||
| ) |
Macro that calls either the scalar or vector piecewise interpolation function based on the type of the fieldPtr parameter (3D array).
Usage example:
// For scalar: PetscReal ***fieldScal; PetscReal outVal; PieceWiseLinearInterpolation(fieldName, fieldScal, i, j, k, &outVal);
// For vector: Cmpnts ***fieldVec; Cmpnts vec; PieceWiseLinearInterpolation(fieldName, fieldVec, i, j, k, &vec);
Definition at line 101 of file interpolation.h.
| #define TrilinearInterpolation | ( | fieldName, | |
| fieldPtr, | |||
| i, | |||
| j, | |||
| k, | |||
| a1, | |||
| a2, | |||
| a3, | |||
| outPtr | |||
| ) |
Macro that calls either the scalar or vector trilinear interpolation function based on the type of the fieldPtr parameter (3D array).
Usage example:
PetscReal result; Cmpnts vec;
// For scalars: TrilinearInterpolation(fieldName, fieldScal, i, j, k, a1, a2, a3, &result);
// For vectors: TrilinearInterpolation(fieldName, fieldVec, i, j, k, a1, a2, a3, &vec);
Definition at line 166 of file interpolation.h.
| PetscErrorCode PieceWiseLinearInterpolation_Scalar | ( | const char * | fieldName, |
| PetscReal *** | fieldScal, | ||
| PetscInt | iCell, | ||
| PetscInt | jCell, | ||
| PetscInt | kCell, | ||
| PetscReal * | val | ||
| ) |
Returns the scalar value from the input cell index without blending.
This is a first-order piecewise-constant lookup helper used in workflows that intentionally disable trilinear blending.
| [in] | fieldName | Field label used for diagnostics. |
| [in] | fieldScal | Scalar field array indexed as [k][j][i]. |
| [in] | iCell | Cell i-index. |
| [in] | jCell | Cell j-index. |
| [in] | kCell | Cell k-index. |
| [out] | val | Output scalar value. |
Returns the scalar value from the input cell index without blending.
Local to this translation unit.
Definition at line 462 of file interpolation.c.
| PetscErrorCode PieceWiseLinearInterpolation_Vector | ( | const char * | fieldName, |
| Cmpnts *** | fieldVec, | ||
| PetscInt | iCell, | ||
| PetscInt | jCell, | ||
| PetscInt | kCell, | ||
| Cmpnts * | vec | ||
| ) |
Returns the vector value from the input cell index without blending.
This is a first-order piecewise-constant lookup helper used in workflows that intentionally disable trilinear blending.
| [in] | fieldName | Field label used for diagnostics. |
| [in] | fieldVec | Vector field array indexed as [k][j][i]. |
| [in] | iCell | Cell i-index. |
| [in] | jCell | Cell j-index. |
| [in] | kCell | Cell k-index. |
| [out] | vec | Output vector value. |
Returns the vector value from the input cell index without blending.
Local to this translation unit.
Definition at line 487 of file interpolation.c.
| PetscErrorCode TrilinearInterpolation_Scalar | ( | const char * | fieldName, |
| PetscReal *** | fieldScal, | ||
| PetscInt | i, | ||
| PetscInt | j, | ||
| PetscInt | k, | ||
| PetscReal | a1, | ||
| PetscReal | a2, | ||
| PetscReal | a3, | ||
| PetscReal * | val | ||
| ) |
Computes the trilinear interpolated scalar at a given point.
each cell a PetscReal. This function uses the standard 8-corner trilinear formula via ComputeTrilinearWeights(). If a different scheme is desired, implement a new function with the same interface.
| fieldName | A |
| fieldScal | 3D |
| i | Integral |
| j | Integral |
| k | Integral |
| a1 | Normalized |
| a2 | Normalized |
| a3 | Normalized |
| val | Pointer |
Computes the trilinear interpolated scalar at a given point.
Local to this translation unit.
Definition at line 577 of file interpolation.c.
| PetscErrorCode TrilinearInterpolation_Vector | ( | const char * | fieldName, |
| Cmpnts *** | fieldVec, | ||
| PetscInt | i, | ||
| PetscInt | j, | ||
| PetscInt | k, | ||
| PetscReal | a1, | ||
| PetscReal | a2, | ||
| PetscReal | a3, | ||
| Cmpnts * | vec | ||
| ) |
Computes the trilinear interpolated vector (e.g., velocity) at a given point.
each cell of type Cmpnts. This function uses the standard 8-corner trilinear formula via ComputeTrilinearWeights(). If a different scheme is desired, implement a new function with the same interface.
| fieldName | A |
| fieldVec | 3D |
| i | Integral |
| j | Integral |
| k | Integral |
| a1 | Normalized |
| a2 | Normalized |
| a3 | Normalized |
| vec | Pointer |
Computes the trilinear interpolated vector (e.g., velocity) at a given point.
Local to this translation unit.
Definition at line 641 of file interpolation.c.
| PetscErrorCode InterpolateEulerFieldToSwarm | ( | UserCtx * | user, |
| Vec | fieldLocal_cellCentered, | ||
| const char * | fieldName, | ||
| const char * | swarmOutFieldName | ||
| ) |
Interpolates a cell-centered field (scalar or vector) onto DMSwarm particles, using a robust, PETSc-idiomatic two-stage process.
This function first converts the cell-centered input data to corner-node data, storing this intermediate result in a PETSc Vec to correctly handle the communication of ghost-point information across parallel ranks. It then performs a final trilinear interpolation from the ghosted corner data to each particle's location.
Workflow:
cornerGlobal, cornerLocal) to manage the intermediate corner-node data, using the existing nodal DMDA (user->fda).bs) to select the correct underlying center-to-corner function, writing results into cornerGlobal.DMGlobalToLocal) to transfer the boundary data from cornerGlobal into the ghost regions of cornerLocal.cornerLocal data.| [in] | user | User context with DMDA, DMSwarm, etc. |
| [in] | fieldLocal_cellCentered | Ghosted local Vec with cell-centered data (e.g., Ucat). |
| [in] | fieldName | Human-readable field name for logging (e.g., "Ucat"). |
| [in] | swarmOutFieldName | Name of the DMSwarm field where interpolation results go. |
Interpolates a cell-centered field (scalar or vector) onto DMSwarm particles, using a robust, PETSc-idiomatic two-stage process.
Routes to InterpolateEulerFieldFromCenterToSwarm (direct trilinear, second-order) or InterpolateEulerFieldFromCornerToSwarm (corner-averaged, legacy) based on user->simCtx->interpolationMethod.
Definition at line 1371 of file interpolation.c.
| PetscErrorCode InterpolateAllFieldsToSwarm | ( | UserCtx * | user | ) |
Interpolates all relevant fields from the DMDA to the DMSwarm.
Currently, it interpolates:
To add more fields, duplicate the call to InterpolateOneFieldOverSwarm and provide:
| [in,out] | user | Pointer to a UserCtx containing:
|
Interpolates all relevant fields from the DMDA to the DMSwarm.
Local to this translation unit.
Definition at line 1397 of file interpolation.c.
| PetscErrorCode InterpolateParticleVelocities | ( | UserCtx * | user | ) |
Interpolates particle velocities using trilinear interpolation.
| [in] | user | Pointer to the user-defined context containing grid and swarm information. |
| PetscErrorCode InterpolateFieldFromCornerToCenter_Scalar | ( | PetscReal *** | field_arr, |
| PetscReal *** | centfield_arr, | ||
| UserCtx * | user | ||
| ) |
Safely interpolate a scalar field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx.
For each cell center in the physical region of the cell-centered DM (fda), this function averages the up to 8 surrounding scalar values from the coordinate DM (da). On boundaries, where fewer corners are available, a partial average is computed.
The coordinate DM (da) is built on corners (IM+1 x JM+1 x KM+1) while the cell-centered DM (fda) covers the physical cells (IM x JM x KM). Index offsets are adjusted via DMDAGetLocalInfo.
| [in] | field_arr | 3D array of corner-based scalar data (from user->da). |
| [out] | centfield_arr | 3D array for interpolated cell-center scalar data (for user->fda). |
| [in] | user | User context containing:
|
Safely interpolate a scalar field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx.
Local to this translation unit.
Definition at line 86 of file interpolation.c.
| PetscErrorCode InterpolateFieldFromCornerToCenter_Vector | ( | Cmpnts *** | field_arr, |
| Cmpnts *** | centfield_arr, | ||
| UserCtx * | user | ||
| ) |
Safely interpolate a vector field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx.
For each cell center in the physical region of the cell-centered DM (fda), this function averages the 8 surrounding corner values from the coordinate DM (da). The coordinate DM (da) is built on corners (IM+1 x JM+1 x KM+1) while the cell-centered DM (fda) covers the physical cells (IM x JM x KM). Index offsets are adjusted using DMDAGetLocalInfo.
| [in] | field_arr | 3D array of corner-based vector data (from user->da). |
| [out] | centfield_arr | 3D array for interpolated cell-center vector data (for user->fda). |
| [in] | user | User context containing:
|
Safely interpolate a vector field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx.
Local to this translation unit.
Definition at line 25 of file interpolation.c.
| PetscErrorCode TestCornerToCenterInterpolation | ( | UserCtx * | user | ) |
Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector.
This function serves as a unit test. It performs the following steps:
InterpolateFieldFromCornerToCenter macro to interpolate these coordinates to the cell centers, storing the result in a new temporary vector.user->Cent vector, which is assumed to have been computed by ComputeCellCentersAndSpacing and serves as the ground truth.ComputeCellCentersAndSpacing has been successfully executed.| user | The UserCtx for a specific grid level. |
Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector.
Local to this translation unit.
Definition at line 137 of file interpolation.c.
| PetscErrorCode InterpolateFieldFromCenterToCorner_Vector | ( | Cmpnts *** | centfield_arr, |
| Cmpnts *** | corner_arr, | ||
| UserCtx * | user | ||
| ) |
Interpolates a vector field from cell centers to corner nodes.
This version is adapted to write directly into a ghosted local array obtained from DMDAVecGetArray(), which allows using GLOBAL indices for writing to the OWNED portion of the array.
| [in] | centfield_arr | Input: 3D array (ghosted) of cell-centered data, accessed via GLOBAL indices. |
| [out] | corner_arr | Output: 3D array (ghosted) where interpolated node values are stored, also accessed via GLOBAL indices for the owned part. |
| [in] | user | User context containing DMDA information. |
Interpolates a vector field from cell centers to corner nodes.
Local to this translation unit.
Definition at line 196 of file interpolation.c.
| PetscErrorCode InterpolateFieldFromCenterToCorner_Scalar | ( | PetscReal *** | centfield_arr, |
| PetscReal *** | corner_arr, | ||
| UserCtx * | user | ||
| ) |
Interpolates a scalar field from cell centers to corner nodes.
This version is adapted to write directly into a ghosted local array obtained from DMDAVecGetArray(), which allows using GLOBAL indices for writing to the OWNED portion of the array.
| [in] | centfield_arr | Input: 3D array (ghosted) of scalar data at cell centers, accessed via GLOBAL indices. |
| [out] | corner_arr | Output: 3D array (ghosted) where interpolated node values are stored, also accessed via GLOBAL indices for the owned part. |
| [in] | user | User context containing DMDA information. |
Interpolates a scalar field from cell centers to corner nodes.
Local to this translation unit.
Definition at line 331 of file interpolation.c.
| PetscErrorCode InterpolateCornerToFaceCenter_Scalar | ( | PetscReal *** | corner_arr, |
| PetscReal *** | faceX_arr, | ||
| PetscReal *** | faceY_arr, | ||
| PetscReal *** | faceZ_arr, | ||
| UserCtx * | user | ||
| ) |
Interpolates a scalar field from corner nodes to all face centers.
This routine computes the average of the four corner-node values defining each face of a hexahedral cell:
| [in] | corner_arr | Ghosted node-centered array (global indexing) from user->fda. |
| [out] | faceX_arr | Local array for X-faces sized [zm][ym][xm+1]. |
| [out] | faceY_arr | Local array for Y-faces sized [zm][ym+1][xm]. |
| [out] | faceZ_arr | Local array for Z-faces sized [zm+1][ym][xm]. |
| [in] | user | User context containing DMDA 'fda' and GetOwnedCellRange. |
Interpolates a scalar field from corner nodes to all face centers.
Local to this translation unit.
Definition at line 2069 of file interpolation.c.
| PetscErrorCode InterpolateCornerToFaceCenter_Vector | ( | Cmpnts *** | corner_arr, |
| Cmpnts *** | faceX_arr, | ||
| Cmpnts *** | faceY_arr, | ||
| Cmpnts *** | faceZ_arr, | ||
| UserCtx * | user | ||
| ) |
Interpolates a vector field from corner nodes to all face centers.
Identical to the scalar version, except it averages each component of the Cmpnts struct at the four corner-nodes per face.
| [in] | corner_arr | Ghosted 3-component array (global node indices). |
| [out] | faceX_arr | Local array of Cmpnts for X-faces sized [zm][ym][xm+1]. |
| [out] | faceY_arr | Local array of Cmpnts for Y-faces sized [zm][ym+1][xm]. |
| [out] | faceZ_arr | Local array of Cmpnts for Z-faces sized [zm+1][ym][xm]. |
| [in] | user | User context containing DMDA 'fda'. |
Interpolates a vector field from corner nodes to all face centers.
Local to this translation unit.
Definition at line 2159 of file interpolation.c.