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

C unit tests for post-processing kernel functions. More...

#include "test_support.h"
#include "postprocessing_kernels.h"
Include dependency graph for test_postprocessing.c:

Go to the source code of this file.

Functions

static PetscErrorCode TestComputeSpecificKE (void)
 Tests specific turbulent kinetic energy computation on particle data.
 
static PetscErrorCode TestComputeDisplacement (void)
 Tests particle displacement computation against reference positions.
 
static PetscErrorCode TestComputeNodalAverageScalar (void)
 Tests nodal averaging of scalar cell-centered data.
 
static PetscErrorCode TestNormalizeRelativeField (void)
 Tests normalization of one field relative to a reference field.
 
static PetscErrorCode TestDimensionalizePressureField (void)
 Tests dimensionalization of pressure data from nondimensional values.
 
static PetscErrorCode TestComputeQCriterionZeroFlow (void)
 Tests Q-criterion computation for a quiescent velocity field.
 
int main (int argc, char **argv)
 Runs the unit-post PETSc test binary.
 

Detailed Description

C unit tests for post-processing kernel functions.

Definition in file test_postprocessing.c.

Function Documentation

◆ TestComputeSpecificKE()

static PetscErrorCode TestComputeSpecificKE ( void  )
static

Tests specific turbulent kinetic energy computation on particle data.

Definition at line 13 of file test_postprocessing.c.

14{
15 SimCtx *simCtx = NULL;
16 UserCtx *user = NULL;
17 PetscReal (*vel_arr)[3] = NULL;
18 PetscScalar *ske_arr = NULL;
19
20 PetscFunctionBeginUser;
21 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
22 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
23
24 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
25 vel_arr[0][0] = 1.0;
26 vel_arr[0][1] = 2.0;
27 vel_arr[0][2] = 2.0;
28 vel_arr[1][0] = 0.0;
29 vel_arr[1][1] = 3.0;
30 vel_arr[1][2] = 4.0;
31 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
32
33 PetscCall(ComputeSpecificKE(user, "velocity", "ske"));
34 PetscCall(DMSwarmGetField(user->post_swarm, "ske", NULL, NULL, (void *)&ske_arr));
35 PetscCall(PicurvAssertRealNear(4.5, ske_arr[0], 1.0e-12, "ComputeSpecificKE first particle"));
36 PetscCall(PicurvAssertRealNear(12.5, ske_arr[1], 1.0e-12, "ComputeSpecificKE second particle"));
37 PetscCall(DMSwarmRestoreField(user->post_swarm, "ske", NULL, NULL, (void *)&ske_arr));
38
39 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
40 PetscFunctionReturn(0);
41}
PetscErrorCode ComputeSpecificKE(UserCtx *user, const char *velocity_field, const char *ske_field)
Computes the specific kinetic energy (KE per unit mass) for each particle.
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 PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Asserts that two real values agree within tolerance.
PetscErrorCode PicurvDestroyMinimalContexts(SimCtx **simCtx_ptr, UserCtx **user_ptr)
Destroys minimal SimCtx/UserCtx fixtures and all owned PETSc objects.
PetscErrorCode PicurvCreateSwarmPair(UserCtx *user, PetscInt nlocal, const char *post_field_name)
Creates matched solver and post-processing swarms for tests.
DM post_swarm
Definition variables.h:886
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:

◆ TestComputeDisplacement()

static PetscErrorCode TestComputeDisplacement ( void  )
static

Tests particle displacement computation against reference positions.

Definition at line 46 of file test_postprocessing.c.

47{
48 SimCtx *simCtx = NULL;
49 UserCtx *user = NULL;
50 PetscReal (*pos_arr)[3] = NULL;
51 PetscScalar *disp_arr = NULL;
52
53 PetscFunctionBeginUser;
54 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
55 PetscCall(PicurvCreateSwarmPair(user, 2, "disp"));
56
57 simCtx->psrc_x = 1.0;
58 simCtx->psrc_y = 2.0;
59 simCtx->psrc_z = 3.0;
60
61 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
62 pos_arr[0][0] = 1.0; pos_arr[0][1] = 2.0; pos_arr[0][2] = 3.0;
63 pos_arr[1][0] = 4.0; pos_arr[1][1] = 6.0; pos_arr[1][2] = 3.0;
64 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
65
66 PetscCall(ComputeDisplacement(user, "disp"));
67 PetscCall(DMSwarmGetField(user->post_swarm, "disp", NULL, NULL, (void *)&disp_arr));
68 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(disp_arr[0]), 1.0e-12, "ComputeDisplacement first particle"));
69 PetscCall(PicurvAssertRealNear(5.0, PetscRealPart(disp_arr[1]), 1.0e-12, "ComputeDisplacement second particle"));
70 PetscCall(DMSwarmRestoreField(user->post_swarm, "disp", NULL, NULL, (void *)&disp_arr));
71
72 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
73 PetscFunctionReturn(0);
74}
PetscErrorCode ComputeDisplacement(UserCtx *user, const char *disp_field)
Computes the displacement magnitude |r_i - r_0| for each particle (per-particle VTK kernel).
PetscReal psrc_x
Definition variables.h:706
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:

◆ TestComputeNodalAverageScalar()

static PetscErrorCode TestComputeNodalAverageScalar ( void  )
static

Tests nodal averaging of scalar cell-centered data.

Definition at line 79 of file test_postprocessing.c.

80{
81 SimCtx *simCtx = NULL;
82 UserCtx *user = NULL;
83 const PetscScalar ***p_nodal_arr = NULL;
84
85 PetscFunctionBeginUser;
86 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
87 PetscCall(VecSet(user->P, 7.0));
88 PetscCall(VecSet(user->P_nodal, -1.0));
89
90 PetscCall(ComputeNodalAverage(user, "P", "P_nodal"));
91
92 PetscCall(DMDAVecGetArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
93 PetscCall(PicurvAssertRealNear(7.0, PetscRealPart(p_nodal_arr[0][0][0]), 1.0e-12, "ComputeNodalAverage interior node"));
94 PetscCall(PicurvAssertRealNear(-1.0, PetscRealPart(p_nodal_arr[user->KM][user->JM][user->IM]), 1.0e-12,
95 "ComputeNodalAverage untouched non-physical boundary node"));
96 PetscCall(DMDAVecRestoreArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
97
98 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
99 PetscFunctionReturn(0);
100}
PetscErrorCode ComputeNodalAverage(UserCtx *user, const char *in_field_name, const char *out_field_name)
Interpolates a cell-centered field to nodal locations using local stencil averaging.
Vec P_nodal
Definition variables.h:887
PetscInt KM
Definition variables.h:820
PetscInt JM
Definition variables.h:820
PetscInt IM
Definition variables.h:820
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestNormalizeRelativeField()

static PetscErrorCode TestNormalizeRelativeField ( void  )
static

Tests normalization of one field relative to a reference field.

Definition at line 105 of file test_postprocessing.c.

106{
107 SimCtx *simCtx = NULL;
108 UserCtx *user = NULL;
109 PetscInt ref_idx = 0;
110 PetscScalar ref_value = 0.0;
111 PetscReal vmin = 0.0, vmax = 0.0;
112
113 PetscFunctionBeginUser;
114 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
115 PetscCall(PetscCalloc1(1, &simCtx->pps));
116 simCtx->pps->reference[0] = 1;
117 simCtx->pps->reference[1] = 1;
118 simCtx->pps->reference[2] = 1;
119
120 PetscCall(VecSet(user->P, 10.0));
121 ref_idx = simCtx->pps->reference[2] * (user->IM * user->JM) +
122 simCtx->pps->reference[1] * user->IM +
123 simCtx->pps->reference[0];
124 PetscCall(VecSetValue(user->P, ref_idx, 4.0, INSERT_VALUES));
125 PetscCall(VecAssemblyBegin(user->P));
126 PetscCall(VecAssemblyEnd(user->P));
127
128 PetscCall(NormalizeRelativeField(user, "P"));
129
130 PetscCall(VecGetValues(user->P, 1, &ref_idx, &ref_value));
131 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(ref_value), 1.0e-12, "NormalizeRelativeField reference value"));
132 PetscCall(VecMin(user->P, NULL, &vmin));
133 PetscCall(VecMax(user->P, NULL, &vmax));
134 PetscCall(PicurvAssertRealNear(0.0, vmin, 1.0e-12, "NormalizeRelativeField minimum"));
135 PetscCall(PicurvAssertRealNear(6.0, vmax, 1.0e-12, "NormalizeRelativeField maximum"));
136
137 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
138 PetscFunctionReturn(0);
139}
PetscErrorCode NormalizeRelativeField(UserCtx *user, const char *relative_field_name)
Normalizes a relative scalar field using the configured reference pressure scale.
PetscInt reference[3]
Definition variables.h:582
PostProcessParams * pps
Definition variables.h:798
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestDimensionalizePressureField()

static PetscErrorCode TestDimensionalizePressureField ( void  )
static

Tests dimensionalization of pressure data from nondimensional values.

Definition at line 144 of file test_postprocessing.c.

145{
146 SimCtx *simCtx = NULL;
147 UserCtx *user = NULL;
148
149 PetscFunctionBeginUser;
150 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
151 simCtx->scaling.P_ref = 3.0;
152
153 PetscCall(VecSet(user->P, 2.0));
154 PetscCall(DimensionalizeField(user, "P"));
155 PetscCall(PicurvAssertVecConstant(user->P, 6.0, 1.0e-12, "DimensionalizeField should scale pressure by P_ref"));
156
157 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
158 PetscFunctionReturn(0);
159}
PetscErrorCode DimensionalizeField(UserCtx *user, const char *field_name)
Scales a specified field from non-dimensional to dimensional units in-place.
PetscErrorCode PicurvAssertVecConstant(Vec vec, PetscScalar expected, PetscReal tol, const char *context)
Asserts that a PETSc vector is spatially constant within tolerance.
ScalingCtx scaling
Definition variables.h:707
PetscReal P_ref
Definition variables.h:628
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestComputeQCriterionZeroFlow()

static PetscErrorCode TestComputeQCriterionZeroFlow ( void  )
static

Tests Q-criterion computation for a quiescent velocity field.

Definition at line 164 of file test_postprocessing.c.

165{
166 SimCtx *simCtx = NULL;
167 UserCtx *user = NULL;
168
169 PetscFunctionBeginUser;
170 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
171
172 PetscCall(VecSet(user->Aj, 1.0));
173 PetscCall(DMGlobalToLocalBegin(user->da, user->Aj, INSERT_VALUES, user->lAj));
174 PetscCall(DMGlobalToLocalEnd(user->da, user->Aj, INSERT_VALUES, user->lAj));
175
176 PetscCall(VecSet(user->Ucat, 0.0));
177 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
178 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
179 PetscCall(VecSet(user->Nvert, 0.0));
180 PetscCall(DMGlobalToLocalBegin(user->da, user->Nvert, INSERT_VALUES, user->lNvert));
181 PetscCall(DMGlobalToLocalEnd(user->da, user->Nvert, INSERT_VALUES, user->lNvert));
182
183 PetscCall(ComputeQCriterion(user));
184 PetscCall(PicurvAssertVecConstant(user->Qcrit, 0.0, 1.0e-12, "ComputeQCriterion should be zero for uniform flow"));
185
186 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
187 PetscFunctionReturn(0);
188}
PetscErrorCode ComputeQCriterion(UserCtx *user)
Computes the Q-criterion diagnostic from the local velocity-gradient tensor.
Vec lNvert
Definition variables.h:837
Vec Qcrit
Definition variables.h:889
Vec Ucat
Definition variables.h:837
Vec lAj
Definition variables.h:858
Vec lUcat
Definition variables.h:837
Vec Nvert
Definition variables.h:837
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-post PETSc test binary.

Definition at line 193 of file test_postprocessing.c.

194{
195 PetscErrorCode ierr;
196 const PicurvTestCase cases[] = {
197 {"compute-specific-ke", TestComputeSpecificKE},
198 {"compute-displacement", TestComputeDisplacement},
199 {"compute-nodal-average-scalar", TestComputeNodalAverageScalar},
200 {"normalize-relative-field", TestNormalizeRelativeField},
201 {"dimensionalize-pressure-field", TestDimensionalizePressureField},
202 {"compute-qcriterion-zero-flow", TestComputeQCriterionZeroFlow},
203 };
204
205 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv post-processing tests");
206 if (ierr) {
207 return (int)ierr;
208 }
209
210 ierr = PicurvRunTests("unit-post", cases, sizeof(cases) / sizeof(cases[0]));
211 if (ierr) {
212 PetscFinalize();
213 return (int)ierr;
214 }
215
216 ierr = PetscFinalize();
217 return (int)ierr;
218}
static PetscErrorCode TestComputeQCriterionZeroFlow(void)
Tests Q-criterion computation for a quiescent velocity field.
static PetscErrorCode TestComputeSpecificKE(void)
Tests specific turbulent kinetic energy computation on particle data.
static PetscErrorCode TestComputeDisplacement(void)
Tests particle displacement computation against reference positions.
static PetscErrorCode TestDimensionalizePressureField(void)
Tests dimensionalization of pressure data from nondimensional values.
static PetscErrorCode TestComputeNodalAverageScalar(void)
Tests nodal averaging of scalar cell-centered data.
static PetscErrorCode TestNormalizeRelativeField(void)
Tests normalization of one field relative to a reference field.
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: