PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_postprocessing.c
Go to the documentation of this file.
1/**
2 * @file test_postprocessing.c
3 * @brief C test module for PICurv.
4 */
5
6#include "test_support.h"
7
9/**
10 * @brief Test-local routine.
11 */
12
13static PetscErrorCode TestComputeSpecificKE(void)
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}
42/**
43 * @brief Test-local routine.
44 */
45
46static PetscErrorCode TestComputeDisplacement(void)
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}
75/**
76 * @brief Test-local routine.
77 */
78
79static PetscErrorCode TestComputeNodalAverageScalar(void)
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[3][3][3]), 1.0e-12, "ComputeNodalAverage untouched boundary node"));
95 PetscCall(DMDAVecRestoreArrayRead(user->da, user->P_nodal, (void *)&p_nodal_arr));
96
97 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
98 PetscFunctionReturn(0);
99}
100/**
101 * @brief Test-local routine.
102 */
103
104static PetscErrorCode TestNormalizeRelativeField(void)
105{
106 SimCtx *simCtx = NULL;
107 UserCtx *user = NULL;
108 PetscInt ref_idx = 0;
109 PetscScalar ref_value = 0.0;
110 PetscReal vmin = 0.0, vmax = 0.0;
111
112 PetscFunctionBeginUser;
113 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
114 PetscCall(PetscCalloc1(1, &simCtx->pps));
115 simCtx->pps->reference[0] = 1;
116 simCtx->pps->reference[1] = 1;
117 simCtx->pps->reference[2] = 1;
118
119 PetscCall(VecSet(user->P, 10.0));
120 ref_idx = simCtx->pps->reference[2] * (user->IM * user->JM) +
121 simCtx->pps->reference[1] * user->IM +
122 simCtx->pps->reference[0];
123 PetscCall(VecSetValue(user->P, ref_idx, 4.0, INSERT_VALUES));
124 PetscCall(VecAssemblyBegin(user->P));
125 PetscCall(VecAssemblyEnd(user->P));
126
127 PetscCall(NormalizeRelativeField(user, "P"));
128
129 PetscCall(VecGetValues(user->P, 1, &ref_idx, &ref_value));
130 PetscCall(PicurvAssertRealNear(0.0, PetscRealPart(ref_value), 1.0e-12, "NormalizeRelativeField reference value"));
131 PetscCall(VecMin(user->P, NULL, &vmin));
132 PetscCall(VecMax(user->P, NULL, &vmax));
133 PetscCall(PicurvAssertRealNear(0.0, vmin, 1.0e-12, "NormalizeRelativeField minimum"));
134 PetscCall(PicurvAssertRealNear(6.0, vmax, 1.0e-12, "NormalizeRelativeField maximum"));
135
136 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
137 PetscFunctionReturn(0);
138}
139/**
140 * @brief Test-local routine.
141 */
142
143static PetscErrorCode TestDimensionalizePressureField(void)
144{
145 SimCtx *simCtx = NULL;
146 UserCtx *user = NULL;
147
148 PetscFunctionBeginUser;
149 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
150 simCtx->scaling.P_ref = 3.0;
151
152 PetscCall(VecSet(user->P, 2.0));
153 PetscCall(DimensionalizeField(user, "P"));
154 PetscCall(PicurvAssertVecConstant(user->P, 6.0, 1.0e-12, "DimensionalizeField should scale pressure by P_ref"));
155
156 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
157 PetscFunctionReturn(0);
158}
159/**
160 * @brief Test-local routine.
161 */
162
163static PetscErrorCode TestComputeQCriterionZeroFlow(void)
164{
165 SimCtx *simCtx = NULL;
166 UserCtx *user = NULL;
167
168 PetscFunctionBeginUser;
169 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
170
171 PetscCall(VecDuplicate(user->P, &user->Aj));
172 PetscCall(VecDuplicate(user->lP, &user->lAj));
173 PetscCall(VecSet(user->Aj, 1.0));
174 PetscCall(DMGlobalToLocalBegin(user->da, user->Aj, INSERT_VALUES, user->lAj));
175 PetscCall(DMGlobalToLocalEnd(user->da, user->Aj, INSERT_VALUES, user->lAj));
176
177 PetscCall(VecSet(user->Ucat, 0.0));
178 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
179 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
180 PetscCall(VecSet(user->Nvert, 0.0));
181 PetscCall(DMGlobalToLocalBegin(user->da, user->Nvert, INSERT_VALUES, user->lNvert));
182 PetscCall(DMGlobalToLocalEnd(user->da, user->Nvert, INSERT_VALUES, user->lNvert));
183
184 PetscCall(ComputeQCriterion(user));
185 PetscCall(PicurvAssertVecConstant(user->Qcrit, 0.0, 1.0e-12, "ComputeQCriterion should be zero for uniform flow"));
186
187 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
188 PetscFunctionReturn(0);
189}
190/**
191 * @brief Entry point for this unit-test binary.
192 */
193
194int main(int argc, char **argv)
195{
196 PetscErrorCode ierr;
197 const PicurvTestCase cases[] = {
198 {"compute-specific-ke", TestComputeSpecificKE},
199 {"compute-displacement", TestComputeDisplacement},
200 {"compute-nodal-average-scalar", TestComputeNodalAverageScalar},
201 {"normalize-relative-field", TestNormalizeRelativeField},
202 {"dimensionalize-pressure-field", TestDimensionalizePressureField},
203 {"compute-qcriterion-zero-flow", TestComputeQCriterionZeroFlow},
204 };
205
206 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv post-processing tests");
207 if (ierr) {
208 return (int)ierr;
209 }
210
211 ierr = PicurvRunTests("unit-post", cases, sizeof(cases) / sizeof(cases[0]));
212 if (ierr) {
213 PetscFinalize();
214 return (int)ierr;
215 }
216
217 ierr = PetscFinalize();
218 return (int)ierr;
219}
PetscErrorCode ComputeQCriterion(UserCtx *user)
Computes the Q-criterion diagnostic from the local velocity-gradient tensor.
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 ComputeDisplacement(UserCtx *user, const char *disp_field)
Computes the displacement magnitude |r_i - r_0| for each particle (per-particle VTK kernel).
PetscErrorCode NormalizeRelativeField(UserCtx *user, const char *relative_field_name)
Normalizes a relative scalar field using the configured reference pressure scale.
PetscErrorCode DimensionalizeField(UserCtx *user, const char *field_name)
Scales a specified field from non-dimensional to dimensional units in-place.
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.
int main(int argc, char **argv)
Entry point for this unit-test binary.
static PetscErrorCode TestComputeQCriterionZeroFlow(void)
Test-local routine.
static PetscErrorCode TestComputeSpecificKE(void)
Test-local routine.
static PetscErrorCode TestComputeDisplacement(void)
Test-local routine.
static PetscErrorCode TestDimensionalizePressureField(void)
Test-local routine.
static PetscErrorCode TestComputeNodalAverageScalar(void)
Test-local routine.
static PetscErrorCode TestNormalizeRelativeField(void)
Test-local 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 PicurvAssertVecConstant(Vec vec, PetscScalar expected, PetscReal tol, const char *context)
Shared test-support routine.
C test module for PICurv.
Named test case descriptor consumed by PicurvRunTests.
Vec P_nodal
Definition variables.h:814
Vec lNvert
Definition variables.h:764
PetscReal psrc_x
Definition variables.h:650
DM post_swarm
Definition variables.h:813
PetscInt reference[3]
Definition variables.h:530
Vec Qcrit
Definition variables.h:816
PetscReal psrc_z
Point source location for PARTICLE_INIT_POINT_SOURCE.
Definition variables.h:650
Vec Ucat
Definition variables.h:764
ScalingCtx scaling
Definition variables.h:651
PetscInt JM
Definition variables.h:747
PetscReal psrc_y
Definition variables.h:650
Vec lAj
Definition variables.h:785
Vec lUcat
Definition variables.h:764
PostProcessParams * pps
Definition variables.h:725
PetscInt IM
Definition variables.h:747
Vec Nvert
Definition variables.h:764
PetscReal P_ref
Definition variables.h:576
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