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

C test module for PICurv. More...

#include "test_support.h"
#include "vtk_io.h"
#include <stdio.h>
#include <string.h>
Include dependency graph for test_vtk_io.c:

Go to the source code of this file.

Functions

static PetscErrorCode TestPrepareOutputCoordinatesSubsamplesInteriorGrid (void)
 Test-local routine.
 
static PetscErrorCode TestPrepareOutputEulerianFieldDataSubsamplesScalar (void)
 Test-local routine.
 
static PetscErrorCode TestPrepareOutputParticleDataSubsampling (void)
 Test-local routine.
 
static PetscErrorCode TestCreateVTKFileFromMetadataWritesStructuredGrid (void)
 Test-local routine.
 
int main (int argc, char **argv)
 Entry point for this unit-test binary.
 

Detailed Description

C test module for PICurv.

Definition in file test_vtk_io.c.

Function Documentation

◆ TestPrepareOutputCoordinatesSubsamplesInteriorGrid()

static PetscErrorCode TestPrepareOutputCoordinatesSubsamplesInteriorGrid ( void  )
static

Test-local routine.

Definition at line 16 of file test_vtk_io.c.

17{
18 SimCtx *simCtx = NULL;
19 UserCtx *user = NULL;
20 PetscScalar *coords = NULL;
21 PetscInt nx = 0, ny = 0, nz = 0, npoints = 0;
22
23 PetscFunctionBeginUser;
24 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
25
26 PetscCall(PrepareOutputCoordinates(user, &coords, &nx, &ny, &nz, &npoints));
27 PetscCall(PicurvAssertIntEqual(3, nx, "PrepareOutputCoordinates should output nx=IM-1"));
28 PetscCall(PicurvAssertIntEqual(3, ny, "PrepareOutputCoordinates should output ny=JM-1"));
29 PetscCall(PicurvAssertIntEqual(3, nz, "PrepareOutputCoordinates should output nz=KM-1"));
30 PetscCall(PicurvAssertIntEqual(27, npoints, "PrepareOutputCoordinates should output (IM-1)*(JM-1)*(KM-1) points"));
31
32 if (simCtx->rank == 0) {
33 PetscCall(PicurvAssertBool((PetscBool)(coords != NULL), "PrepareOutputCoordinates should allocate coordinates on rank 0"));
34 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(coords[0]), 1.0e-12, "First coordinate x should be 0"));
35 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(coords[1]), 1.0e-12, "First coordinate y should be 0"));
36 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(coords[2]), 1.0e-12, "First coordinate z should be 0"));
37
38 PetscCall(PicurvAssertRealNear(2.0, PetscRealPart(coords[3 * (npoints - 1) + 0]), 1.0e-12, "Last subsampled coordinate x should be IM-2"));
39 PetscCall(PicurvAssertRealNear(2.0, PetscRealPart(coords[3 * (npoints - 1) + 1]), 1.0e-12, "Last subsampled coordinate y should be JM-2"));
40 PetscCall(PicurvAssertRealNear(2.0, PetscRealPart(coords[3 * (npoints - 1) + 2]), 1.0e-12, "Last subsampled coordinate z should be KM-2"));
41
42 PetscCall(PetscFree(coords));
43 }
44
45 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
46 PetscFunctionReturn(0);
47}
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 PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Shared test-support routine.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Shared test-support routine.
PetscMPIInt rank
Definition variables.h:594
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
PetscErrorCode PrepareOutputCoordinates(UserCtx *user, PetscScalar **out_coords, PetscInt *out_nx, PetscInt *out_ny, PetscInt *out_nz, PetscInt *out_npoints)
Creates a C array of coordinates corresponding to a subsampled (legacy-style) grid.
Definition vtk_io.c:215
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestPrepareOutputEulerianFieldDataSubsamplesScalar()

static PetscErrorCode TestPrepareOutputEulerianFieldDataSubsamplesScalar ( void  )
static

Test-local routine.

Definition at line 52 of file test_vtk_io.c.

53{
54 SimCtx *simCtx = NULL;
55 UserCtx *user = NULL;
56 PetscScalar *field_out = NULL;
57 PetscInt start = 0, end = 0;
58
59 PetscFunctionBeginUser;
60 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
61
62 PetscCall(VecGetOwnershipRange(user->P_nodal, &start, &end));
63 for (PetscInt idx = start; idx < end; ++idx) {
64 PetscCall(VecSetValue(user->P_nodal, idx, (PetscScalar)idx, INSERT_VALUES));
65 }
66 PetscCall(VecAssemblyBegin(user->P_nodal));
67 PetscCall(VecAssemblyEnd(user->P_nodal));
68
69 PetscCall(PrepareOutputEulerianFieldData(user, user->P_nodal, 1, &field_out));
70 if (simCtx->rank == 0) {
71 PetscCall(PicurvAssertBool((PetscBool)(field_out != NULL),
72 "PrepareOutputEulerianFieldData should allocate subsampled data on rank 0"));
73 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(field_out[0]), 1.0e-12, "Subsampled first value should map to source index 0"));
74 PetscCall(PicurvAssertRealNear(42.0, PetscRealPart(field_out[26]), 1.0e-12,
75 "Subsampled last value should map to source index (2,2,2)=42 in 4x4x4"));
76 PetscCall(PetscFree(field_out));
77 }
78
79 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
80 PetscFunctionReturn(0);
81}
Vec P_nodal
Definition variables.h:814
PetscErrorCode PrepareOutputEulerianFieldData(UserCtx *user, Vec field_vec, PetscInt num_components, PetscScalar **out_data)
Creates a C array of field data corresponding to a subsampled (legacy-style) grid.
Definition vtk_io.c:283
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestPrepareOutputParticleDataSubsampling()

static PetscErrorCode TestPrepareOutputParticleDataSubsampling ( void  )
static

Test-local routine.

Definition at line 86 of file test_vtk_io.c.

87{
88 SimCtx *simCtx = NULL;
89 UserCtx *user = NULL;
91 VTKMetaData meta;
92 PetscInt n_total = 0;
93 PetscReal (*position)[3] = NULL;
94 PetscReal (*velocity)[3] = NULL;
95
96 PetscFunctionBeginUser;
97 PetscCall(PetscMemzero(&pps, sizeof(pps)));
98 PetscCall(PetscMemzero(&meta, sizeof(meta)));
99 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
100 PetscCall(PicurvCreateSwarmPair(user, 5, "ske"));
101
102 pps.particle_output_freq = 2;
103 PetscCall(PetscStrncpy(pps.particle_fields, "position,velocity", sizeof(pps.particle_fields)));
104
105 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&position));
106 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&velocity));
107 for (PetscInt p = 0; p < 5; ++p) {
108 position[p][0] = (PetscReal)p;
109 position[p][1] = 0.0;
110 position[p][2] = 0.0;
111 velocity[p][0] = 10.0 + (PetscReal)p;
112 velocity[p][1] = 20.0 + (PetscReal)p;
113 velocity[p][2] = 30.0 + (PetscReal)p;
114 }
115 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&position));
116 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&velocity));
117
118 PetscCall(PrepareOutputParticleData(user, &pps, &meta, &n_total));
119 if (simCtx->rank == 0) {
120 PetscCall(PicurvAssertIntEqual(5, n_total, "PrepareOutputParticleData should report total particles before subsampling"));
121 PetscCall(PicurvAssertIntEqual(3, meta.npoints, "Stride-2 subsampling over 5 particles should keep 3 output particles"));
122 PetscCall(PicurvAssertIntEqual(VTK_POLYDATA, meta.fileType, "Particle output metadata should be PolyData"));
123 PetscCall(PicurvAssertIntEqual(1, meta.num_point_data_fields, "Velocity field should be exported as one point-data field"));
124 PetscCall(PicurvAssertBool((PetscBool)(strcmp(meta.point_data_fields[0].name, "velocity") == 0),
125 "Particle field name should propagate to VTK metadata"));
126
127 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(meta.coords[0]), 1.0e-12, "First exported particle should be index 0"));
128 PetscCall(PicurvAssertRealNear(2.0, PetscRealPart(meta.coords[3]), 1.0e-12, "Second exported particle should be index 2"));
129 PetscCall(PicurvAssertRealNear(4.0, PetscRealPart(meta.coords[6]), 1.0e-12, "Third exported particle should be index 4"));
130
131 PetscCall(PetscFree(meta.coords));
132 for (PetscInt f = 0; f < meta.num_point_data_fields; ++f) {
133 PetscCall(PetscFree(meta.point_data_fields[f].data));
134 }
135 PetscCall(PetscFree(meta.connectivity));
136 PetscCall(PetscFree(meta.offsets));
137 }
138
139 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
140 PetscFunctionReturn(0);
141}
PetscErrorCode PicurvCreateSwarmPair(UserCtx *user, PetscInt nlocal, const char *post_field_name)
Shared test-support routine.
PetscInt npoints
Definition variables.h:552
PetscInt num_point_data_fields
Definition variables.h:555
PetscInt * connectivity
Definition variables.h:556
PetscInt * offsets
Definition variables.h:557
VTKFileType fileType
Definition variables.h:550
PetscScalar * data
Definition variables.h:540
char name[64]
Definition variables.h:538
PetscScalar * coords
Definition variables.h:553
PetscInt particle_output_freq
Definition variables.h:519
char particle_fields[1024]
Definition variables.h:517
VTKFieldInfo point_data_fields[20]
Definition variables.h:554
@ VTK_POLYDATA
Definition variables.h:546
Holds all configuration parameters for a post-processing run.
Definition variables.h:501
PetscErrorCode PrepareOutputParticleData(UserCtx *user, PostProcessParams *pps, VTKMetaData *meta, PetscInt *p_n_total)
Gathers, subsamples, and prepares all particle data for VTK output.
Definition vtk_io.c:392
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestCreateVTKFileFromMetadataWritesStructuredGrid()

static PetscErrorCode TestCreateVTKFileFromMetadataWritesStructuredGrid ( void  )
static

Test-local routine.

Definition at line 146 of file test_vtk_io.c.

147{
148 char tmpdir[PETSC_MAX_PATH_LEN];
149 char vtk_path[PETSC_MAX_PATH_LEN];
150 VTKMetaData meta;
151 FILE *file = NULL;
152 char probe[512];
153 size_t nread = 0;
154
155 PetscFunctionBeginUser;
156 PetscCall(PetscMemzero(&meta, sizeof(meta)));
157 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
158 PetscCall(PetscSNPrintf(vtk_path, sizeof(vtk_path), "%s/field_00001.vts", tmpdir));
159
161 meta.mx = 2;
162 meta.my = 2;
163 meta.mz = 2;
164 meta.npoints = 8;
165 meta.num_point_data_fields = 1;
166
167 PetscCall(PetscMalloc1(3 * meta.npoints, &meta.coords));
168 PetscCall(PetscMalloc1(meta.npoints, &meta.point_data_fields[0].data));
169 PetscCall(PetscStrncpy(meta.point_data_fields[0].name, "P_nodal", sizeof(meta.point_data_fields[0].name)));
171
172 for (PetscInt p = 0; p < meta.npoints; ++p) {
173 meta.coords[3 * p + 0] = (PetscScalar)(p % 2);
174 meta.coords[3 * p + 1] = (PetscScalar)((p / 2) % 2);
175 meta.coords[3 * p + 2] = (PetscScalar)(p / 4);
176 meta.point_data_fields[0].data[p] = (PetscScalar)p;
177 }
178
179 PetscCall(CreateVTKFileFromMetadata(vtk_path, &meta, PETSC_COMM_WORLD));
180 PetscCall(PicurvAssertFileExists(vtk_path, "CreateVTKFileFromMetadata should emit a .vts file"));
181
182 file = fopen(vtk_path, "rb");
183 PetscCheck(file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open generated VTK file '%s'.", vtk_path);
184 nread = fread(probe, 1, sizeof(probe) - 1, file);
185 fclose(file);
186 probe[nread] = '\0';
187 PetscCall(PicurvAssertBool((PetscBool)(strstr(probe, "StructuredGrid") != NULL),
188 "Generated VTK header should declare StructuredGrid type"));
189
190 PetscCall(PetscFree(meta.coords));
191 PetscCall(PetscFree(meta.point_data_fields[0].data));
192 PetscFunctionReturn(0);
193}
PetscInt CreateVTKFileFromMetadata(const char *filename, const VTKMetaData *meta, MPI_Comm comm)
Implementation of CreateVTKFileFromMetadata().
Definition vtk_io.c:133
PetscErrorCode PicurvMakeTempDir(char *path, size_t path_len)
Shared test-support routine.
PetscErrorCode PicurvAssertFileExists(const char *path, const char *context)
Shared test-support routine.
PetscInt num_components
Definition variables.h:539
PetscInt mz
Definition variables.h:551
PetscInt my
Definition variables.h:551
PetscInt mx
Definition variables.h:551
@ VTK_STRUCTURED
Definition variables.h:545
Here is the call graph for this function:
Here is the caller graph for this function:

◆ main()

int main ( int  argc,
char **  argv 
)

Entry point for this unit-test binary.

Definition at line 198 of file test_vtk_io.c.

199{
200 PetscErrorCode ierr;
201 const PicurvTestCase cases[] = {
202 {"prepare-output-coordinates-subsamples-interior-grid", TestPrepareOutputCoordinatesSubsamplesInteriorGrid},
203 {"prepare-output-eulerian-field-data-subsamples-scalar", TestPrepareOutputEulerianFieldDataSubsamplesScalar},
204 {"prepare-output-particle-data-subsampling", TestPrepareOutputParticleDataSubsampling},
205 {"create-vtk-file-from-metadata-writes-structured-grid", TestCreateVTKFileFromMetadataWritesStructuredGrid},
206 };
207
208 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv VTK I/O tests");
209 if (ierr) {
210 return (int)ierr;
211 }
212
213 ierr = PicurvRunTests("unit-post-vtk", cases, sizeof(cases) / sizeof(cases[0]));
214 if (ierr) {
215 PetscFinalize();
216 return (int)ierr;
217 }
218
219 ierr = PetscFinalize();
220 return (int)ierr;
221}
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Shared test-support routine.
Named test case descriptor consumed by PicurvRunTests.
static PetscErrorCode TestPrepareOutputEulerianFieldDataSubsamplesScalar(void)
Test-local routine.
Definition test_vtk_io.c:52
static PetscErrorCode TestPrepareOutputParticleDataSubsampling(void)
Test-local routine.
Definition test_vtk_io.c:86
static PetscErrorCode TestCreateVTKFileFromMetadataWritesStructuredGrid(void)
Test-local routine.
static PetscErrorCode TestPrepareOutputCoordinatesSubsamplesInteriorGrid(void)
Test-local routine.
Definition test_vtk_io.c:16
Here is the call graph for this function: