xref: /petsc/src/dm/dt/tests/ex9.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1fbdc3dfeSToby Isaac const char help[] = "Tests PetscDTPKDEvalJet()";
2fbdc3dfeSToby Isaac 
3fbdc3dfeSToby Isaac #include <petscdt.h>
4fbdc3dfeSToby Isaac #include <petscblaslapack.h>
5fbdc3dfeSToby Isaac 
testOrthogonality(PetscInt dim,PetscInt deg)6d71ae5a4SJacob Faibussowitsch static PetscErrorCode testOrthogonality(PetscInt dim, PetscInt deg)
7d71ae5a4SJacob Faibussowitsch {
8fbdc3dfeSToby Isaac   PetscQuadrature  q;
9fbdc3dfeSToby Isaac   const PetscReal *points, *weights;
10fbdc3dfeSToby Isaac   PetscInt         Npoly, npoints, i, j, k;
11fbdc3dfeSToby Isaac   PetscReal       *p;
12fbdc3dfeSToby Isaac 
13fbdc3dfeSToby Isaac   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q));
159566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights));
169566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + deg, dim, &Npoly));
179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Npoly * npoints, &p));
189566063dSJacob Faibussowitsch   PetscCall(PetscDTPKDEvalJet(dim, npoints, points, deg, 0, p));
19fbdc3dfeSToby Isaac   for (i = 0; i < Npoly; i++) {
20fbdc3dfeSToby Isaac     for (j = i; j < Npoly; j++) {
21fbdc3dfeSToby Isaac       PetscReal integral = 0.;
22fbdc3dfeSToby Isaac       PetscReal exact    = (i == j) ? 1. : 0.;
23fbdc3dfeSToby Isaac 
24fbdc3dfeSToby Isaac       for (k = 0; k < npoints; k++) integral += weights[k] * p[i * npoints + k] * p[j * npoints + k];
2563a3b9bcSJacob Faibussowitsch       PetscCheck(PetscAbsReal(integral - exact) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "<P[%" PetscInt_FMT "], P[%" PetscInt_FMT "]> = %g != delta_{%" PetscInt_FMT ",%" PetscInt_FMT "}", i, j, (double)integral, i, j);
26fbdc3dfeSToby Isaac     }
27fbdc3dfeSToby Isaac   }
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(p));
299566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&q));
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31fbdc3dfeSToby Isaac }
32fbdc3dfeSToby Isaac 
testDerivativesLegendre(PetscInt dim,PetscInt deg,PetscInt k)33d71ae5a4SJacob Faibussowitsch static PetscErrorCode testDerivativesLegendre(PetscInt dim, PetscInt deg, PetscInt k)
34d71ae5a4SJacob Faibussowitsch {
35fbdc3dfeSToby Isaac   PetscInt         Np, Nk, i, j, l, d, npoints;
36fbdc3dfeSToby Isaac   PetscRandom      rand;
37fbdc3dfeSToby Isaac   PetscReal       *point;
38fbdc3dfeSToby Isaac   PetscReal       *lgndre_coeffs;
39fbdc3dfeSToby Isaac   PetscReal       *pkd_coeffs;
40fbdc3dfeSToby Isaac   PetscReal       *proj;
41fbdc3dfeSToby Isaac   PetscReal      **B;
42fbdc3dfeSToby Isaac   PetscQuadrature  q;
43fbdc3dfeSToby Isaac   PetscReal       *points1d;
44fbdc3dfeSToby Isaac   PetscInt        *degrees;
45fbdc3dfeSToby Isaac   PetscInt        *degtup, *ktup;
46fbdc3dfeSToby Isaac   const PetscReal *points;
47fbdc3dfeSToby Isaac   const PetscReal *weights;
48fbdc3dfeSToby Isaac   PetscReal       *lgndre_jet;
49fbdc3dfeSToby Isaac   PetscReal      **D;
50fbdc3dfeSToby Isaac   PetscReal       *pkd_jet, *pkd_jet_basis;
51fbdc3dfeSToby Isaac 
52fbdc3dfeSToby Isaac   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + deg, dim, &Np));
549566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + k, dim, &Nk));
55fbdc3dfeSToby Isaac 
56fbdc3dfeSToby Isaac   /* create the projector (because it is an orthonormal basis, the projector is the moment integrals) */
579566063dSJacob Faibussowitsch   PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q));
589566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights));
599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints * Np, &proj));
609566063dSJacob Faibussowitsch   PetscCall(PetscDTPKDEvalJet(dim, npoints, points, deg, 0, proj));
619371c9d4SSatish Balay   for (i = 0; i < Np; i++)
629371c9d4SSatish Balay     for (j = 0; j < npoints; j++) proj[i * npoints + j] *= weights[j];
63fbdc3dfeSToby Isaac 
649566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
659566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand, -1., 1.));
66fbdc3dfeSToby Isaac 
67fbdc3dfeSToby Isaac   /* create a random coefficient vector */
689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Np, &lgndre_coeffs, Np, &pkd_coeffs));
6948a46eb9SPierre Jolivet   for (i = 0; i < Np; i++) PetscCall(PetscRandomGetValueReal(rand, &lgndre_coeffs[i]));
70fbdc3dfeSToby Isaac 
719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dim, &degtup, dim, &ktup));
729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(deg + 1, &degrees));
73fbdc3dfeSToby Isaac   for (i = 0; i < deg + 1; i++) degrees[i] = i;
74fbdc3dfeSToby Isaac 
75fbdc3dfeSToby Isaac   /* project the lgndre_coeffs to pkd_coeffs */
769566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pkd_coeffs, Np));
779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &points1d));
789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &B));
79fbdc3dfeSToby Isaac   for (d = 0; d < dim; d++) {
80f4f49eeaSPierre Jolivet     PetscCall(PetscMalloc1((deg + 1) * npoints, &B[d]));
81fbdc3dfeSToby Isaac     /* get this coordinate */
82fbdc3dfeSToby Isaac     for (i = 0; i < npoints; i++) points1d[i] = points[i * dim + d];
839566063dSJacob Faibussowitsch     PetscCall(PetscDTLegendreEval(npoints, points1d, deg + 1, degrees, B[d], NULL, NULL));
84fbdc3dfeSToby Isaac   }
859566063dSJacob Faibussowitsch   PetscCall(PetscFree(points1d));
86fbdc3dfeSToby Isaac   for (i = 0; i < npoints; i++) {
87fbdc3dfeSToby Isaac     PetscReal val = 0.;
88fbdc3dfeSToby Isaac 
89fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
90fbdc3dfeSToby Isaac       PetscReal mul  = lgndre_coeffs[j];
91fbdc3dfeSToby Isaac       PetscReal valj = 1.;
92fbdc3dfeSToby Isaac 
939566063dSJacob Faibussowitsch       PetscCall(PetscDTIndexToGradedOrder(dim, j, degtup));
94ad540459SPierre Jolivet       for (l = 0; l < dim; l++) valj *= B[l][i * (deg + 1) + degtup[l]];
95fbdc3dfeSToby Isaac       val += mul * valj;
96fbdc3dfeSToby Isaac     }
97ad540459SPierre Jolivet     for (j = 0; j < Np; j++) pkd_coeffs[j] += proj[j * npoints + i] * val;
98fbdc3dfeSToby Isaac   }
9948a46eb9SPierre Jolivet   for (i = 0; i < dim; i++) PetscCall(PetscFree(B[i]));
1009566063dSJacob Faibussowitsch   PetscCall(PetscFree(B));
101fbdc3dfeSToby Isaac 
102fbdc3dfeSToby Isaac   /* create a random point in the biunit simplex */
1039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &point));
10448a46eb9SPierre Jolivet   for (i = 0; i < dim; i++) PetscCall(PetscRandomGetValueReal(rand, &point[i]));
105fbdc3dfeSToby Isaac   for (i = dim - 1; i > 0; i--) {
106fbdc3dfeSToby Isaac     PetscReal val = point[i];
107fbdc3dfeSToby Isaac     PetscInt  j;
108fbdc3dfeSToby Isaac 
109ad540459SPierre Jolivet     for (j = 0; j < i; j++) point[j] = (point[j] + 1.) * (1. - val) * 0.5 - 1.;
110fbdc3dfeSToby Isaac   }
111fbdc3dfeSToby Isaac 
1129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nk * Np, &pkd_jet_basis, Nk, &lgndre_jet, Nk, &pkd_jet));
113fbdc3dfeSToby Isaac   /* evaluate the jet at the point with PKD polynomials */
1149566063dSJacob Faibussowitsch   PetscCall(PetscDTPKDEvalJet(dim, 1, point, deg, k, pkd_jet_basis));
115fbdc3dfeSToby Isaac   for (i = 0; i < Nk; i++) {
116fbdc3dfeSToby Isaac     PetscReal val = 0.;
117ad540459SPierre Jolivet     for (j = 0; j < Np; j++) val += pkd_coeffs[j] * pkd_jet_basis[j * Nk + i];
118fbdc3dfeSToby Isaac     pkd_jet[i] = val;
119fbdc3dfeSToby Isaac   }
120fbdc3dfeSToby Isaac 
121fbdc3dfeSToby Isaac   /* evaluate the 1D jets of the Legendre polynomials */
1229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &D));
123fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
124f4f49eeaSPierre Jolivet     PetscCall(PetscMalloc1((deg + 1) * (k + 1), &D[i]));
125f4f49eeaSPierre Jolivet     PetscCall(PetscDTJacobiEvalJet(0., 0., 1, &point[i], deg, k, D[i]));
126fbdc3dfeSToby Isaac   }
127fbdc3dfeSToby Isaac   /* compile the 1D Legendre jets into the tensor Legendre jet */
128fbdc3dfeSToby Isaac   for (j = 0; j < Nk; j++) lgndre_jet[j] = 0.;
129fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) {
130fbdc3dfeSToby Isaac     PetscReal mul = lgndre_coeffs[i];
131fbdc3dfeSToby Isaac 
1329566063dSJacob Faibussowitsch     PetscCall(PetscDTIndexToGradedOrder(dim, i, degtup));
133fbdc3dfeSToby Isaac     for (j = 0; j < Nk; j++) {
134fbdc3dfeSToby Isaac       PetscReal val = 1.;
135fbdc3dfeSToby Isaac 
1369566063dSJacob Faibussowitsch       PetscCall(PetscDTIndexToGradedOrder(dim, j, ktup));
137ad540459SPierre Jolivet       for (l = 0; l < dim; l++) val *= D[l][degtup[l] * (k + 1) + ktup[l]];
138fbdc3dfeSToby Isaac       lgndre_jet[j] += mul * val;
139fbdc3dfeSToby Isaac     }
140fbdc3dfeSToby Isaac   }
14148a46eb9SPierre Jolivet   for (i = 0; i < dim; i++) PetscCall(PetscFree(D[i]));
1429566063dSJacob Faibussowitsch   PetscCall(PetscFree(D));
143fbdc3dfeSToby Isaac 
144fbdc3dfeSToby Isaac   for (i = 0; i < Nk; i++) {
145fbdc3dfeSToby Isaac     PetscReal diff  = lgndre_jet[i] - pkd_jet[i];
146fbdc3dfeSToby Isaac     PetscReal scale = 1. + PetscAbsReal(lgndre_jet[i]) + PetscAbsReal(pkd_jet[i]);
147fbdc3dfeSToby Isaac     PetscReal tol   = 10. * PETSC_SMALL * scale;
148fbdc3dfeSToby Isaac 
14908401ef6SPierre Jolivet     PetscCheck(PetscAbsReal(diff) <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Jet mismatch between PKD and tensor Legendre bases: error %g at tolerance %g", (double)diff, (double)tol);
150fbdc3dfeSToby Isaac   }
151fbdc3dfeSToby Isaac 
1529566063dSJacob Faibussowitsch   PetscCall(PetscFree2(degtup, ktup));
1539566063dSJacob Faibussowitsch   PetscCall(PetscFree(degrees));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree3(pkd_jet_basis, lgndre_jet, pkd_jet));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(point));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lgndre_coeffs, pkd_coeffs));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(proj));
1589566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
1599566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&q));
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161fbdc3dfeSToby Isaac }
162fbdc3dfeSToby Isaac 
main(int argc,char ** argv)163d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
164d71ae5a4SJacob Faibussowitsch {
165fbdc3dfeSToby Isaac   PetscInt dim, deg, k;
166fbdc3dfeSToby Isaac 
167327415f7SBarry Smith   PetscFunctionBeginUser;
1689566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
169fbdc3dfeSToby Isaac   dim = 3;
170fbdc3dfeSToby Isaac   deg = 4;
171fbdc3dfeSToby Isaac   k   = 3;
172d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PetscDTPKDEval() tests", "none");
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dim", "Dimension of the simplex", "ex9.c", dim, &dim, NULL));
1749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-degree", "The degree of the polynomial space", "ex9.c", deg, &deg, NULL));
1759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-k", "The number of derivatives to use in the taylor test", "ex9.c", k, &k, NULL));
176d0609cedSBarry Smith   PetscOptionsEnd();
1779566063dSJacob Faibussowitsch   PetscCall(testOrthogonality(dim, deg));
1789566063dSJacob Faibussowitsch   PetscCall(testDerivativesLegendre(dim, deg, k));
1799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
180b122ec5aSJacob Faibussowitsch   return 0;
181fbdc3dfeSToby Isaac }
182fbdc3dfeSToby Isaac 
183fbdc3dfeSToby Isaac /*TEST
184fbdc3dfeSToby Isaac 
185fbdc3dfeSToby Isaac   test:
186fbdc3dfeSToby Isaac     args: -dim {{1 2 3 4}}
187*3886731fSPierre Jolivet     output_file: output/empty.out
188fbdc3dfeSToby Isaac 
189fbdc3dfeSToby Isaac TEST*/
190