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

C test module for PICurv. More...

#include "test_support.h"
#include "interpolation.h"
#include "walkingsearch.h"
Include dependency graph for test_geometry.c:

Go to the source code of this file.

Functions

static PetscErrorCode TestScalarInterpolationConstantField (void)
 Test-local routine.
 
static PetscErrorCode TestVectorInterpolationConstantField (void)
 Test-local routine.
 
static PetscErrorCode TestSignedDistanceAndClassification (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_geometry.c.

Function Documentation

◆ TestScalarInterpolationConstantField()

static PetscErrorCode TestScalarInterpolationConstantField ( void  )
static

Test-local routine.

Definition at line 14 of file test_geometry.c.

15{
16 SimCtx *simCtx = NULL;
17 UserCtx *user = NULL;
18 Vec input = NULL, output = NULL;
19 PetscReal ***in_arr = NULL, ***out_arr = NULL;
20
21 PetscFunctionBeginUser;
22 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
23 PetscCall(VecDuplicate(user->P, &input));
24 PetscCall(VecDuplicate(user->P, &output));
25 PetscCall(VecSet(input, 8.0));
26 PetscCall(VecSet(output, 0.0));
27
28 PetscCall(DMDAVecGetArrayRead(user->da, input, &in_arr));
29 PetscCall(DMDAVecGetArray(user->da, output, &out_arr));
30 PetscCall(InterpolateFieldFromCornerToCenter_Scalar(in_arr, out_arr, user));
31 PetscCall(PicurvAssertRealNear(8.0, out_arr[1][1][1], 1.0e-12, "scalar corner->center preserves constants"));
32 PetscCall(DMDAVecRestoreArrayRead(user->da, input, &in_arr));
33 PetscCall(DMDAVecRestoreArray(user->da, output, &out_arr));
34
35 PetscCall(VecDestroy(&input));
36 PetscCall(VecDestroy(&output));
37 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
38 PetscFunctionReturn(0);
39}
PetscErrorCode InterpolateFieldFromCornerToCenter_Scalar(PetscReal ***field_arr, PetscReal ***centfield_arr, UserCtx *user)
Safely interpolate a scalar field from corner nodes (from the coordinate DM) to cell centers (from th...
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.
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestVectorInterpolationConstantField()

static PetscErrorCode TestVectorInterpolationConstantField ( void  )
static

Test-local routine.

Definition at line 44 of file test_geometry.c.

45{
46 SimCtx *simCtx = NULL;
47 UserCtx *user = NULL;
48 Vec input = NULL, output = NULL;
49 Cmpnts ***in_arr = NULL, ***out_arr = NULL;
50
51 PetscFunctionBeginUser;
52 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
53 PetscCall(VecDuplicate(user->Ucat, &input));
54 PetscCall(VecDuplicate(user->Ucat, &output));
55 PetscCall(VecSet(input, 3.5));
56 PetscCall(VecSet(output, 0.0));
57
58 PetscCall(DMDAVecGetArrayRead(user->fda, input, &in_arr));
59 PetscCall(DMDAVecGetArray(user->fda, output, &out_arr));
60 PetscCall(InterpolateFieldFromCornerToCenter_Vector(in_arr, out_arr, user));
61 PetscCall(PicurvAssertRealNear(3.5, out_arr[1][1][1].x, 1.0e-12, "vector corner->center x preserves constants"));
62 PetscCall(PicurvAssertRealNear(3.5, out_arr[1][1][1].y, 1.0e-12, "vector corner->center y preserves constants"));
63 PetscCall(PicurvAssertRealNear(3.5, out_arr[1][1][1].z, 1.0e-12, "vector corner->center z preserves constants"));
64 PetscCall(DMDAVecRestoreArrayRead(user->fda, input, &in_arr));
65 PetscCall(DMDAVecRestoreArray(user->fda, output, &out_arr));
66
67 PetscCall(VecDestroy(&input));
68 PetscCall(VecDestroy(&output));
69 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
70 PetscFunctionReturn(0);
71}
PetscErrorCode InterpolateFieldFromCornerToCenter_Vector(Cmpnts ***field_arr, Cmpnts ***centfield_arr, UserCtx *user)
Safely interpolate a vector field from corner nodes (from the coordinate DM) to cell centers (from th...
Vec Ucat
Definition variables.h:764
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestSignedDistanceAndClassification()

static PetscErrorCode TestSignedDistanceAndClassification ( void  )
static

Test-local routine.

Definition at line 76 of file test_geometry.c.

77{
78 const Cmpnts v1 = {0.0, 0.0, 0.0};
79 const Cmpnts v2 = {1.0, 0.0, 0.0};
80 const Cmpnts v3 = {1.0, 1.0, 0.0};
81 const Cmpnts v4 = {0.0, 1.0, 0.0};
82 const Cmpnts centroid = {0.5, 0.5, 0.5};
83 const Cmpnts inside_point = {0.5, 0.5, 0.25};
84 const Cmpnts outside_point = {0.5, 0.5, -1.0};
85 PetscReal d_inside = 0.0, d_outside = 0.0;
86 PetscReal distances[NUM_FACES] = {1.0, 1.0, 0.0, 1.0, 1.0, 1.0};
87 PetscInt result = -99;
88
89 PetscFunctionBeginUser;
90 PetscCall(ComputeSignedDistanceToPlane(v1, v2, v3, v4, centroid, inside_point, &d_inside, 1.0e-12));
91 PetscCall(ComputeSignedDistanceToPlane(v1, v2, v3, v4, centroid, outside_point, &d_outside, 1.0e-12));
92 PetscCall(PicurvAssertBool((PetscBool)(d_inside > 0.0), "inside point should yield positive signed distance"));
93 PetscCall(PicurvAssertBool((PetscBool)(d_outside < 0.0), "outside point should yield negative signed distance"));
94
95 PetscCall(DeterminePointPosition(distances, &result));
96 PetscCall(PicurvAssertIntEqual(1, result, "one zero face distance should classify as on-face"));
97
98 distances[0] = -0.1;
99 PetscCall(DeterminePointPosition(distances, &result));
100 PetscCall(PicurvAssertIntEqual(-1, result, "negative face distance should classify as outside"));
101 PetscFunctionReturn(0);
102}
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Shared test-support routine.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Shared test-support routine.
@ NUM_FACES
Definition variables.h:145
PetscErrorCode ComputeSignedDistanceToPlane(const Cmpnts v1, const Cmpnts v2, const Cmpnts v3, const Cmpnts v4, const Cmpnts cell_centroid, const Cmpnts p_target, PetscReal *d_signed, const PetscReal threshold)
Computes the signed distance from a point to the plane approximating a quadrilateral face.
PetscErrorCode DeterminePointPosition(PetscReal *d, PetscInt *result)
Classifies a point based on precomputed face distances.
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 107 of file test_geometry.c.

108{
109 PetscErrorCode ierr;
110 const PicurvTestCase cases[] = {
111 {"scalar-interpolation-constant", TestScalarInterpolationConstantField},
112 {"vector-interpolation-constant", TestVectorInterpolationConstantField},
113 {"signed-distance-and-classification", TestSignedDistanceAndClassification},
114 };
115
116 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv geometry tests");
117 if (ierr) {
118 return (int)ierr;
119 }
120
121 ierr = PicurvRunTests("unit-geometry", cases, sizeof(cases) / sizeof(cases[0]));
122 if (ierr) {
123 PetscFinalize();
124 return (int)ierr;
125 }
126
127 ierr = PetscFinalize();
128 return (int)ierr;
129}
static PetscErrorCode TestSignedDistanceAndClassification(void)
Test-local routine.
static PetscErrorCode TestScalarInterpolationConstantField(void)
Test-local routine.
static PetscErrorCode TestVectorInterpolationConstantField(void)
Test-local routine.
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Shared test-support routine.
Named test case descriptor consumed by PicurvRunTests.
Here is the call graph for this function: