PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_poisson_rhs.c
Go to the documentation of this file.
1/**
2 * @file test_poisson_rhs.c
3 * @brief C test module for PICurv.
4 */
5
6#include "test_support.h"
7
8#include "poisson.h"
9#include "rhs.h"
10/**
11 * @brief Test-local routine.
12 */
13
14static PetscErrorCode EnsurePoissonAndRhsVectors(UserCtx *user)
15{
16 PetscFunctionBeginUser;
17 if (!user->Phi) PetscCall(DMCreateGlobalVector(user->da, &user->Phi));
18 if (!user->lPhi) PetscCall(DMCreateLocalVector(user->da, &user->lPhi));
19 if (!user->Aj) PetscCall(DMCreateGlobalVector(user->da, &user->Aj));
20 if (!user->lAj) PetscCall(DMCreateLocalVector(user->da, &user->lAj));
21 if (!user->Diffusivity) PetscCall(DMCreateGlobalVector(user->da, &user->Diffusivity));
22 if (!user->lDiffusivity) PetscCall(DMCreateLocalVector(user->da, &user->lDiffusivity));
23 PetscFunctionReturn(0);
24}
25/**
26 * @brief Test-local routine.
27 */
28
29static PetscErrorCode TestUpdatePressureAddsPhi(void)
30{
31 SimCtx *simCtx = NULL;
32 UserCtx *user = NULL;
33
34 PetscFunctionBeginUser;
35 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
36 PetscCall(EnsurePoissonAndRhsVectors(user));
37 PetscCall(VecSet(user->P, 1.0));
38 PetscCall(VecSet(user->Phi, 0.25));
39
40 PetscCall(UpdatePressure(user));
41 PetscCall(PicurvAssertVecConstant(user->P, 1.25, 1.0e-12, "UpdatePressure should add Phi into P"));
42
43 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
44 PetscFunctionReturn(0);
45}
46/**
47 * @brief Test-local routine.
48 */
49
50static PetscErrorCode TestPoissonRHSZeroDivergence(void)
51{
52 SimCtx *simCtx = NULL;
53 UserCtx *user = NULL;
54 Vec B = NULL;
55
56 PetscFunctionBeginUser;
57 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
58 PetscCall(EnsurePoissonAndRhsVectors(user));
59
60 simCtx->dt = 0.5;
61 PetscCall(VecSet(user->Ucont, 0.0));
62 PetscCall(VecSet(user->Nvert, 0.0));
63 PetscCall(VecSet(user->Aj, 1.0));
64 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
65 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
66 PetscCall(DMGlobalToLocalBegin(user->da, user->Nvert, INSERT_VALUES, user->lNvert));
67 PetscCall(DMGlobalToLocalEnd(user->da, user->Nvert, INSERT_VALUES, user->lNvert));
68 PetscCall(DMGlobalToLocalBegin(user->da, user->Aj, INSERT_VALUES, user->lAj));
69 PetscCall(DMGlobalToLocalEnd(user->da, user->Aj, INSERT_VALUES, user->lAj));
70
71 PetscCall(VecDuplicate(user->P, &B));
72 PetscCall(PoissonRHS(user, B));
73 PetscCall(PicurvAssertVecConstant(B, 0.0, 1.0e-12, "zero velocity divergence should produce zero Poisson RHS"));
74 PetscCall(PicurvAssertRealNear(0.0, simCtx->summationRHS, 1.0e-12, "global Poisson RHS sum should be zero"));
75
76 PetscCall(VecDestroy(&B));
77 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
78 PetscFunctionReturn(0);
79}
80/**
81 * @brief Test-local routine.
82 */
83
84static PetscErrorCode TestComputeBodyForcesDispatcher(void)
85{
86 SimCtx *simCtx = NULL;
87 UserCtx *user = NULL;
88 Vec rct = NULL;
89 Cmpnts ***rct_arr = NULL;
90
91 PetscFunctionBeginUser;
92 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
93 PetscCall(VecDuplicate(user->Ucont, &rct));
94 PetscCall(VecZeroEntries(rct));
95
98 simCtx->bulkVelocityCorrection = 2.0;
99 simCtx->dt = 1.0;
100 simCtx->forceScalingFactor = 1.0;
101 simCtx->drivingForceMagnitude = 0.0;
102
103 PetscCall(ComputeBodyForces(user, rct));
104 PetscCall(DMDAVecGetArrayRead(user->fda, rct, &rct_arr));
105 PetscCall(PicurvAssertRealNear(1.5, rct_arr[1][1][1].x, 1.0e-12, "ComputeBodyForces should apply driven-channel source"));
106 PetscCall(PicurvAssertRealNear(0.0, rct_arr[1][1][1].y, 1.0e-12, "ComputeBodyForces should keep y unchanged"));
107 PetscCall(PicurvAssertRealNear(0.0, rct_arr[1][1][1].z, 1.0e-12, "ComputeBodyForces should keep z unchanged"));
108 PetscCall(DMDAVecRestoreArrayRead(user->fda, rct, &rct_arr));
109
110 PetscCall(VecDestroy(&rct));
111 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
112 PetscFunctionReturn(0);
113}
114/**
115 * @brief Test-local routine.
116 */
117
119{
120 SimCtx *simCtx = NULL;
121 UserCtx *user = NULL;
122 PetscReal ***diff_arr = NULL;
123
124 PetscFunctionBeginUser;
125 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
126 PetscCall(EnsurePoissonAndRhsVectors(user));
127
128 simCtx->ren = 2.0;
129 simCtx->schmidt_number = 4.0;
130 simCtx->les = PETSC_FALSE;
131 simCtx->rans = PETSC_FALSE;
132 PetscCall(VecZeroEntries(user->Diffusivity));
133
134 PetscCall(ComputeEulerianDiffusivity(user));
135 PetscCall(DMDAVecGetArrayRead(user->da, user->Diffusivity, &diff_arr));
136 PetscCall(PicurvAssertRealNear(0.125, diff_arr[1][1][1], 1.0e-12, "molecular diffusivity interior value"));
137 PetscCall(DMDAVecRestoreArrayRead(user->da, user->Diffusivity, &diff_arr));
138
139 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
140 PetscFunctionReturn(0);
141}
142/**
143 * @brief Entry point for this unit-test binary.
144 */
145
146int main(int argc, char **argv)
147{
148 PetscErrorCode ierr;
149 const PicurvTestCase cases[] = {
150 {"update-pressure-adds-phi", TestUpdatePressureAddsPhi},
151 {"poisson-rhs-zero-divergence", TestPoissonRHSZeroDivergence},
152 {"compute-body-forces-dispatcher", TestComputeBodyForcesDispatcher},
153 {"compute-eulerian-diffusivity-molecular-only", TestComputeEulerianDiffusivityMolecularOnly},
154 };
155
156 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv Poisson/RHS tests");
157 if (ierr) {
158 return (int)ierr;
159 }
160
161 ierr = PicurvRunTests("unit-poisson-rhs", cases, sizeof(cases) / sizeof(cases[0]));
162 if (ierr) {
163 PetscFinalize();
164 return (int)ierr;
165 }
166
167 ierr = PetscFinalize();
168 return (int)ierr;
169}
PetscErrorCode UpdatePressure(UserCtx *user)
Updates the pressure field P with the pressure correction Phi computed by the Poisson solver.
Definition poisson.c:846
PetscErrorCode PoissonRHS(UserCtx *user, Vec B)
Computes the Right-Hand-Side (RHS) of the Poisson equation, which is the divergence of the intermedia...
Definition poisson.c:2129
PetscErrorCode ComputeBodyForces(UserCtx *user, Vec Rct)
General dispatcher for applying all active body forces (momentum sources).
Definition rhs.c:1161
PetscErrorCode ComputeEulerianDiffusivity(UserCtx *user)
Computes the effective diffusivity scalar field (Gamma_eff) on the Eulerian grid.
Definition rhs.c:1956
static PetscErrorCode TestComputeEulerianDiffusivityMolecularOnly(void)
Test-local routine.
int main(int argc, char **argv)
Entry point for this unit-test binary.
static PetscErrorCode TestComputeBodyForcesDispatcher(void)
Test-local routine.
static PetscErrorCode TestUpdatePressureAddsPhi(void)
Test-local routine.
static PetscErrorCode TestPoissonRHSZeroDivergence(void)
Test-local routine.
static PetscErrorCode EnsurePoissonAndRhsVectors(UserCtx *user)
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 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.
PetscReal schmidt_number
Definition variables.h:653
BoundaryFaceConfig boundary_faces[6]
Definition variables.h:756
Vec lNvert
Definition variables.h:764
Vec Phi
Definition variables.h:764
PetscReal forceScalingFactor
Definition variables.h:667
PetscInt rans
Definition variables.h:676
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
Definition variables.h:245
PetscReal ren
Definition variables.h:638
BCHandlerType handler_type
Definition variables.h:296
PetscReal dt
Definition variables.h:606
PetscReal bulkVelocityCorrection
Definition variables.h:669
Vec Ucont
Definition variables.h:764
Vec lPhi
Definition variables.h:764
Vec lUcont
Definition variables.h:764
Vec Diffusivity
Definition variables.h:767
Vec lAj
Definition variables.h:785
PetscInt les
Definition variables.h:676
Vec Nvert
Definition variables.h:764
Vec lDiffusivity
Definition variables.h:767
PetscReal summationRHS
Definition variables.h:711
PetscReal drivingForceMagnitude
Definition variables.h:667
@ BC_FACE_NEG_X
Definition variables.h:204
A 3D point or vector with PetscScalar components.
Definition variables.h:100
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