xref: /petsc/src/dm/dt/tests/ex9.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1fbdc3dfeSToby Isaac const char help[] = "Tests PetscDTPKDEvalJet()";
2fbdc3dfeSToby Isaac 
3fbdc3dfeSToby Isaac #include <petscdt.h>
4fbdc3dfeSToby Isaac #include <petscblaslapack.h>
5fbdc3dfeSToby Isaac 
6fbdc3dfeSToby Isaac static PetscErrorCode testOrthogonality(PetscInt dim, PetscInt deg)
7fbdc3dfeSToby Isaac {
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;
145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim + deg, dim, &Npoly));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Npoly * npoints, &p));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(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];
252c71b3e2SJacob Faibussowitsch       PetscCheckFalse(PetscAbsReal(integral - exact) > PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_PLIB, "<P[%D], P[%D]> = %g != delta_{%D,%D}", i, j, (double) integral, i, j);
26fbdc3dfeSToby Isaac     }
27fbdc3dfeSToby Isaac   }
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(p));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&q));
30fbdc3dfeSToby Isaac   PetscFunctionReturn(0);
31fbdc3dfeSToby Isaac }
32fbdc3dfeSToby Isaac 
33fbdc3dfeSToby Isaac static PetscErrorCode testDerivativesLegendre(PetscInt dim, PetscInt deg, PetscInt k)
34fbdc3dfeSToby Isaac {
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;
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim + deg, dim, &Np));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim + k, dim, &Nk));
55fbdc3dfeSToby Isaac 
56fbdc3dfeSToby Isaac   /* create the projector (because it is an orthonormal basis, the projector is the moment integrals) */
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints * Np, &proj));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTPKDEvalJet(dim, npoints, points, deg, 0, proj));
61fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) for (j = 0; j < npoints; j++) proj[i * npoints + j] *= weights[j];
62fbdc3dfeSToby Isaac 
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rand));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rand, -1., 1.));
65fbdc3dfeSToby Isaac 
66fbdc3dfeSToby Isaac   /* create a random coefficient vector */
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(Np, &lgndre_coeffs, Np, &pkd_coeffs));
68fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) {
695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValueReal(rand, &lgndre_coeffs[i]));
70fbdc3dfeSToby Isaac   }
71fbdc3dfeSToby Isaac 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(dim, &degtup, dim, &ktup));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(deg + 1, &degrees));
74fbdc3dfeSToby Isaac   for (i = 0; i < deg + 1; i++) degrees[i] = i;
75fbdc3dfeSToby Isaac 
76fbdc3dfeSToby Isaac   /* project the lgndre_coeffs to pkd_coeffs */
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(pkd_coeffs, Np));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(npoints, &points1d));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dim, &B));
80fbdc3dfeSToby Isaac   for (d = 0; d < dim; d++) {
815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1((deg + 1)*npoints, &(B[d])));
82fbdc3dfeSToby Isaac     /* get this coordinate */
83fbdc3dfeSToby Isaac     for (i = 0; i < npoints; i++) points1d[i] = points[i * dim + d];
845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTLegendreEval(npoints, points1d, deg + 1, degrees, B[d], NULL, NULL));
85fbdc3dfeSToby Isaac   }
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(points1d));
87fbdc3dfeSToby Isaac   for (i = 0; i < npoints; i++) {
88fbdc3dfeSToby Isaac     PetscReal val = 0.;
89fbdc3dfeSToby Isaac 
90fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
91fbdc3dfeSToby Isaac       PetscReal mul = lgndre_coeffs[j];
92fbdc3dfeSToby Isaac       PetscReal valj = 1.;
93fbdc3dfeSToby Isaac 
945f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDTIndexToGradedOrder(dim, j, degtup));
95fbdc3dfeSToby Isaac       for (l = 0; l < dim; l++) {
96fbdc3dfeSToby Isaac         valj *= B[l][i * (deg + 1) + degtup[l]];
97fbdc3dfeSToby Isaac       }
98fbdc3dfeSToby Isaac       val += mul * valj;
99fbdc3dfeSToby Isaac     }
100fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
101fbdc3dfeSToby Isaac       pkd_coeffs[j] += proj[j * npoints + i] * val;
102fbdc3dfeSToby Isaac     }
103fbdc3dfeSToby Isaac   }
104fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
1055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(B[i]));
106fbdc3dfeSToby Isaac   }
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(B));
108fbdc3dfeSToby Isaac 
109fbdc3dfeSToby Isaac   /* create a random point in the biunit simplex */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dim, &point));
111fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValueReal(rand, &point[i]));
113fbdc3dfeSToby Isaac   }
114fbdc3dfeSToby Isaac   for (i = dim - 1; i > 0; i--) {
115fbdc3dfeSToby Isaac     PetscReal val = point[i];
116fbdc3dfeSToby Isaac     PetscInt  j;
117fbdc3dfeSToby Isaac 
118fbdc3dfeSToby Isaac     for (j = 0; j < i; j++) {
119fbdc3dfeSToby Isaac       point[j] = (point[j] + 1.)*(1. - val)*0.5 - 1.;
120fbdc3dfeSToby Isaac     }
121fbdc3dfeSToby Isaac   }
122fbdc3dfeSToby Isaac 
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(Nk*Np, &pkd_jet_basis, Nk, &lgndre_jet, Nk, &pkd_jet));
124fbdc3dfeSToby Isaac   /* evaluate the jet at the point with PKD polynomials */
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTPKDEvalJet(dim, 1, point, deg, k, pkd_jet_basis));
126fbdc3dfeSToby Isaac   for (i = 0; i < Nk; i++) {
127fbdc3dfeSToby Isaac     PetscReal val = 0.;
128fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
129fbdc3dfeSToby Isaac       val += pkd_coeffs[j] * pkd_jet_basis[j * Nk + i];
130fbdc3dfeSToby Isaac     }
131fbdc3dfeSToby Isaac     pkd_jet[i] = val;
132fbdc3dfeSToby Isaac   }
133fbdc3dfeSToby Isaac 
134fbdc3dfeSToby Isaac   /* evaluate the 1D jets of the Legendre polynomials */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(dim, &D));
136fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1((deg + 1) * (k+1), &(D[i])));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTJacobiEvalJet(0., 0., 1, &(point[i]), deg, k, D[i]));
139fbdc3dfeSToby Isaac   }
140fbdc3dfeSToby Isaac   /* compile the 1D Legendre jets into the tensor Legendre jet */
141fbdc3dfeSToby Isaac   for (j = 0; j < Nk; j++) lgndre_jet[j] = 0.;
142fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) {
143fbdc3dfeSToby Isaac     PetscReal mul = lgndre_coeffs[i];
144fbdc3dfeSToby Isaac 
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTIndexToGradedOrder(dim, i, degtup));
146fbdc3dfeSToby Isaac     for (j = 0; j < Nk; j++) {
147fbdc3dfeSToby Isaac       PetscReal val = 1.;
148fbdc3dfeSToby Isaac 
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDTIndexToGradedOrder(dim, j, ktup));
150fbdc3dfeSToby Isaac       for (l = 0; l < dim; l++) {
151fbdc3dfeSToby Isaac         val *= D[l][degtup[l]*(k+1) + ktup[l]];
152fbdc3dfeSToby Isaac       }
153fbdc3dfeSToby Isaac       lgndre_jet[j] += mul * val;
154fbdc3dfeSToby Isaac     }
155fbdc3dfeSToby Isaac   }
156fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(D[i]));
158fbdc3dfeSToby Isaac   }
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(D));
160fbdc3dfeSToby Isaac 
161fbdc3dfeSToby Isaac   for (i = 0; i < Nk; i++) {
162fbdc3dfeSToby Isaac     PetscReal diff = lgndre_jet[i] - pkd_jet[i];
163fbdc3dfeSToby Isaac     PetscReal scale = 1. + PetscAbsReal(lgndre_jet[i]) + PetscAbsReal(pkd_jet[i]);
164fbdc3dfeSToby Isaac     PetscReal tol = 10. * PETSC_SMALL * scale;
165fbdc3dfeSToby Isaac 
1662c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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);
167fbdc3dfeSToby Isaac   }
168fbdc3dfeSToby Isaac 
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(degtup,ktup));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(degrees));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(pkd_jet_basis, lgndre_jet, pkd_jet));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(point));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(lgndre_coeffs, pkd_coeffs));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(proj));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&q));
177fbdc3dfeSToby Isaac   PetscFunctionReturn(0);
178fbdc3dfeSToby Isaac }
179fbdc3dfeSToby Isaac 
180fbdc3dfeSToby Isaac int main(int argc, char **argv)
181fbdc3dfeSToby Isaac {
182fbdc3dfeSToby Isaac   PetscInt       dim, deg, k;
183fbdc3dfeSToby Isaac   PetscErrorCode ierr;
184fbdc3dfeSToby Isaac 
185*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
186fbdc3dfeSToby Isaac   dim = 3;
187fbdc3dfeSToby Isaac   deg = 4;
188fbdc3dfeSToby Isaac   k = 3;
189fbdc3dfeSToby Isaac   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PetscDTPKDEval() tests","none");CHKERRQ(ierr);
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-dim", "Dimension of the simplex","ex9.c",dim,&dim,NULL));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-degree", "The degree of the polynomial space","ex9.c",deg,&deg,NULL));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-k", "The number of derivatives to use in the taylor test","ex9.c",k,&k,NULL));
1931e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(testOrthogonality(dim, deg));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(testDerivativesLegendre(dim, deg, k));
196*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
197*b122ec5aSJacob Faibussowitsch   return 0;
198fbdc3dfeSToby Isaac }
199fbdc3dfeSToby Isaac 
200fbdc3dfeSToby Isaac /*TEST
201fbdc3dfeSToby Isaac 
202fbdc3dfeSToby Isaac   test:
203fbdc3dfeSToby Isaac     args: -dim {{1 2 3 4}}
204fbdc3dfeSToby Isaac 
205fbdc3dfeSToby Isaac TEST*/
206