PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
|
Go to the source code of this file.
Macros | |
#define | BBOX_TOLERANCE 1e-6 |
#define | __FUNCT__ "InitializeSingleGridDM" |
#define | __FUNCT__ "InitializeAllGridDMs" |
#define | __FUNCT__ "AssignAllGridCoordinates" |
#define | __FUNCT__ "SetFinestLevelCoordinates" |
#define | __FUNCT__ "GenerateAndSetCoordinates" |
#define | __FUNCT__ "ReadAndSetCoordinates" |
#define | __FUNCT__ "RestrictCoordinates" |
Functions | |
static PetscErrorCode | ParseAndSetGridInputs (UserCtx *user) |
Determines the grid source and calls the appropriate parsing routine. | |
PetscErrorCode | DefineAllGridDimensions (SimCtx *simCtx) |
Orchestrates the parsing and setting of grid dimensions for all blocks. | |
static PetscErrorCode | InitializeSingleGridDM (UserCtx *user, UserCtx *coarse_user) |
Creates the DMDA objects (da and fda) for a single UserCtx. | |
PetscErrorCode | InitializeAllGridDMs (SimCtx *simCtx) |
Orchestrates the creation of DMDA objects for every block and multigrid level. | |
static PetscErrorCode | SetFinestLevelCoordinates (UserCtx *user) |
A router that populates the coordinates for a single finest-level DMDA. | |
static PetscErrorCode | GenerateAndSetCoordinates (UserCtx *user) |
Programmatically generates and sets grid coordinates based on user parameters. | |
static PetscErrorCode | ReadAndSetCoordinates (UserCtx *user, FILE *fd) |
Reads physical coordinates from a file and populates the DMDA for a specific block. | |
static PetscErrorCode | RestrictCoordinates (UserCtx *coarse_user, UserCtx *fine_user) |
Populates coarse grid coordinates by restricting from a fine grid. | |
PetscErrorCode | AssignAllGridCoordinates (SimCtx *simCtx) |
Orchestrates the assignment of physical coordinates to all DMDA objects. | |
static PetscReal | ComputeStretchedCoord (PetscInt i, PetscInt N, PetscReal L, PetscReal r) |
Computes a stretched coordinate along one dimension. | |
PetscErrorCode | ComputeLocalBoundingBox (UserCtx *user, BoundingBox *localBBox) |
Computes the local bounding box of the grid on the current process. | |
PetscErrorCode | GatherAllBoundingBoxes (UserCtx *user, BoundingBox **allBBoxes) |
Gathers local bounding boxes from all MPI processes to rank 0. | |
PetscErrorCode | BroadcastAllBoundingBoxes (UserCtx *user, BoundingBox **bboxlist) |
Broadcasts the bounding box list from rank 0 to all ranks. | |
|
static |
Determines the grid source and calls the appropriate parsing routine.
This function acts as a router. It checks the global simCtx->generate_grid
flag (accessed via the user->simCtx
back-pointer) to decide whether to call the parser for a programmatically generated grid or for a grid defined in a file.
user | Pointer to the UserCtx for a specific block. The function will populate the geometric fields within this struct. |
Definition at line 20 of file grid.c.
PetscErrorCode DefineAllGridDimensions | ( | SimCtx * | simCtx | ) |
Orchestrates the parsing and setting of grid dimensions for all blocks.
This function serves as the high-level entry point for defining the geometric properties of each grid block in the simulation. It iterates through every block defined by simCtx->block_number
.
For each block, it performs two key actions:
_this
) in the corresponding UserCtx
struct for the finest multigrid level. This makes the context "self-aware".ParseAndSetGridInputs
) to handle the detailed work of parsing options or files to populate the rest of the geometric properties for that specific block (e.g., IM
, Min_X
, rx
).simCtx | The master SimCtx, which contains the number of blocks and the UserCtx hierarchy to be configured. |
Definition at line 57 of file grid.c.
Creates the DMDA objects (da and fda) for a single UserCtx.
This function is a direct adaptation of the core logic in MGDACreate
. It creates the scalar (da
) and vector (fda
) DMs for a single grid level.
If a coarse_user
context is provided, it performs the critical processor alignment calculation from the legacy code. This ensures the new (fine) DM aligns with the coarse DM for multigrid efficiency. If coarse_user
is NULL, it creates the DM with a default PETSc decomposition, intended for the coarsest grid level.
user | The UserCtx for which the DMs will be created. Its IM, JM, KM fields must be pre-populated. |
coarse_user | The UserCtx of the next-coarser grid level, or NULL if user is the coarsest level. |
Definition at line 107 of file grid.c.
PetscErrorCode InitializeAllGridDMs | ( | SimCtx * | simCtx | ) |
Orchestrates the creation of DMDA objects for every block and multigrid level.
This function systematically builds the entire DMDA hierarchy. It first calculates the dimensions (IM, JM, KM) for all coarse grids based on the finest grid's dimensions and the semi-coarsening flags. It then iterates from the coarsest to the finest level, calling a powerful helper function (InitializeSingleGridDM
) to create the DMs for each block, ensuring that finer grids are properly aligned with their coarser parents for multigrid efficiency.
Definition at line 230 of file grid.c.
|
static |
A router that populates the coordinates for a single finest-level DMDA.
This function orchestrates the coordinate setting for one block. It checks the global generate_grid
flag and calls the appropriate helper for either programmatic generation or reading from a file.
After the local coordinate vector is populated by a helper, this function performs the necessary DMLocalToGlobal and DMGlobalToLocal scatters to ensure that the ghost node coordinate values are correctly communicated and updated across all MPI ranks.
user | The UserCtx for a specific block on the finest level. |
Definition at line 367 of file grid.c.
|
static |
Programmatically generates and sets grid coordinates based on user parameters.
This function populates the local coordinate vector of the provided UserCtx
using the geometric properties (Min_X
, Max_X
, rx
, etc.) that were parsed from command-line options. It supports non-linear grid stretching.
user | The UserCtx for a specific block. |
Definition at line 476 of file grid.c.
|
static |
Reads physical coordinates from a file and populates the DMDA for a specific block.
This function handles the collective read of an interleaved (X Y Z per line) multi-block grid file. It assumes the file header (nblk, dimensions) has already been processed by ReadGridFile.
The process is robust for parallel execution:
user->_this
).gc
.gc
is broadcast to all other MPI ranks.gc
array, and copies the (x,y,z) coordinates into its local PETSc coordinate vector.user | The UserCtx for a specific block. Its _this field must be set, and its IM, JM, KM fields must be correctly pre-populated. |
Definition at line 535 of file grid.c.
Populates coarse grid coordinates by restricting from a fine grid.
This function is a direct adaptation of the coordinate restriction loop from the legacy MG_Initial
function. It ensures that the physical location of a coarse grid node is identical to its corresponding parent node on the fine grid. The mapping from coarse to fine index (i
-> ih
) is determined by the semi-coarsening flags (isc
, jsc
, ksc
) stored in the UserCtx
.
coarse_user | The UserCtx for the coarse grid (destination). |
fine_user | The UserCtx for the fine grid (source). |
Definition at line 617 of file grid.c.
PetscErrorCode AssignAllGridCoordinates | ( | SimCtx * | simCtx | ) |
Orchestrates the assignment of physical coordinates to all DMDA objects.
This function manages the entire process of populating the coordinate vectors for every DMDA across all multigrid levels and blocks. It follows a two-part strategy that is essential for multigrid methods:
SetFinestLevelCoordinates
) to set the physical coordinates for the highest-resolution grid (the finest multigrid level).RestrictCoordinates
) to copy the coordinate values from the fine grid nodes to their corresponding parent nodes on the coarser grids. This ensures all levels represent the exact same geometry.Definition at line 317 of file grid.c.
|
inlinestatic |
Computes a stretched coordinate along one dimension.
This function computes a coordinate based on a geometric stretching ratio. If the ratio (r) is 1.0, a uniform distribution is used: x(i) = L * (i/N)
If r != 1.0, a geometric stretching is applied: x(i) = L * [ (r^(i/N) - 1 ) / (r - 1) ]
Here:
[in] | i | Current index (0 <= i <= N). |
[in] | N | Number of segments along the dimension. |
[in] | L | Total length of the domain. |
[in] | r | Stretching ratio. |
Definition at line 453 of file grid.c.
PetscErrorCode ComputeLocalBoundingBox | ( | UserCtx * | user, |
BoundingBox * | localBBox | ||
) |
Computes the local bounding box of the grid on the current process.
This function calculates the minimum and maximum coordinates (x, y, z) of the local grid points owned by the current MPI process. It iterates over the local portion of the grid, examines each grid point's coordinates, and updates the minimum and maximum values accordingly.
The computed bounding box is stored in the provided localBBox
structure, and the user->bbox
field is also updated with this bounding box for consistency within the user context.
[in] | user | Pointer to the user-defined context containing grid information. This context must be properly initialized before calling this function. |
[out] | localBBox | Pointer to the BoundingBox structure where the computed local bounding box will be stored. The structure should be allocated by the caller. |
0
on success, non-zero on failure. Definition at line 702 of file grid.c.
PetscErrorCode GatherAllBoundingBoxes | ( | UserCtx * | user, |
BoundingBox ** | allBBoxes | ||
) |
Gathers local bounding boxes from all MPI processes to rank 0.
Each rank computes its local bounding box, then all ranks participate in an MPI_Gather to send their BoundingBox to rank 0. Rank 0 allocates the result array and returns it via allBBoxes.
[in] | user | Pointer to UserCtx (must be non-NULL). |
[out] | allBBoxes | On rank 0, receives malloc’d array of size size . On other ranks, set to NULL. |
Definition at line 837 of file grid.c.
PetscErrorCode BroadcastAllBoundingBoxes | ( | UserCtx * | user, |
BoundingBox ** | bboxlist | ||
) |
Broadcasts the bounding box list from rank 0 to all ranks.
Broadcasts the bounding box information collected on rank 0 to all other ranks.
After GatherAllBoundingBoxes, rank 0 has an array of size
boxes. This routine makes sure every rank ends up with its own malloc’d copy.
[in] | user | Pointer to UserCtx (unused here, but kept for signature). |
[in,out] | bboxlist | On entry: rank 0’s array; on exit: every rank’s array. |
Definition at line 896 of file grid.c.