PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_postprocessor.c
Go to the documentation of this file.
1/**
2 * @file test_postprocessor.c
3 * @brief C test module for PICurv.
4 */
5
6#include "test_support.h"
7
8#include "postprocessor.h"
9
10#include <stdio.h>
11#include <string.h>
12/**
13 * @brief Test-local routine.
14 */
15
17{
18 SimCtx *simCtx = NULL;
19 UserCtx *user = NULL;
21 PetscInt bs = 0;
22 void *field_ptr = NULL;
23
24 PetscFunctionBeginUser;
25 PetscCall(PetscMemzero(&pps, sizeof(pps)));
26 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
27 PetscCall(PetscStrncpy(pps.particle_pipeline, "ComputeSpecificKE:velocity>ske", sizeof(pps.particle_pipeline)));
28
29 PetscCall(SetupPostProcessSwarm(user, &pps));
30 PetscCall(PicurvAssertBool((PetscBool)(user->post_swarm != NULL),
31 "SetupPostProcessSwarm should allocate a dedicated post swarm"));
32
33 PetscCall(DMSwarmSetLocalSizes(user->post_swarm, 0, 0));
34 PetscCall(DMSwarmGetField(user->post_swarm, "ske", &bs, NULL, &field_ptr));
35 PetscCall(PicurvAssertIntEqual(1, bs, "SetupPostProcessSwarm should register scalar output field 'ske'"));
36 PetscCall(DMSwarmRestoreField(user->post_swarm, "ske", &bs, NULL, &field_ptr));
37
38 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
39 PetscFunctionReturn(0);
40}
41/**
42 * @brief Test-local routine.
43 */
44
46{
47 SimCtx *simCtx = NULL;
48 UserCtx *user = NULL;
50 const PetscScalar ***p_nodal_arr = NULL;
51
52 PetscFunctionBeginUser;
53 PetscCall(PetscMemzero(&pps, sizeof(pps)));
54 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
55 PetscCall(PetscStrncpy(pps.process_pipeline, "UnknownStep;CellToNodeAverage:P>P_nodal", sizeof(pps.process_pipeline)));
56
57 PetscCall(VecSet(user->P, 9.0));
58 PetscCall(VecSet(user->P_nodal, -1.0));
59 PetscCall(EulerianDataProcessingPipeline(user, &pps));
60
61 PetscCall(DMDAVecGetArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
62 PetscCall(PicurvAssertRealNear(9.0, PetscRealPart(p_nodal_arr[0][0][0]), 1.0e-12,
63 "Eulerian pipeline should execute CellToNodeAverage for P->P_nodal"));
64 PetscCall(PicurvAssertRealNear(-1.0, PetscRealPart(p_nodal_arr[3][3][3]), 1.0e-12,
65 "Eulerian pipeline should leave non-computed boundary nodes unchanged"));
66 PetscCall(DMDAVecRestoreArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
67
68 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
69 PetscFunctionReturn(0);
70}
71/**
72 * @brief Test-local routine.
73 */
74
76{
77 SimCtx *simCtx = NULL;
78 UserCtx *user = NULL;
80 PetscReal (*vel_arr)[3] = NULL;
81 PetscScalar *ske_arr = NULL;
82
83 PetscFunctionBeginUser;
84 PetscCall(PetscMemzero(&pps, sizeof(pps)));
85 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
86 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
87 PetscCall(PetscStrncpy(pps.particle_pipeline, "ComputeSpecificKE:velocity>ske", sizeof(pps.particle_pipeline)));
88
89 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
90 vel_arr[0][0] = 1.0;
91 vel_arr[0][1] = 2.0;
92 vel_arr[0][2] = 2.0;
93 vel_arr[1][0] = 0.0;
94 vel_arr[1][1] = 3.0;
95 vel_arr[1][2] = 4.0;
96 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
97
98 PetscCall(ParticleDataProcessingPipeline(user, &pps));
99 PetscCall(DMSwarmGetField(user->post_swarm, "ske", NULL, NULL, (void *)&ske_arr));
100 PetscCall(PicurvAssertRealNear(4.5, PetscRealPart(ske_arr[0]), 1.0e-12,
101 "Particle pipeline should write first specific kinetic energy value"));
102 PetscCall(PicurvAssertRealNear(12.5, PetscRealPart(ske_arr[1]), 1.0e-12,
103 "Particle pipeline should write second specific kinetic energy value"));
104 PetscCall(DMSwarmRestoreField(user->post_swarm, "ske", NULL, NULL, (void *)&ske_arr));
105
106 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
107 PetscFunctionReturn(0);
108}
109/**
110 * @brief Test-local routine.
111 */
112
114{
115 SimCtx *simCtx = NULL;
116 UserCtx *user = NULL;
118 char tmpdir[PETSC_MAX_PATH_LEN];
119 char csv_prefix[PETSC_MAX_PATH_LEN];
120 char csv_path[PETSC_MAX_PATH_LEN];
121 PetscReal (*pos_arr)[3] = NULL;
122
123 PetscFunctionBeginUser;
124 PetscCall(PetscMemzero(&pps, sizeof(pps)));
125 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
126 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
127
128 simCtx->dt = 1.0;
129 simCtx->ren = 1.0;
130 simCtx->schmidt_number = 1.0;
131 simCtx->psrc_x = 0.0;
132 simCtx->psrc_y = 0.0;
133 simCtx->psrc_z = 0.0;
134
135 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
136 PetscCall(PetscSNPrintf(csv_prefix, sizeof(csv_prefix), "%s/stats", tmpdir));
137 PetscCall(PetscSNPrintf(csv_path, sizeof(csv_path), "%s_msd.csv", csv_prefix));
138
139 PetscCall(PetscStrncpy(pps.statistics_pipeline, "ComputeMSD", sizeof(pps.statistics_pipeline)));
140 PetscCall(PetscStrncpy(pps.statistics_output_prefix, csv_prefix, sizeof(pps.statistics_output_prefix)));
141
142 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
143 pos_arr[0][0] = 1.0;
144 pos_arr[0][1] = 0.0;
145 pos_arr[0][2] = 0.0;
146 pos_arr[1][0] = -1.0;
147 pos_arr[1][1] = 0.0;
148 pos_arr[1][2] = 0.0;
149 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
150
151 PetscCall(GlobalStatisticsPipeline(user, &pps, 1));
152 PetscCall(PicurvAssertFileExists(csv_path, "GlobalStatisticsPipeline should dispatch ComputeMSD and emit CSV output"));
153
154 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
155 PetscFunctionReturn(0);
156}
157/**
158 * @brief Test-local routine.
159 */
160
161static PetscErrorCode TestWriteEulerianFileWritesVTS(void)
162{
163 SimCtx *simCtx = NULL;
164 UserCtx *user = NULL;
166 char tmpdir[PETSC_MAX_PATH_LEN];
167 char vtk_path[PETSC_MAX_PATH_LEN];
168
169 PetscFunctionBeginUser;
170 PetscCall(PetscMemzero(&pps, sizeof(pps)));
171 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
172 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
173 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/field_00003.vts", tmpdir));
174
175 PetscCall(PetscStrncpy(pps.output_fields_instantaneous, "P_nodal", sizeof(pps.output_fields_instantaneous)));
176 PetscCall(PetscSNPrintf(pps.output_prefix, sizeof(pps.output_prefix), "%s/field", tmpdir));
177 PetscCall(VecSet(user->P_nodal, 7.0));
178
179 PetscCall(WriteEulerianFile(user, &pps, 3));
180 PetscCall(PicurvAssertFileExists(vtk_path, "WriteEulerianFile should emit a .vts file for requested output fields"));
181
182 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
183 PetscFunctionReturn(0);
184}
185/**
186 * @brief Test-local routine.
187 */
188
189static PetscErrorCode TestWriteParticleFileWritesVTP(void)
190{
191 SimCtx *simCtx = NULL;
192 UserCtx *user = NULL;
194 char tmpdir[PETSC_MAX_PATH_LEN];
195 char vtk_path[PETSC_MAX_PATH_LEN];
196 PetscReal (*pos_arr)[3] = NULL;
197 PetscReal (*vel_arr)[3] = NULL;
198
199 PetscFunctionBeginUser;
200 PetscCall(PetscMemzero(&pps, sizeof(pps)));
201 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
202 PetscCall(PicurvCreateSwarmPair(user, 3, "ske"));
203 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
204 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/particles_00004.vtp", tmpdir));
205
206 pps.outputParticles = PETSC_TRUE;
207 pps.particle_output_freq = 1;
208 PetscCall(PetscStrncpy(pps.particle_fields, "position,velocity", sizeof(pps.particle_fields)));
209 PetscCall(PetscSNPrintf(pps.particle_output_prefix, sizeof(pps.particle_output_prefix), "%s/particles", tmpdir));
210
211 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
212 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
213 for (PetscInt p = 0; p < 3; ++p) {
214 pos_arr[p][0] = (PetscReal)p;
215 pos_arr[p][1] = (PetscReal)(p + 1);
216 pos_arr[p][2] = (PetscReal)(p + 2);
217 vel_arr[p][0] = 2.0 * (PetscReal)p;
218 vel_arr[p][1] = 3.0 * (PetscReal)p;
219 vel_arr[p][2] = 4.0 * (PetscReal)p;
220 }
221 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
222 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
223
224 PetscCall(WriteParticleFile(user, &pps, 4));
225 PetscCall(PicurvAssertFileExists(vtk_path, "WriteParticleFile should emit a .vtp file for requested particle fields"));
226
227 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
228 PetscFunctionReturn(0);
229}
230/**
231 * @brief Entry point for this unit-test binary.
232 */
233
234int main(int argc, char **argv)
235{
236 PetscErrorCode ierr;
237 const PicurvTestCase cases[] = {
238 {"setup-postprocess-swarm-registers-pipeline-fields", TestSetupPostProcessSwarmRegistersPipelineFields},
239 {"eulerian-data-processing-pipeline-runs-configured-kernels", TestEulerianDataProcessingPipelineRunsConfiguredKernels},
240 {"particle-data-processing-pipeline-computes-specific-ke", TestParticleDataProcessingPipelineComputesSpecificKE},
241 {"global-statistics-pipeline-writes-msd-csv", TestGlobalStatisticsPipelineWritesMSDCSV},
242 {"write-eulerian-file-writes-vts", TestWriteEulerianFileWritesVTS},
243 {"write-particle-file-writes-vtp", TestWriteParticleFileWritesVTP},
244 };
245
246 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv postprocessor tests");
247 if (ierr) {
248 return (int)ierr;
249 }
250
251 ierr = PicurvRunTests("unit-postprocessor", cases, sizeof(cases) / sizeof(cases[0]));
252 if (ierr) {
253 PetscFinalize();
254 return (int)ierr;
255 }
256
257 ierr = PetscFinalize();
258 return (int)ierr;
259}
PetscErrorCode EulerianDataProcessingPipeline(UserCtx *user, PostProcessParams *pps)
Parses the processing pipeline string and executes the requested kernels.
PetscErrorCode WriteEulerianFile(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Orchestrates the writing of a combined, multi-field VTK file for a single time step.
PetscErrorCode GlobalStatisticsPipeline(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Executes the global statistics pipeline, computing aggregate reductions over all particles.
PetscErrorCode ParticleDataProcessingPipeline(UserCtx *user, PostProcessParams *pps)
Parses and executes the particle pipeline using a robust two-pass approach.
PetscErrorCode WriteParticleFile(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Writes particle data to a VTP file using the Prepare-Write-Cleanup pattern.
PetscErrorCode SetupPostProcessSwarm(UserCtx *user, PostProcessParams *pps)
Creates a new, dedicated DMSwarm for post-processing tasks.
static PetscErrorCode TestWriteParticleFileWritesVTP(void)
Test-local routine.
int main(int argc, char **argv)
Entry point for this unit-test binary.
static PetscErrorCode TestParticleDataProcessingPipelineComputesSpecificKE(void)
Test-local routine.
static PetscErrorCode TestGlobalStatisticsPipelineWritesMSDCSV(void)
Test-local routine.
static PetscErrorCode TestWriteEulerianFileWritesVTS(void)
Test-local routine.
static PetscErrorCode TestEulerianDataProcessingPipelineRunsConfiguredKernels(void)
Test-local routine.
static PetscErrorCode TestSetupPostProcessSwarmRegistersPipelineFields(void)
Test-local routine.
PetscErrorCode PicurvMakeTempDir(char *path, size_t path_len)
Shared test-support routine.
PetscErrorCode PicurvCreateMinimalContexts(SimCtx **simCtx_out, UserCtx **user_out, PetscInt mx, PetscInt my, PetscInt mz)
Shared test-support routine.
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Shared test-support routine.
PetscErrorCode PicurvDestroyMinimalContexts(SimCtx **simCtx_ptr, UserCtx **user_ptr)
Shared test-support routine.
PetscErrorCode PicurvCreateSwarmPair(UserCtx *user, PetscInt nlocal, const char *post_field_name)
Shared test-support routine.
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Shared test-support routine.
PetscErrorCode PicurvAssertFileExists(const char *path, const char *context)
Shared test-support routine.
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Shared test-support routine.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Shared test-support routine.
C test module for PICurv.
Named test case descriptor consumed by PicurvRunTests.
char statistics_output_prefix[256]
basename for CSV output, e.g.
Definition variables.h:523
Vec P_nodal
Definition variables.h:814
char particle_output_prefix[256]
Definition variables.h:518
PetscReal schmidt_number
Definition variables.h:653
char output_prefix[256]
Definition variables.h:515
PetscReal psrc_x
Definition variables.h:650
PetscReal ren
Definition variables.h:638
DM post_swarm
Definition variables.h:813
PetscReal dt
Definition variables.h:606
char statistics_pipeline[1024]
e.g.
Definition variables.h:522
PetscReal psrc_z
Point source location for PARTICLE_INIT_POINT_SOURCE.
Definition variables.h:650
char output_fields_instantaneous[1024]
Definition variables.h:513
char particle_pipeline[1024]
Definition variables.h:516
PetscReal psrc_y
Definition variables.h:650
char process_pipeline[1024]
Definition variables.h:512
PetscInt particle_output_freq
Definition variables.h:519
char particle_fields[1024]
Definition variables.h:517
PetscBool outputParticles
Definition variables.h:509
Holds all configuration parameters for a post-processing run.
Definition variables.h:501
The master context for the entire simulation.
Definition variables.h:591
User-defined context containing data specific to a single computational grid level.
Definition variables.h:738