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

C test module for PICurv. More...

#include "test_support.h"
#include "Metric.h"
Include dependency graph for test_metric.c:

Go to the source code of this file.

Functions

static PetscErrorCode TestInvertCovariantMetricTensorDiagonal (void)
 Test-local routine.
 
static PetscErrorCode TestMetricVelocityContravariantIdentity (void)
 Test-local routine.
 
static PetscErrorCode TestMetricVelocityContravariantScaledAxes (void)
 Test-local routine.
 
static PetscErrorCode TestCalculateFaceNormalAndAreaAxisAligned (void)
 Test-local routine.
 
static PetscErrorCode TestComputeCellCharacteristicLengthScaleAxisAligned (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_metric.c.

Function Documentation

◆ TestInvertCovariantMetricTensorDiagonal()

static PetscErrorCode TestInvertCovariantMetricTensorDiagonal ( void  )
static

Test-local routine.

Definition at line 13 of file test_metric.c.

14{
15 double covariant[3][3] = {
16 {4.0, 0.0, 0.0},
17 {0.0, 9.0, 0.0},
18 {0.0, 0.0, 16.0},
19 };
20 double contravariant[3][3] = {{0.0}};
21
22 PetscFunctionBeginUser;
23 PetscCall(InvertCovariantMetricTensor(covariant, contravariant));
24 PetscCall(PicurvAssertRealNear(0.25, contravariant[0][0], 1.0e-12, "inverse(0,0)"));
25 PetscCall(PicurvAssertRealNear(1.0 / 9.0, contravariant[1][1], 1.0e-12, "inverse(1,1)"));
26 PetscCall(PicurvAssertRealNear(0.0625, contravariant[2][2], 1.0e-12, "inverse(2,2)"));
27 PetscFunctionReturn(0);
28}
PetscErrorCode InvertCovariantMetricTensor(double covariantTensor[3][3], double contravariantTensor[3][3])
Inverts the 3x3 covariant metric tensor to obtain the contravariant metric tensor.
Definition Metric.c:203
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Shared test-support routine.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestMetricVelocityContravariantIdentity()

static PetscErrorCode TestMetricVelocityContravariantIdentity ( void  )
static

Test-local routine.

Definition at line 33 of file test_metric.c.

34{
35 const PetscReal J[3][3] = {
36 {1.0, 0.0, 0.0},
37 {0.0, 1.0, 0.0},
38 {0.0, 0.0, 1.0},
39 };
40 const PetscReal u[3] = {1.5, -2.0, 4.0};
41 PetscReal uc[3] = {0.0, 0.0, 0.0};
42
43 PetscFunctionBeginUser;
44 PetscCall(MetricVelocityContravariant(J, 1.0, u, uc));
45 PetscCall(PicurvAssertRealNear(1.5, uc[0], 1.0e-12, "identity map uc[0]"));
46 PetscCall(PicurvAssertRealNear(-2.0, uc[1], 1.0e-12, "identity map uc[1]"));
47 PetscCall(PicurvAssertRealNear(4.0, uc[2], 1.0e-12, "identity map uc[2]"));
48 PetscFunctionReturn(0);
49}
PetscErrorCode MetricVelocityContravariant(const PetscReal J[3][3], PetscReal detJ, const PetscReal u[3], PetscReal uc[3])
Implementation of MetricVelocityContravariant().
Definition Metric.c:166
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestMetricVelocityContravariantScaledAxes()

static PetscErrorCode TestMetricVelocityContravariantScaledAxes ( void  )
static

Test-local routine.

Definition at line 54 of file test_metric.c.

55{
56 const PetscReal J[3][3] = {
57 {2.0, 0.0, 0.0},
58 {0.0, 3.0, 0.0},
59 {0.0, 0.0, 4.0},
60 };
61 const PetscReal u[3] = {2.0, 3.0, 4.0};
62 PetscReal uc[3] = {0.0, 0.0, 0.0};
63
64 PetscFunctionBeginUser;
65 PetscCall(MetricVelocityContravariant(J, 24.0, u, uc));
66 PetscCall(PicurvAssertRealNear(1.0, uc[0], 1.0e-12, "scaled map uc[0]"));
67 PetscCall(PicurvAssertRealNear(1.0, uc[1], 1.0e-12, "scaled map uc[1]"));
68 PetscCall(PicurvAssertRealNear(1.0, uc[2], 1.0e-12, "scaled map uc[2]"));
69 PetscFunctionReturn(0);
70}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestCalculateFaceNormalAndAreaAxisAligned()

static PetscErrorCode TestCalculateFaceNormalAndAreaAxisAligned ( void  )
static

Test-local routine.

Definition at line 75 of file test_metric.c.

76{
77 const Cmpnts csi = {2.0, 0.0, 0.0};
78 const Cmpnts eta = {0.0, 3.0, 0.0};
79 const Cmpnts zet = {0.0, 0.0, 4.0};
80 double ni[3] = {0.0, 0.0, 0.0};
81 double nj[3] = {0.0, 0.0, 0.0};
82 double nk[3] = {0.0, 0.0, 0.0};
83 double Ai = 0.0, Aj = 0.0, Ak = 0.0;
84
85 PetscFunctionBeginUser;
86 PetscCall(CalculateFaceNormalAndArea(csi, eta, zet, ni, nj, nk, &Ai, &Aj, &Ak));
87 PetscCall(PicurvAssertRealNear(2.0, Ai, 1.0e-12, "Ai"));
88 PetscCall(PicurvAssertRealNear(3.0, Aj, 1.0e-12, "Aj"));
89 PetscCall(PicurvAssertRealNear(4.0, Ak, 1.0e-12, "Ak"));
90 PetscCall(PicurvAssertRealNear(1.0, ni[0], 1.0e-12, "ni.x"));
91 PetscCall(PicurvAssertRealNear(0.0, ni[1], 1.0e-12, "ni.y"));
92 PetscCall(PicurvAssertRealNear(0.0, ni[2], 1.0e-12, "ni.z"));
93 PetscCall(PicurvAssertRealNear(0.0, nj[0], 1.0e-12, "nj.x"));
94 PetscCall(PicurvAssertRealNear(1.0, nj[1], 1.0e-12, "nj.y"));
95 PetscCall(PicurvAssertRealNear(0.0, nj[2], 1.0e-12, "nj.z"));
96 PetscCall(PicurvAssertRealNear(0.0, nk[0], 1.0e-12, "nk.x"));
97 PetscCall(PicurvAssertRealNear(0.0, nk[1], 1.0e-12, "nk.y"));
98 PetscCall(PicurvAssertRealNear(1.0, nk[2], 1.0e-12, "nk.z"));
99 PetscFunctionReturn(0);
100}
PetscErrorCode CalculateFaceNormalAndArea(Cmpnts csi, Cmpnts eta, Cmpnts zet, double ni[3], double nj[3], double nk[3], double *Ai, double *Aj, double *Ak)
Computes the unit normal vectors and areas of the three faces of a computational cell.
Definition Metric.c:237
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:

◆ TestComputeCellCharacteristicLengthScaleAxisAligned()

static PetscErrorCode TestComputeCellCharacteristicLengthScaleAxisAligned ( void  )
static

Test-local routine.

Definition at line 105 of file test_metric.c.

106{
107 const Cmpnts csi = {2.0, 0.0, 0.0};
108 const Cmpnts eta = {0.0, 3.0, 0.0};
109 const Cmpnts zet = {0.0, 0.0, 4.0};
110 double dx = 0.0, dy = 0.0, dz = 0.0;
111
112 PetscFunctionBeginUser;
113 PetscCall(ComputeCellCharacteristicLengthScale(1.0, csi, eta, zet, &dx, &dy, &dz));
114 PetscCall(PicurvAssertRealNear(0.5, dx, 1.0e-12, "dx"));
115 PetscCall(PicurvAssertRealNear(1.0 / 3.0, dy, 1.0e-12, "dy"));
116 PetscCall(PicurvAssertRealNear(0.25, dz, 1.0e-12, "dz"));
117 PetscFunctionReturn(0);
118}
PetscErrorCode ComputeCellCharacteristicLengthScale(PetscReal ajc, Cmpnts csi, Cmpnts eta, Cmpnts zet, double *dx, double *dy, double *dz)
Computes characteristic length scales (dx, dy, dz) for a curvilinear cell.
Definition Metric.c:283
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 123 of file test_metric.c.

124{
125 PetscErrorCode ierr;
126 const PicurvTestCase cases[] = {
127 {"invert-covariant-metric-diagonal", TestInvertCovariantMetricTensorDiagonal},
128 {"metric-velocity-contravariant-identity", TestMetricVelocityContravariantIdentity},
129 {"metric-velocity-contravariant-scaled-axes", TestMetricVelocityContravariantScaledAxes},
130 {"face-normal-and-area-axis-aligned", TestCalculateFaceNormalAndAreaAxisAligned},
131 {"characteristic-length-scale-axis-aligned", TestComputeCellCharacteristicLengthScaleAxisAligned},
132 };
133
134 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv metric tests");
135 if (ierr) {
136 return (int)ierr;
137 }
138
139 ierr = PicurvRunTests("unit-metric", cases, sizeof(cases) / sizeof(cases[0]));
140 if (ierr) {
141 PetscFinalize();
142 return (int)ierr;
143 }
144
145 ierr = PetscFinalize();
146 return (int)ierr;
147}
static PetscErrorCode TestMetricVelocityContravariantScaledAxes(void)
Test-local routine.
Definition test_metric.c:54
static PetscErrorCode TestInvertCovariantMetricTensorDiagonal(void)
Test-local routine.
Definition test_metric.c:13
static PetscErrorCode TestMetricVelocityContravariantIdentity(void)
Test-local routine.
Definition test_metric.c:33
static PetscErrorCode TestCalculateFaceNormalAndAreaAxisAligned(void)
Test-local routine.
Definition test_metric.c:75
static PetscErrorCode TestComputeCellCharacteristicLengthScaleAxisAligned(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: