PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
particle_statistics.h
Go to the documentation of this file.
1#ifndef PARTICLE_STATISTICS_H
2#define PARTICLE_STATISTICS_H
3
4#include "variables.h"
5#include "logging.h"
6
7/**
8 * @file particle_statistics.h
9 * @brief Global statistics kernels for the Statistics Pipeline.
10 *
11 * These kernels compute global aggregate quantities via MPI reduction — they do NOT
12 * produce per-particle VTK data. Each kernel appends one CSV row per call and logs
13 * a one-line summary via LOG_INFO.
14 *
15 * Output file convention: {stats_prefix}_{kernel_name}.csv
16 * All MPI_Allreduce operations are internal; file I/O is rank-0-only.
17 *
18 * To add a new statistic:
19 * 1. Implement PetscErrorCode ComputeXxx(UserCtx*, const char*, PetscInt) in
20 * src/particle_statistics.c
21 * 2. Declare it below.
22 * 3. Add one else-if in GlobalStatisticsPipeline() in src/postprocessor.c
23 * No other files need to change.
24 */
25
26/**
27 * @brief Computes the mean-squared displacement (MSD) of a particle cloud.
28 *
29 * Reference point r0 = (simCtx->psrc_x, psrc_y, psrc_z).
30 * Diffusivity D = 1 / (Re * Sc), time t = ti * dt.
31 * Computes MSD_x, MSD_y, MSD_z (isotropy check), MSD_total, r_rms,
32 * centre-of-mass drift, and fractions inside 1σ/2σ/3σ theoretical shells.
33 * Appends one row to {stats_prefix}_msd.csv; writes header on first call.
34 *
35 * @param user The UserCtx containing the DMSwarm (user->swarm).
36 * @param stats_prefix Base filename prefix (e.g. "brownian_stats").
37 * @param ti Current time-step index.
38 * @return PetscErrorCode
39 */
40PetscErrorCode ComputeParticleMSD(UserCtx *user, const char *stats_prefix, PetscInt ti);
41
42/* -----------------------------------------------------------------------
43 * Future kernels — add declarations here:
44 *
45 * PetscErrorCode ComputeVelocityPDF(UserCtx *user, const char *stats_prefix, PetscInt ti);
46 * PetscErrorCode ComputeConcentrationField(UserCtx *user, const char *stats_prefix, PetscInt ti);
47 * PetscErrorCode ComputeResidenceTimeDistribution(UserCtx *user, const char *stats_prefix, PetscInt ti);
48 * ----------------------------------------------------------------------- */
49
50#endif /* PARTICLE_STATISTICS_H */
Logging utilities and macros for PETSc-based applications.
PetscErrorCode ComputeParticleMSD(UserCtx *user, const char *stats_prefix, PetscInt ti)
Computes the mean-squared displacement (MSD) of a particle cloud.
Main header file for a complex fluid dynamics solver.
User-defined context containing data specific to a single computational grid level.
Definition variables.h:732