PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
|
Lower-level functions used by the main scattering routines. More...
Functions | |
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 | 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. | |
static PetscErrorCode | ScatterParticleFieldToEulerField_Internal (UserCtx *user, const char *particleFieldName, DM targetDM, PetscInt expected_dof, Vec eulerFieldAverageVec) |
Internal helper function to orchestrate the scatter operation (accumulate + normalize). | |
Lower-level functions used by the main scattering routines.
PetscErrorCode AccumulateParticleField | ( | DM | swarm, |
const char * | particleFieldName, | ||
DM | gridSumDM, | ||
Vec | gridSumVec | ||
) |
Accumulates a particle field (scalar or vector) into a target grid sum vector.
This function iterates through local particles, identifies their cell using the "DMSwarm_CellID" field, and adds the particle's field value (particleFieldName
) to the corresponding cell location in the gridSumVec
. It handles both scalar (DOF=1) and vector (DOF=3) fields automatically based on the DOF of gridSumDM
.
IMPORTANT: The caller must ensure gridSumVec
is zeroed before calling this function if a fresh sum calculation is desired.
[in] | swarm | The DMSwarm containing particles. |
[in] | particleFieldName | Name of the field on the particles (must match DOF). |
[in] | gridSumDM | The DMDA associated with gridSumVec . Its DOF determines how many components are accumulated. |
[in,out] | gridSumVec | The Vec (associated with gridSumDM ) to accumulate sums into. |
This function iterates through local particles, identifies their cell using the "DMSwarm_CellID"
field, and adds the particle's field value (particleFieldName
) to the corresponding cell location in the gridSumVec
. It handles both scalar (DOF=1) and vector (DOF=3) fields automatically based on the DOF of gridSumDM
.
IMPORTANT: The caller must ensure gridSumVec
is zeroed before calling this function if a fresh sum calculation is desired.
[in] | swarm | The DMSwarm containing particles. |
[in] | particleFieldName | Name of the field on the particles (must match DOF of gridSumDM ). |
[in] | gridSumDM | The DMDA associated with gridSumVec . Its DOF determines how many components are accumulated. |
[in,out] | gridSumVec | The Vec (associated with gridSumDM ) to accumulate sums into. |
Definition at line 1849 of file interpolation.c.
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.
Calculates avgVec[i] = sumVec[i] / countVec[i] for each component of each OWNED cell where countVec[i] > 0. Sets avgVec[i] = 0 otherwise. Handles both scalar (DOF=1) and vector (DOF=3) data fields based on dataDM
. Uses basic VecGetArray
/VecGetArrayRead
and manual index calculation.
[in] | countDM | The DMDA associated with countVec (must have DOF=1). |
[in] | countVec | The Vec containing particle counts per cell (read-only). |
[in] | dataDM | The DMDA associated with sumVec and avgVec (must have DOF=1 or DOF=3). |
[in] | sumVec | The Vec containing the accumulated sums per cell (read-only). |
[in,out] | avgVec | The Vec where the calculated averages will be stored (overwritten). |
Calculates avgVec[i] = sumVec[i] / countVec[i]
for each component of each cell owned by the current process, provided countVec[i] > 0
. Otherwise, sets avgVec[i] = 0
. Handles both scalar (DOF=1) and vector (DOF=3) data fields based on the DOF of dataDM
. Uses DMDA multi-dimensional array accessors (DMDAVecGetArray...
) for safe and convenient indexing.
[in] | countDM | The DMDA associated with countVec (must have DOF=1). |
[in] | countVec | The Vec containing particle counts per cell (read-only). |
[in] | dataDM | The DMDA associated with sumVec and avgVec (must have DOF=1 or DOF=3). |
[in] | sumVec | The Vec containing the accumulated sums per cell (read-only). |
[in,out] | avgVec | The Vec where the calculated averages will be stored (overwritten). Must be associated with dataDM . |
Definition at line 1954 of file interpolation.c.
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.
Based on hardcoded rules mapping particle field names ("P", "Nvert", "Ucat", "Ucont") to user context DMs (user->da
or user->fda
). This function encapsulates the policy of where different fields should be scattered. Modify this function to add rules for custom fields.
[in] | user | Pointer to the UserCtx containing required DMs (da , fda ). |
[in] | particleFieldName | Name of the particle field. |
[out] | targetDM | Pointer to store the determined target DM (user->da or user->fda ). |
[out] | expected_dof | Pointer to store the expected DOF (1 or 3) for this field. |
PETSC_ERR_ARG_NULL
if required inputs are NULL.PETSC_ERR_ARG_WRONG
if particleFieldName
is not recognized. Definition at line 1784 of file interpolation.c.
|
static |
Internal helper function to orchestrate the scatter operation (accumulate + normalize).
Manages the temporary sum vector and calls the accumulation and normalization functions. Assumes caller determined target DM and DOF. Checks for particle field existence. NOTE: If DMSwarmGetField fails, the subsequent SETERRQ will overwrite the original error.
[in] | user | Pointer to UserCtx containing swarm, ParticleCount, da. |
[in] | particleFieldName | Name of the field in the DMSwarm. |
[in] | targetDM | The DMDA where the final average and intermediate sum reside. |
[in] | expected_dof | The expected DOF (1 or 3) for the targetDM and field. |
[in,out] | eulerFieldAverageVec | The pre-created Vec associated with targetDM to store the result. |
Definition at line 2074 of file interpolation.c.