C unit tests for metric tensors, face geometry, and length-scale helpers.
More...
#include "test_support.h"
#include "Metric.h"
Go to the source code of this file.
C unit tests for metric tensors, face geometry, and length-scale helpers.
Definition in file test_metric.c.
◆ TestInvertCovariantMetricTensorDiagonal()
| static PetscErrorCode TestInvertCovariantMetricTensorDiagonal |
( |
void |
| ) |
|
|
static |
Tests inversion of a diagonal covariant metric tensor.
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;
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.
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Asserts that two real values agree within tolerance.
◆ TestMetricVelocityContravariantIdentity()
| static PetscErrorCode TestMetricVelocityContravariantIdentity |
( |
void |
| ) |
|
|
static |
Tests contravariant velocity computation for identity metrics.
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;
48 PetscFunctionReturn(0);
49}
PetscErrorCode MetricVelocityContravariant(const PetscReal J[3][3], PetscReal detJ, const PetscReal u[3], PetscReal uc[3])
Public interface for MetricVelocityContravariant().
◆ TestMetricVelocityContravariantScaledAxes()
| static PetscErrorCode TestMetricVelocityContravariantScaledAxes |
( |
void |
| ) |
|
|
static |
Tests contravariant velocity computation for scaled axes.
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;
69 PetscFunctionReturn(0);
70}
◆ TestCalculateFaceNormalAndAreaAxisAligned()
| static PetscErrorCode TestCalculateFaceNormalAndAreaAxisAligned |
( |
void |
| ) |
|
|
static |
Tests face-normal and face-area calculation on an axis-aligned cell.
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;
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.
A 3D point or vector with PetscScalar components.
◆ TestComputeCellCharacteristicLengthScaleAxisAligned()
| static PetscErrorCode TestComputeCellCharacteristicLengthScaleAxisAligned |
( |
void |
| ) |
|
|
static |
Tests characteristic length-scale computation on an axis-aligned cell.
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;
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.
◆ main()
| int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Runs the unit-metric PETSc test binary.
Definition at line 123 of file test_metric.c.
124{
125 PetscErrorCode ierr;
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)
Tests contravariant velocity computation for scaled axes.
static PetscErrorCode TestInvertCovariantMetricTensorDiagonal(void)
Tests inversion of a diagonal covariant metric tensor.
static PetscErrorCode TestMetricVelocityContravariantIdentity(void)
Tests contravariant velocity computation for identity metrics.
static PetscErrorCode TestCalculateFaceNormalAndAreaAxisAligned(void)
Tests face-normal and face-area calculation on an axis-aligned cell.
static PetscErrorCode TestComputeCellCharacteristicLengthScaleAxisAligned(void)
Tests characteristic length-scale computation on an axis-aligned cell.
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.