PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Functions
test_postprocessor.c File Reference

C unit tests for postprocessor orchestration and output paths. More...

#include "test_support.h"
#include "postprocessor.h"
#include "setup.h"
#include <stdio.h>
#include <string.h>
Include dependency graph for test_postprocessor.c:

Go to the source code of this file.

Functions

static PetscErrorCode PicurvAssertFileMissing (const char *path, const char *context)
 Assert that a file path does not exist.
 
static PetscErrorCode CountLinesWithPrefix (const char *path, const char *prefix, PetscInt *count_out)
 Count file lines that begin with a given prefix.
 
static PetscErrorCode TestSetupPostProcessSwarmRegistersPipelineFields (void)
 Tests post-processing swarm setup and pipeline-field registration.
 
static PetscErrorCode TestEulerianDataProcessingPipelineRunsConfiguredKernels (void)
 Tests Eulerian post-processing pipeline dispatch for configured kernels.
 
static PetscErrorCode TestParticleDataProcessingPipelineComputesSpecificKE (void)
 Tests particle post-processing pipeline dispatch for specific-KE output.
 
static PetscErrorCode TestGlobalStatisticsPipelineWritesMSDCSV (void)
 Tests global statistics pipeline generation of MSD CSV output.
 
static PetscErrorCode TestGlobalStatisticsPipelineRewritesExistingMSDStep (void)
 Tests that reprocessing the same MSD step rewrites the CSV row instead of duplicating it.
 
static PetscErrorCode TestSetupSimulationEnvironmentCreatesStatisticsDirectoryWithoutClobbering (void)
 Tests that postprocessor setup creates nested statistics directories without clobbering existing contents.
 
static PetscErrorCode TestWriteEulerianFileWritesVTS (void)
 Tests VTS emission for Eulerian post-processing output.
 
static PetscErrorCode TestWriteParticleFileWritesVTP (void)
 Tests VTP emission for particle post-processing output.
 
static PetscErrorCode TestWriteEulerianFileRewritesSameStepCleanly (void)
 Tests repeated same-step VTS emission leaves only the final artifact.
 
static PetscErrorCode TestWriteParticleFileRewritesSameStepCleanly (void)
 Tests repeated same-step VTP emission leaves only the final artifact.
 
int main (int argc, char **argv)
 Runs the unit-postprocessor PETSc test binary.
 

Detailed Description

C unit tests for postprocessor orchestration and output paths.

Definition in file test_postprocessor.c.

Function Documentation

◆ PicurvAssertFileMissing()

static PetscErrorCode PicurvAssertFileMissing ( const char *  path,
const char *  context 
)
static

Assert that a file path does not exist.

Parameters
[in]pathFile path to check.
[in]contextFailure context message.
Returns
Petsc error code.

Definition at line 20 of file test_postprocessor.c.

21{
22 PetscBool exists = PETSC_FALSE;
23
24 PetscFunctionBeginUser;
25 PetscCall(PetscTestFile(path, 'r', &exists));
26 if (exists) {
27 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[FAIL] %s | unexpected file: %s\n", context, path));
28 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Unexpected file is present.");
29 }
30 PetscFunctionReturn(0);
31}
Here is the caller graph for this function:

◆ CountLinesWithPrefix()

static PetscErrorCode CountLinesWithPrefix ( const char *  path,
const char *  prefix,
PetscInt *  count_out 
)
static

Count file lines that begin with a given prefix.

Parameters
[in]pathFile path to scan.
[in]prefixPrefix to match at the beginning of each line.
[out]count_outCount of matching lines.
Returns
Petsc error code.

Definition at line 40 of file test_postprocessor.c.

41{
42 FILE *file = NULL;
43 char line[4096];
44 PetscInt count = 0;
45
46 PetscFunctionBeginUser;
47 PetscCheck(count_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "count_out must not be NULL.");
48 file = fopen(path, "r");
49 PetscCheck(file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open %s.", path);
50 while (fgets(line, sizeof(line), file)) {
51 if (strncmp(line, prefix, strlen(prefix)) == 0) count++;
52 }
53 fclose(file);
54 *count_out = count;
55 PetscFunctionReturn(0);
56}
Here is the caller graph for this function:

◆ TestSetupPostProcessSwarmRegistersPipelineFields()

static PetscErrorCode TestSetupPostProcessSwarmRegistersPipelineFields ( void  )
static

Tests post-processing swarm setup and pipeline-field registration.

Definition at line 61 of file test_postprocessor.c.

62{
63 SimCtx *simCtx = NULL;
64 UserCtx *user = NULL;
66 PetscInt bs = 0;
67 void *field_ptr = NULL;
68
69 PetscFunctionBeginUser;
70 PetscCall(PetscMemzero(&pps, sizeof(pps)));
71 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
72 PetscCall(PetscStrncpy(pps.particle_pipeline, "ComputeSpecificKE:velocity>ske", sizeof(pps.particle_pipeline)));
73
74 PetscCall(SetupPostProcessSwarm(user, &pps));
75 PetscCall(PicurvAssertBool((PetscBool)(user->post_swarm != NULL),
76 "SetupPostProcessSwarm should allocate a dedicated post swarm"));
77
78 PetscCall(DMSwarmSetLocalSizes(user->post_swarm, 0, 0));
79 PetscCall(DMSwarmGetField(user->post_swarm, "ske", &bs, NULL, &field_ptr));
80 PetscCall(PicurvAssertIntEqual(1, bs, "SetupPostProcessSwarm should register scalar output field 'ske'"));
81 PetscCall(DMSwarmRestoreField(user->post_swarm, "ske", &bs, NULL, &field_ptr));
82
83 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
84 PetscFunctionReturn(0);
85}
PetscErrorCode SetupPostProcessSwarm(UserCtx *user, PostProcessParams *pps)
Creates a new, dedicated DMSwarm for post-processing tasks.
PetscErrorCode PicurvCreateMinimalContexts(SimCtx **simCtx_out, UserCtx **user_out, PetscInt mx, PetscInt my, PetscInt mz)
Builds minimal SimCtx and UserCtx fixtures for C unit tests.
PetscErrorCode PicurvDestroyMinimalContexts(SimCtx **simCtx_ptr, UserCtx **user_ptr)
Destroys minimal SimCtx/UserCtx fixtures and all owned PETSc objects.
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Asserts that two integer values are equal.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Asserts that one boolean condition is true.
DM post_swarm
Definition variables.h:886
char particle_pipeline[1024]
Definition variables.h:568
Holds all configuration parameters for a post-processing run.
Definition variables.h:553
The master context for the entire simulation.
Definition variables.h:643
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestEulerianDataProcessingPipelineRunsConfiguredKernels()

static PetscErrorCode TestEulerianDataProcessingPipelineRunsConfiguredKernels ( void  )
static

Tests Eulerian post-processing pipeline dispatch for configured kernels.

Definition at line 90 of file test_postprocessor.c.

91{
92 SimCtx *simCtx = NULL;
93 UserCtx *user = NULL;
95 const PetscScalar ***p_nodal_arr = NULL;
96
97 PetscFunctionBeginUser;
98 PetscCall(PetscMemzero(&pps, sizeof(pps)));
99 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
100 PetscCall(PetscStrncpy(pps.process_pipeline, "UnknownStep;CellToNodeAverage:P>P_nodal", sizeof(pps.process_pipeline)));
101
102 PetscCall(VecSet(user->P, 9.0));
103 PetscCall(VecSet(user->P_nodal, -1.0));
104 PetscCall(EulerianDataProcessingPipeline(user, &pps));
105
106 PetscCall(DMDAVecGetArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
107 PetscCall(PicurvAssertRealNear(9.0, PetscRealPart(p_nodal_arr[0][0][0]), 1.0e-12,
108 "Eulerian pipeline should execute CellToNodeAverage for P->P_nodal"));
109 PetscCall(PicurvAssertRealNear(-1.0, PetscRealPart(p_nodal_arr[user->KM][user->JM][user->IM]), 1.0e-12,
110 "Eulerian pipeline should leave the extra non-physical boundary node unchanged"));
111 PetscCall(DMDAVecRestoreArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
112
113 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
114 PetscFunctionReturn(0);
115}
PetscErrorCode EulerianDataProcessingPipeline(UserCtx *user, PostProcessParams *pps)
Parses the processing pipeline string and executes the requested kernels.
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Asserts that two real values agree within tolerance.
Vec P_nodal
Definition variables.h:887
PetscInt KM
Definition variables.h:820
PetscInt JM
Definition variables.h:820
char process_pipeline[1024]
Definition variables.h:564
PetscInt IM
Definition variables.h:820
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestParticleDataProcessingPipelineComputesSpecificKE()

static PetscErrorCode TestParticleDataProcessingPipelineComputesSpecificKE ( void  )
static

Tests particle post-processing pipeline dispatch for specific-KE output.

Definition at line 120 of file test_postprocessor.c.

121{
122 SimCtx *simCtx = NULL;
123 UserCtx *user = NULL;
125 PetscReal (*vel_arr)[3] = NULL;
126 PetscScalar *ske_arr = NULL;
127
128 PetscFunctionBeginUser;
129 PetscCall(PetscMemzero(&pps, sizeof(pps)));
130 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
131 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
132 PetscCall(PetscStrncpy(pps.particle_pipeline, "ComputeSpecificKE:velocity>ske", sizeof(pps.particle_pipeline)));
133
134 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
135 vel_arr[0][0] = 1.0;
136 vel_arr[0][1] = 2.0;
137 vel_arr[0][2] = 2.0;
138 vel_arr[1][0] = 0.0;
139 vel_arr[1][1] = 3.0;
140 vel_arr[1][2] = 4.0;
141 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
142
143 PetscCall(ParticleDataProcessingPipeline(user, &pps));
144 PetscCall(DMSwarmGetField(user->post_swarm, "ske", NULL, NULL, (void *)&ske_arr));
145 PetscCall(PicurvAssertRealNear(4.5, PetscRealPart(ske_arr[0]), 1.0e-12,
146 "Particle pipeline should write first specific kinetic energy value"));
147 PetscCall(PicurvAssertRealNear(12.5, PetscRealPart(ske_arr[1]), 1.0e-12,
148 "Particle pipeline should write second specific kinetic energy value"));
149 PetscCall(DMSwarmRestoreField(user->post_swarm, "ske", NULL, NULL, (void *)&ske_arr));
150
151 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
152 PetscFunctionReturn(0);
153}
PetscErrorCode ParticleDataProcessingPipeline(UserCtx *user, PostProcessParams *pps)
Parses and executes the particle pipeline using a robust two-pass approach.
PetscErrorCode PicurvCreateSwarmPair(UserCtx *user, PetscInt nlocal, const char *post_field_name)
Creates matched solver and post-processing swarms for tests.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestGlobalStatisticsPipelineWritesMSDCSV()

static PetscErrorCode TestGlobalStatisticsPipelineWritesMSDCSV ( void  )
static

Tests global statistics pipeline generation of MSD CSV output.

Definition at line 158 of file test_postprocessor.c.

159{
160 SimCtx *simCtx = NULL;
161 UserCtx *user = NULL;
163 char tmpdir[PETSC_MAX_PATH_LEN];
164 char csv_prefix[PETSC_MAX_PATH_LEN];
165 char csv_path[PETSC_MAX_PATH_LEN];
166 PetscReal (*pos_arr)[3] = NULL;
167
168 PetscFunctionBeginUser;
169 PetscCall(PetscMemzero(&pps, sizeof(pps)));
170 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
171 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
172
173 simCtx->dt = 1.0;
174 simCtx->ren = 1.0;
175 simCtx->schmidt_number = 1.0;
176 simCtx->psrc_x = 0.0;
177 simCtx->psrc_y = 0.0;
178 simCtx->psrc_z = 0.0;
179
180 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
181 PetscCall(PetscSNPrintf(csv_prefix, sizeof(csv_prefix), "%s/stats", tmpdir));
182 PetscCall(PetscSNPrintf(csv_path, sizeof(csv_path), "%s_msd.csv", csv_prefix));
183
184 PetscCall(PetscStrncpy(pps.statistics_pipeline, "ComputeMSD", sizeof(pps.statistics_pipeline)));
185 PetscCall(PetscStrncpy(pps.statistics_output_prefix, csv_prefix, sizeof(pps.statistics_output_prefix)));
186
187 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
188 pos_arr[0][0] = 1.0;
189 pos_arr[0][1] = 0.0;
190 pos_arr[0][2] = 0.0;
191 pos_arr[1][0] = -1.0;
192 pos_arr[1][1] = 0.0;
193 pos_arr[1][2] = 0.0;
194 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
195
196 PetscCall(GlobalStatisticsPipeline(user, &pps, 1));
197 PetscCall(PicurvAssertFileExists(csv_path, "GlobalStatisticsPipeline should dispatch ComputeMSD and emit CSV output"));
198
199 PetscCall(PicurvRemoveTempDir(tmpdir));
200 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
201 PetscFunctionReturn(0);
202}
PetscErrorCode GlobalStatisticsPipeline(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Executes the global statistics pipeline, computing aggregate reductions over all particles.
PetscErrorCode PicurvMakeTempDir(char *path, size_t path_len)
Creates a unique temporary directory for one test case.
PetscErrorCode PicurvAssertFileExists(const char *path, const char *context)
Asserts that a filesystem path exists as a readable file.
PetscErrorCode PicurvRemoveTempDir(const char *path)
Recursively removes a temporary directory created by PicurvMakeTempDir.
char statistics_output_prefix[256]
basename for CSV output, e.g.
Definition variables.h:575
PetscReal schmidt_number
Definition variables.h:709
PetscReal psrc_x
Definition variables.h:706
PetscReal ren
Definition variables.h:691
PetscReal dt
Definition variables.h:658
char statistics_pipeline[1024]
e.g.
Definition variables.h:574
PetscReal psrc_z
Point source location for PARTICLE_INIT_POINT_SOURCE.
Definition variables.h:706
PetscReal psrc_y
Definition variables.h:706
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestGlobalStatisticsPipelineRewritesExistingMSDStep()

static PetscErrorCode TestGlobalStatisticsPipelineRewritesExistingMSDStep ( void  )
static

Tests that reprocessing the same MSD step rewrites the CSV row instead of duplicating it.

Definition at line 207 of file test_postprocessor.c.

208{
209 SimCtx *simCtx = NULL;
210 UserCtx *user = NULL;
212 char tmpdir[PETSC_MAX_PATH_LEN];
213 char csv_prefix[PETSC_MAX_PATH_LEN];
214 char csv_path[PETSC_MAX_PATH_LEN];
215 char csv_tmp_path[PETSC_MAX_PATH_LEN];
216 PetscReal (*pos_arr)[3] = NULL;
217 PetscInt row_count = 0;
218
219 PetscFunctionBeginUser;
220 PetscCall(PetscMemzero(&pps, sizeof(pps)));
221 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
222 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
223
224 simCtx->dt = 1.0;
225 simCtx->ren = 1.0;
226 simCtx->schmidt_number = 1.0;
227 simCtx->psrc_x = 0.0;
228 simCtx->psrc_y = 0.0;
229 simCtx->psrc_z = 0.0;
230
231 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
232 PetscCall(PetscSNPrintf(csv_prefix, sizeof(csv_prefix), "%s/stats", tmpdir));
233 PetscCall(PetscSNPrintf(csv_path, sizeof(csv_path), "%s_msd.csv", csv_prefix));
234 PetscCall(PetscSNPrintf(csv_tmp_path, sizeof(csv_tmp_path), "%s.tmp", csv_path));
235
236 PetscCall(PetscStrncpy(pps.statistics_pipeline, "ComputeMSD", sizeof(pps.statistics_pipeline)));
237 PetscCall(PetscStrncpy(pps.statistics_output_prefix, csv_prefix, sizeof(pps.statistics_output_prefix)));
238
239 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
240 pos_arr[0][0] = 1.0;
241 pos_arr[0][1] = 0.0;
242 pos_arr[0][2] = 0.0;
243 pos_arr[1][0] = -1.0;
244 pos_arr[1][1] = 0.0;
245 pos_arr[1][2] = 0.0;
246 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
247 PetscCall(GlobalStatisticsPipeline(user, &pps, 1));
248
249 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
250 pos_arr[0][0] = 2.0;
251 pos_arr[0][1] = 0.0;
252 pos_arr[0][2] = 0.0;
253 pos_arr[1][0] = -2.0;
254 pos_arr[1][1] = 0.0;
255 pos_arr[1][2] = 0.0;
256 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
257 PetscCall(GlobalStatisticsPipeline(user, &pps, 1));
258
259 PetscCall(CountLinesWithPrefix(csv_path, "1,", &row_count));
260 PetscCall(PicurvAssertIntEqual(1, row_count,
261 "Reprocessing the same MSD step should keep one CSV row for that step"));
262 PetscCall(PicurvAssertFileMissing(csv_tmp_path,
263 "Reprocessing the same MSD step should not leave a temporary CSV artifact behind"));
264
265 PetscCall(PicurvRemoveTempDir(tmpdir));
266 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
267 PetscFunctionReturn(0);
268}
static PetscErrorCode CountLinesWithPrefix(const char *path, const char *prefix, PetscInt *count_out)
Count file lines that begin with a given prefix.
static PetscErrorCode PicurvAssertFileMissing(const char *path, const char *context)
Assert that a file path does not exist.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestSetupSimulationEnvironmentCreatesStatisticsDirectoryWithoutClobbering()

static PetscErrorCode TestSetupSimulationEnvironmentCreatesStatisticsDirectoryWithoutClobbering ( void  )
static

Tests that postprocessor setup creates nested statistics directories without clobbering existing contents.

Definition at line 273 of file test_postprocessor.c.

274{
275 SimCtx *simCtx = NULL;
276 UserCtx *user = NULL;
277 PostProcessParams *pps = NULL;
278 char tmpdir[PETSC_MAX_PATH_LEN];
279 char post_path[PETSC_MAX_PATH_LEN];
280 char source_dir[PETSC_MAX_PATH_LEN];
281 char stats_dir[PETSC_MAX_PATH_LEN];
282 char stats_prefix[PETSC_MAX_PATH_LEN];
283 char sentinel_path[PETSC_MAX_PATH_LEN];
284 PetscBool exists = PETSC_FALSE;
285 FILE *file = NULL;
286
287 PetscFunctionBeginUser;
288 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
289 PetscCall(PetscNew(&simCtx->pps));
290 pps = simCtx->pps;
291 PetscCall(PetscMemzero(pps, sizeof(*pps)));
292 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
293
294 PetscCall(PetscSNPrintf(post_path, sizeof(post_path), "%s/post.run", tmpdir));
295 file = fopen(post_path, "w");
296 PetscCheck(file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to create post config %s.", post_path);
297 fputs("startTime = 0\nendTime = 0\ntimeStep = 1\n", file);
298 fclose(file);
299 file = NULL;
300
301 PetscCall(PetscSNPrintf(source_dir, sizeof(source_dir), "%s/output", tmpdir));
302 PetscCall(PicurvEnsureDir(source_dir));
303 PetscCall(PetscSNPrintf(stats_dir, sizeof(stats_dir), "%s/output/statistics", tmpdir));
304 PetscCall(PetscSNPrintf(stats_prefix, sizeof(stats_prefix), "%s/BrownianStats", stats_dir));
305 PetscCall(PetscSNPrintf(sentinel_path, sizeof(sentinel_path), "%s/sentinel.keep", stats_dir));
306
308 simCtx->generate_grid = PETSC_TRUE;
309 PetscCall(PetscStrncpy(simCtx->PostprocessingControlFile, post_path, sizeof(simCtx->PostprocessingControlFile)));
310 PetscCall(PetscStrncpy(pps->source_dir, source_dir, sizeof(pps->source_dir)));
311 PetscCall(PetscSNPrintf(pps->output_prefix, sizeof(pps->output_prefix), "%s/viz/Field", tmpdir));
312 PetscCall(PetscStrncpy(pps->statistics_pipeline, "ComputeMSD", sizeof(pps->statistics_pipeline)));
313 PetscCall(PetscStrncpy(pps->statistics_output_prefix, stats_prefix, sizeof(pps->statistics_output_prefix)));
314
315 PetscCall(SetupSimulationEnvironment(simCtx));
316 PetscCall(PetscTestDirectory(stats_dir, 'r', &exists));
317 PetscCall(PicurvAssertBool(exists,
318 "SetupSimulationEnvironment should create the statistics directory when it is missing"));
319
320 file = fopen(sentinel_path, "w");
321 PetscCheck(file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to create statistics sentinel %s.", sentinel_path);
322 fputs("keep\n", file);
323 fclose(file);
324 file = NULL;
325
326 PetscCall(SetupSimulationEnvironment(simCtx));
327 PetscCall(PicurvAssertFileExists(sentinel_path,
328 "SetupSimulationEnvironment should preserve existing statistics directory contents on repeat runs"));
329
330 PetscCall(PicurvRemoveTempDir(tmpdir));
331 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
332 PetscFunctionReturn(0);
333}
PetscErrorCode SetupSimulationEnvironment(SimCtx *simCtx)
Verifies and prepares the complete I/O environment for a simulation run.
Definition setup.c:813
PetscErrorCode PicurvEnsureDir(const char *path)
Ensures a directory exists for test output.
char output_prefix[256]
Definition variables.h:567
char source_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:555
PetscBool generate_grid
Definition variables.h:714
PostProcessParams * pps
Definition variables.h:798
@ EXEC_MODE_POSTPROCESSOR
Definition variables.h:617
ExecutionMode exec_mode
Definition variables.h:662
char PostprocessingControlFile[PETSC_MAX_PATH_LEN]
Definition variables.h:797
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestWriteEulerianFileWritesVTS()

static PetscErrorCode TestWriteEulerianFileWritesVTS ( void  )
static

Tests VTS emission for Eulerian post-processing output.

Definition at line 338 of file test_postprocessor.c.

339{
340 SimCtx *simCtx = NULL;
341 UserCtx *user = NULL;
343 char tmpdir[PETSC_MAX_PATH_LEN];
344 char vtk_path[PETSC_MAX_PATH_LEN];
345
346 PetscFunctionBeginUser;
347 PetscCall(PetscMemzero(&pps, sizeof(pps)));
348 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
349 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
350 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/field_00003.vts", tmpdir));
351
352 PetscCall(PetscStrncpy(pps.output_fields_instantaneous, "P_nodal", sizeof(pps.output_fields_instantaneous)));
353 PetscCall(PetscSNPrintf(pps.output_prefix, sizeof(pps.output_prefix), "%s/field", tmpdir));
354 PetscCall(VecSet(user->P_nodal, 7.0));
355
356 PetscCall(WriteEulerianFile(user, &pps, 3));
357 PetscCall(PicurvAssertFileExists(vtk_path, "WriteEulerianFile should emit a .vts file for requested output fields"));
358
359 PetscCall(PicurvRemoveTempDir(tmpdir));
360 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
361 PetscFunctionReturn(0);
362}
PetscErrorCode WriteEulerianFile(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Orchestrates the writing of a combined, multi-field VTK file for a single time step.
char output_fields_instantaneous[1024]
Definition variables.h:565
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestWriteParticleFileWritesVTP()

static PetscErrorCode TestWriteParticleFileWritesVTP ( void  )
static

Tests VTP emission for particle post-processing output.

Definition at line 367 of file test_postprocessor.c.

368{
369 SimCtx *simCtx = NULL;
370 UserCtx *user = NULL;
372 char tmpdir[PETSC_MAX_PATH_LEN];
373 char vtk_path[PETSC_MAX_PATH_LEN];
374 PetscReal (*pos_arr)[3] = NULL;
375 PetscReal (*vel_arr)[3] = NULL;
376
377 PetscFunctionBeginUser;
378 PetscCall(PetscMemzero(&pps, sizeof(pps)));
379 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
380 PetscCall(PicurvCreateSwarmPair(user, 3, "ske"));
381 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
382 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/particles_00004.vtp", tmpdir));
383
384 pps.outputParticles = PETSC_TRUE;
385 pps.particle_output_freq = 1;
386 PetscCall(PetscStrncpy(pps.particle_fields, "position,velocity", sizeof(pps.particle_fields)));
387 PetscCall(PetscSNPrintf(pps.particle_output_prefix, sizeof(pps.particle_output_prefix), "%s/particles", tmpdir));
388
389 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
390 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
391 for (PetscInt p = 0; p < 3; ++p) {
392 pos_arr[p][0] = (PetscReal)p;
393 pos_arr[p][1] = (PetscReal)(p + 1);
394 pos_arr[p][2] = (PetscReal)(p + 2);
395 vel_arr[p][0] = 2.0 * (PetscReal)p;
396 vel_arr[p][1] = 3.0 * (PetscReal)p;
397 vel_arr[p][2] = 4.0 * (PetscReal)p;
398 }
399 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
400 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
401
402 PetscCall(WriteParticleFile(user, &pps, 4));
403 PetscCall(PicurvAssertFileExists(vtk_path, "WriteParticleFile should emit a .vtp file for requested particle fields"));
404
405 PetscCall(PicurvRemoveTempDir(tmpdir));
406 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
407 PetscFunctionReturn(0);
408}
PetscErrorCode WriteParticleFile(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Writes particle data to a VTP file using the Prepare-Write-Cleanup pattern.
char particle_output_prefix[256]
Definition variables.h:570
PetscInt particle_output_freq
Definition variables.h:571
char particle_fields[1024]
Definition variables.h:569
PetscBool outputParticles
Definition variables.h:561
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestWriteEulerianFileRewritesSameStepCleanly()

static PetscErrorCode TestWriteEulerianFileRewritesSameStepCleanly ( void  )
static

Tests repeated same-step VTS emission leaves only the final artifact.

Definition at line 412 of file test_postprocessor.c.

413{
414 SimCtx *simCtx = NULL;
415 UserCtx *user = NULL;
417 char tmpdir[PETSC_MAX_PATH_LEN];
418 char vtk_path[PETSC_MAX_PATH_LEN];
419 char vtk_tmp_path[PETSC_MAX_PATH_LEN];
420
421 PetscFunctionBeginUser;
422 PetscCall(PetscMemzero(&pps, sizeof(pps)));
423 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
424 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
425 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/field_00003.vts", tmpdir));
426 PetscCall(PetscSNPrintf(vtk_tmp_path, sizeof(vtk_tmp_path), "%s.tmp", vtk_path));
427
428 PetscCall(PetscStrncpy(pps.output_fields_instantaneous, "P_nodal", sizeof(pps.output_fields_instantaneous)));
429 PetscCall(PetscSNPrintf(pps.output_prefix, sizeof(pps.output_prefix), "%s/field", tmpdir));
430 PetscCall(VecSet(user->P_nodal, 7.0));
431
432 PetscCall(WriteEulerianFile(user, &pps, 3));
433 PetscCall(VecSet(user->P_nodal, 9.0));
434 PetscCall(WriteEulerianFile(user, &pps, 3));
435
436 PetscCall(PicurvAssertFileExists(vtk_path,
437 "Rewriting the same Eulerian step should leave the final .vts file present"));
438 PetscCall(PicurvAssertFileMissing(vtk_tmp_path,
439 "Rewriting the same Eulerian step should not leave a temporary .vts file behind"));
440
441 PetscCall(PicurvRemoveTempDir(tmpdir));
442 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
443 PetscFunctionReturn(0);
444}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestWriteParticleFileRewritesSameStepCleanly()

static PetscErrorCode TestWriteParticleFileRewritesSameStepCleanly ( void  )
static

Tests repeated same-step VTP emission leaves only the final artifact.

Tests repeated same-step VTP emission leaves only the final artifact.

Definition at line 452 of file test_postprocessor.c.

453{
454 SimCtx *simCtx = NULL;
455 UserCtx *user = NULL;
457 char tmpdir[PETSC_MAX_PATH_LEN];
458 char vtk_path[PETSC_MAX_PATH_LEN];
459 char vtk_tmp_path[PETSC_MAX_PATH_LEN];
460 PetscReal (*pos_arr)[3] = NULL;
461 PetscReal (*vel_arr)[3] = NULL;
462
463 PetscFunctionBeginUser;
464 PetscCall(PetscMemzero(&pps, sizeof(pps)));
465 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
466 PetscCall(PicurvCreateSwarmPair(user, 3, "ske"));
467 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
468 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/particles_00004.vtp", tmpdir));
469 PetscCall(PetscSNPrintf(vtk_tmp_path, sizeof(vtk_tmp_path), "%s.tmp", vtk_path));
470
471 pps.outputParticles = PETSC_TRUE;
472 pps.particle_output_freq = 1;
473 PetscCall(PetscStrncpy(pps.particle_fields, "position,velocity", sizeof(pps.particle_fields)));
474 PetscCall(PetscSNPrintf(pps.particle_output_prefix, sizeof(pps.particle_output_prefix), "%s/particles", tmpdir));
475
476 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
477 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
478 for (PetscInt p = 0; p < 3; ++p) {
479 pos_arr[p][0] = (PetscReal)p;
480 pos_arr[p][1] = (PetscReal)(p + 1);
481 pos_arr[p][2] = (PetscReal)(p + 2);
482 vel_arr[p][0] = 2.0 * (PetscReal)p;
483 vel_arr[p][1] = 3.0 * (PetscReal)p;
484 vel_arr[p][2] = 4.0 * (PetscReal)p;
485 }
486 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
487 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
488
489 PetscCall(WriteParticleFile(user, &pps, 4));
490 PetscCall(WriteParticleFile(user, &pps, 4));
491
492 PetscCall(PicurvAssertFileExists(vtk_path,
493 "Rewriting the same particle step should leave the final .vtp file present"));
494 PetscCall(PicurvAssertFileMissing(vtk_tmp_path,
495 "Rewriting the same particle step should not leave a temporary .vtp file behind"));
496
497 PetscCall(PicurvRemoveTempDir(tmpdir));
498 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
499 PetscFunctionReturn(0);
500}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ main()

int main ( int  argc,
char **  argv 
)

Runs the unit-postprocessor PETSc test binary.

Definition at line 506 of file test_postprocessor.c.

507{
508 PetscErrorCode ierr;
509 const PicurvTestCase cases[] = {
510 {"setup-postprocess-swarm-registers-pipeline-fields", TestSetupPostProcessSwarmRegistersPipelineFields},
511 {"eulerian-data-processing-pipeline-runs-configured-kernels", TestEulerianDataProcessingPipelineRunsConfiguredKernels},
512 {"particle-data-processing-pipeline-computes-specific-ke", TestParticleDataProcessingPipelineComputesSpecificKE},
513 {"global-statistics-pipeline-writes-msd-csv", TestGlobalStatisticsPipelineWritesMSDCSV},
514 {"global-statistics-pipeline-rewrites-existing-msd-step", TestGlobalStatisticsPipelineRewritesExistingMSDStep},
515 {"setup-simulation-environment-creates-statistics-directory-without-clobbering", TestSetupSimulationEnvironmentCreatesStatisticsDirectoryWithoutClobbering},
516 {"write-eulerian-file-writes-vts", TestWriteEulerianFileWritesVTS},
517 {"write-particle-file-writes-vtp", TestWriteParticleFileWritesVTP},
518 {"write-eulerian-file-rewrites-same-step-cleanly", TestWriteEulerianFileRewritesSameStepCleanly},
519 {"write-particle-file-rewrites-same-step-cleanly", TestWriteParticleFileRewritesSameStepCleanly},
520 };
521
522 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv postprocessor tests");
523 if (ierr) {
524 return (int)ierr;
525 }
526
527 ierr = PicurvRunTests("unit-postprocessor", cases, sizeof(cases) / sizeof(cases[0]));
528 if (ierr) {
529 PetscFinalize();
530 return (int)ierr;
531 }
532
533 ierr = PetscFinalize();
534 return (int)ierr;
535}
static PetscErrorCode TestWriteParticleFileWritesVTP(void)
Tests VTP emission for particle post-processing output.
static PetscErrorCode TestSetupSimulationEnvironmentCreatesStatisticsDirectoryWithoutClobbering(void)
Tests that postprocessor setup creates nested statistics directories without clobbering existing cont...
static PetscErrorCode TestParticleDataProcessingPipelineComputesSpecificKE(void)
Tests particle post-processing pipeline dispatch for specific-KE output.
static PetscErrorCode TestGlobalStatisticsPipelineRewritesExistingMSDStep(void)
Tests that reprocessing the same MSD step rewrites the CSV row instead of duplicating it.
static PetscErrorCode TestWriteParticleFileRewritesSameStepCleanly(void)
Tests repeated same-step VTP emission leaves only the final artifact.
static PetscErrorCode TestGlobalStatisticsPipelineWritesMSDCSV(void)
Tests global statistics pipeline generation of MSD CSV output.
static PetscErrorCode TestWriteEulerianFileWritesVTS(void)
Tests VTS emission for Eulerian post-processing output.
static PetscErrorCode TestWriteEulerianFileRewritesSameStepCleanly(void)
Tests repeated same-step VTS emission leaves only the final artifact.
static PetscErrorCode TestEulerianDataProcessingPipelineRunsConfiguredKernels(void)
Tests Eulerian post-processing pipeline dispatch for configured kernels.
static PetscErrorCode TestSetupPostProcessSwarmRegistersPipelineFields(void)
Tests post-processing swarm setup and pipeline-field registration.
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Runs a named C test suite and prints pass/fail progress markers.
Named test case descriptor consumed by PicurvRunTests.
Here is the call graph for this function: