xref: /petsc/src/dm/dt/tests/ex9.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode  ierr;
13fbdc3dfeSToby Isaac 
14fbdc3dfeSToby Isaac   PetscFunctionBegin;
15fbdc3dfeSToby Isaac   ierr = PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q);CHKERRQ(ierr);
16fbdc3dfeSToby Isaac   ierr = PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights);CHKERRQ(ierr);
17fbdc3dfeSToby Isaac   ierr = PetscDTBinomialInt(dim + deg, dim, &Npoly);CHKERRQ(ierr);
18fbdc3dfeSToby Isaac   ierr = PetscMalloc1(Npoly * npoints, &p);CHKERRQ(ierr);
19fbdc3dfeSToby Isaac   ierr = PetscDTPKDEvalJet(dim, npoints, points, deg, 0, p);CHKERRQ(ierr);
20fbdc3dfeSToby Isaac   for (i = 0; i < Npoly; i++) {
21fbdc3dfeSToby Isaac     for (j = i; j < Npoly; j++) {
22fbdc3dfeSToby Isaac       PetscReal integral = 0.;
23fbdc3dfeSToby Isaac       PetscReal exact = (i == j) ? 1. : 0.;
24fbdc3dfeSToby Isaac 
25fbdc3dfeSToby Isaac       for (k = 0; k < npoints; k++) integral += weights[k] * p[i * npoints + k] * p[j * npoints + k];
26*2c71b3e2SJacob 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);
27fbdc3dfeSToby Isaac     }
28fbdc3dfeSToby Isaac   }
29fbdc3dfeSToby Isaac   ierr = PetscFree(p);CHKERRQ(ierr);
30fbdc3dfeSToby Isaac   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
31fbdc3dfeSToby Isaac   PetscFunctionReturn(0);
32fbdc3dfeSToby Isaac }
33fbdc3dfeSToby Isaac 
34fbdc3dfeSToby Isaac static PetscErrorCode testDerivativesLegendre(PetscInt dim, PetscInt deg, PetscInt k)
35fbdc3dfeSToby Isaac {
36fbdc3dfeSToby Isaac   PetscInt       Np, Nk, i, j, l, d, npoints;
37fbdc3dfeSToby Isaac   PetscRandom    rand;
38fbdc3dfeSToby Isaac   PetscReal      *point;
39fbdc3dfeSToby Isaac   PetscReal      *lgndre_coeffs;
40fbdc3dfeSToby Isaac   PetscReal      *pkd_coeffs;
41fbdc3dfeSToby Isaac   PetscReal      *proj;
42fbdc3dfeSToby Isaac   PetscReal     **B;
43fbdc3dfeSToby Isaac   PetscQuadrature q;
44fbdc3dfeSToby Isaac   PetscReal       *points1d;
45fbdc3dfeSToby Isaac   PetscInt        *degrees;
46fbdc3dfeSToby Isaac   PetscInt        *degtup, *ktup;
47fbdc3dfeSToby Isaac   const PetscReal *points;
48fbdc3dfeSToby Isaac   const PetscReal *weights;
49fbdc3dfeSToby Isaac   PetscReal      *lgndre_jet;
50fbdc3dfeSToby Isaac   PetscReal     **D;
51fbdc3dfeSToby Isaac   PetscReal      *pkd_jet, *pkd_jet_basis;
52fbdc3dfeSToby Isaac   PetscErrorCode ierr;
53fbdc3dfeSToby Isaac 
54fbdc3dfeSToby Isaac   PetscFunctionBegin;
55fbdc3dfeSToby Isaac   ierr = PetscDTBinomialInt(dim + deg, dim, &Np);CHKERRQ(ierr);
56fbdc3dfeSToby Isaac   ierr = PetscDTBinomialInt(dim + k, dim, &Nk);CHKERRQ(ierr);
57fbdc3dfeSToby Isaac 
58fbdc3dfeSToby Isaac   /* create the projector (because it is an orthonormal basis, the projector is the moment integrals) */
59fbdc3dfeSToby Isaac   ierr = PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q);CHKERRQ(ierr);
60fbdc3dfeSToby Isaac   ierr = PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights);CHKERRQ(ierr);
61fbdc3dfeSToby Isaac   ierr = PetscMalloc1(npoints * Np, &proj);CHKERRQ(ierr);
62fbdc3dfeSToby Isaac   ierr = PetscDTPKDEvalJet(dim, npoints, points, deg, 0, proj);CHKERRQ(ierr);
63fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) for (j = 0; j < npoints; j++) proj[i * npoints + j] *= weights[j];
64fbdc3dfeSToby Isaac 
65fbdc3dfeSToby Isaac   ierr = PetscRandomCreate(PETSC_COMM_SELF, &rand);CHKERRQ(ierr);
66fbdc3dfeSToby Isaac   ierr = PetscRandomSetInterval(rand, -1., 1.);CHKERRQ(ierr);
67fbdc3dfeSToby Isaac 
68fbdc3dfeSToby Isaac   /* create a random coefficient vector */
69fbdc3dfeSToby Isaac   ierr = PetscMalloc2(Np, &lgndre_coeffs, Np, &pkd_coeffs);CHKERRQ(ierr);
70fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) {
71fbdc3dfeSToby Isaac     ierr = PetscRandomGetValueReal(rand, &lgndre_coeffs[i]);CHKERRQ(ierr);
72fbdc3dfeSToby Isaac   }
73fbdc3dfeSToby Isaac 
74fbdc3dfeSToby Isaac   ierr = PetscMalloc2(dim, &degtup, dim, &ktup);CHKERRQ(ierr);
75fbdc3dfeSToby Isaac   ierr = PetscMalloc1(deg + 1, &degrees);CHKERRQ(ierr);
76fbdc3dfeSToby Isaac   for (i = 0; i < deg + 1; i++) degrees[i] = i;
77fbdc3dfeSToby Isaac 
78fbdc3dfeSToby Isaac   /* project the lgndre_coeffs to pkd_coeffs */
79fbdc3dfeSToby Isaac   ierr = PetscArrayzero(pkd_coeffs, Np);CHKERRQ(ierr);
80fbdc3dfeSToby Isaac   ierr = PetscMalloc1(npoints, &points1d);CHKERRQ(ierr);
81fbdc3dfeSToby Isaac   ierr = PetscMalloc1(dim, &B);CHKERRQ(ierr);
82fbdc3dfeSToby Isaac   for (d = 0; d < dim; d++) {
83fbdc3dfeSToby Isaac     ierr = PetscMalloc1((deg + 1)*npoints, &(B[d]));CHKERRQ(ierr);
84fbdc3dfeSToby Isaac     /* get this coordinate */
85fbdc3dfeSToby Isaac     for (i = 0; i < npoints; i++) points1d[i] = points[i * dim + d];
86fbdc3dfeSToby Isaac     ierr = PetscDTLegendreEval(npoints, points1d, deg + 1, degrees, B[d], NULL, NULL);CHKERRQ(ierr);
87fbdc3dfeSToby Isaac   }
88fbdc3dfeSToby Isaac   ierr = PetscFree(points1d);CHKERRQ(ierr);
89fbdc3dfeSToby Isaac   for (i = 0; i < npoints; i++) {
90fbdc3dfeSToby Isaac     PetscReal val = 0.;
91fbdc3dfeSToby Isaac 
92fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
93fbdc3dfeSToby Isaac       PetscReal mul = lgndre_coeffs[j];
94fbdc3dfeSToby Isaac       PetscReal valj = 1.;
95fbdc3dfeSToby Isaac 
96fbdc3dfeSToby Isaac       ierr = PetscDTIndexToGradedOrder(dim, j, degtup);CHKERRQ(ierr);
97fbdc3dfeSToby Isaac       for (l = 0; l < dim; l++) {
98fbdc3dfeSToby Isaac         valj *= B[l][i * (deg + 1) + degtup[l]];
99fbdc3dfeSToby Isaac       }
100fbdc3dfeSToby Isaac       val += mul * valj;
101fbdc3dfeSToby Isaac     }
102fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
103fbdc3dfeSToby Isaac       pkd_coeffs[j] += proj[j * npoints + i] * val;
104fbdc3dfeSToby Isaac     }
105fbdc3dfeSToby Isaac   }
106fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
107fbdc3dfeSToby Isaac     ierr = PetscFree(B[i]);CHKERRQ(ierr);
108fbdc3dfeSToby Isaac   }
109fbdc3dfeSToby Isaac   ierr = PetscFree(B);CHKERRQ(ierr);
110fbdc3dfeSToby Isaac 
111fbdc3dfeSToby Isaac   /* create a random point in the biunit simplex */
112fbdc3dfeSToby Isaac   ierr = PetscMalloc1(dim, &point);CHKERRQ(ierr);
113fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
114fbdc3dfeSToby Isaac     ierr = PetscRandomGetValueReal(rand, &point[i]);CHKERRQ(ierr);
115fbdc3dfeSToby Isaac   }
116fbdc3dfeSToby Isaac   for (i = dim - 1; i > 0; i--) {
117fbdc3dfeSToby Isaac     PetscReal val = point[i];
118fbdc3dfeSToby Isaac     PetscInt  j;
119fbdc3dfeSToby Isaac 
120fbdc3dfeSToby Isaac     for (j = 0; j < i; j++) {
121fbdc3dfeSToby Isaac       point[j] = (point[j] + 1.)*(1. - val)*0.5 - 1.;
122fbdc3dfeSToby Isaac     }
123fbdc3dfeSToby Isaac   }
124fbdc3dfeSToby Isaac 
125fbdc3dfeSToby Isaac   ierr = PetscMalloc3(Nk*Np, &pkd_jet_basis, Nk, &lgndre_jet, Nk, &pkd_jet);CHKERRQ(ierr);
126fbdc3dfeSToby Isaac   /* evaluate the jet at the point with PKD polynomials */
127fbdc3dfeSToby Isaac   ierr = PetscDTPKDEvalJet(dim, 1, point, deg, k, pkd_jet_basis);CHKERRQ(ierr);
128fbdc3dfeSToby Isaac   for (i = 0; i < Nk; i++) {
129fbdc3dfeSToby Isaac     PetscReal val = 0.;
130fbdc3dfeSToby Isaac     for (j = 0; j < Np; j++) {
131fbdc3dfeSToby Isaac       val += pkd_coeffs[j] * pkd_jet_basis[j * Nk + i];
132fbdc3dfeSToby Isaac     }
133fbdc3dfeSToby Isaac     pkd_jet[i] = val;
134fbdc3dfeSToby Isaac   }
135fbdc3dfeSToby Isaac 
136fbdc3dfeSToby Isaac   /* evaluate the 1D jets of the Legendre polynomials */
137fbdc3dfeSToby Isaac   ierr = PetscMalloc1(dim, &D);CHKERRQ(ierr);
138fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
139fbdc3dfeSToby Isaac     ierr = PetscMalloc1((deg + 1) * (k+1), &(D[i]));CHKERRQ(ierr);
140fbdc3dfeSToby Isaac     ierr = PetscDTJacobiEvalJet(0., 0., 1, &(point[i]), deg, k, D[i]);CHKERRQ(ierr);
141fbdc3dfeSToby Isaac   }
142fbdc3dfeSToby Isaac   /* compile the 1D Legendre jets into the tensor Legendre jet */
143fbdc3dfeSToby Isaac   for (j = 0; j < Nk; j++) lgndre_jet[j] = 0.;
144fbdc3dfeSToby Isaac   for (i = 0; i < Np; i++) {
145fbdc3dfeSToby Isaac     PetscReal mul = lgndre_coeffs[i];
146fbdc3dfeSToby Isaac 
147fbdc3dfeSToby Isaac     ierr = PetscDTIndexToGradedOrder(dim, i, degtup);CHKERRQ(ierr);
148fbdc3dfeSToby Isaac     for (j = 0; j < Nk; j++) {
149fbdc3dfeSToby Isaac       PetscReal val = 1.;
150fbdc3dfeSToby Isaac 
151fbdc3dfeSToby Isaac       ierr = PetscDTIndexToGradedOrder(dim, j, ktup);CHKERRQ(ierr);
152fbdc3dfeSToby Isaac       for (l = 0; l < dim; l++) {
153fbdc3dfeSToby Isaac         val *= D[l][degtup[l]*(k+1) + ktup[l]];
154fbdc3dfeSToby Isaac       }
155fbdc3dfeSToby Isaac       lgndre_jet[j] += mul * val;
156fbdc3dfeSToby Isaac     }
157fbdc3dfeSToby Isaac   }
158fbdc3dfeSToby Isaac   for (i = 0; i < dim; i++) {
159fbdc3dfeSToby Isaac     ierr = PetscFree(D[i]);CHKERRQ(ierr);
160fbdc3dfeSToby Isaac   }
161fbdc3dfeSToby Isaac   ierr = PetscFree(D);CHKERRQ(ierr);
162fbdc3dfeSToby Isaac 
163fbdc3dfeSToby Isaac   for (i = 0; i < Nk; i++) {
164fbdc3dfeSToby Isaac     PetscReal diff = lgndre_jet[i] - pkd_jet[i];
165fbdc3dfeSToby Isaac     PetscReal scale = 1. + PetscAbsReal(lgndre_jet[i]) + PetscAbsReal(pkd_jet[i]);
166fbdc3dfeSToby Isaac     PetscReal tol = 10. * PETSC_SMALL * scale;
167fbdc3dfeSToby Isaac 
168*2c71b3e2SJacob 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);
169fbdc3dfeSToby Isaac   }
170fbdc3dfeSToby Isaac 
171fbdc3dfeSToby Isaac   ierr = PetscFree2(degtup,ktup);CHKERRQ(ierr);
172fbdc3dfeSToby Isaac   ierr = PetscFree(degrees);CHKERRQ(ierr);
173fbdc3dfeSToby Isaac   ierr = PetscFree3(pkd_jet_basis, lgndre_jet, pkd_jet);CHKERRQ(ierr);
174fbdc3dfeSToby Isaac   ierr = PetscFree(point);CHKERRQ(ierr);
175fbdc3dfeSToby Isaac   ierr = PetscFree2(lgndre_coeffs, pkd_coeffs);CHKERRQ(ierr);
176fbdc3dfeSToby Isaac   ierr = PetscFree(proj);CHKERRQ(ierr);
177fbdc3dfeSToby Isaac   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
178fbdc3dfeSToby Isaac   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
179fbdc3dfeSToby Isaac   PetscFunctionReturn(0);
180fbdc3dfeSToby Isaac }
181fbdc3dfeSToby Isaac 
182fbdc3dfeSToby Isaac int main(int argc, char **argv)
183fbdc3dfeSToby Isaac {
184fbdc3dfeSToby Isaac   PetscInt       dim, deg, k;
185fbdc3dfeSToby Isaac   PetscErrorCode ierr;
186fbdc3dfeSToby Isaac 
187fbdc3dfeSToby Isaac   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
188fbdc3dfeSToby Isaac   dim = 3;
189fbdc3dfeSToby Isaac   deg = 4;
190fbdc3dfeSToby Isaac   k = 3;
191fbdc3dfeSToby Isaac   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PetscDTPKDEval() tests","none");CHKERRQ(ierr);
192fbdc3dfeSToby Isaac   ierr = PetscOptionsInt("-dim", "Dimension of the simplex","ex9.c",dim,&dim,NULL);CHKERRQ(ierr);
193fbdc3dfeSToby Isaac   ierr = PetscOptionsInt("-degree", "The degree of the polynomial space","ex9.c",deg,&deg,NULL);CHKERRQ(ierr);
194fbdc3dfeSToby Isaac   ierr = PetscOptionsInt("-k", "The number of derivatives to use in the taylor test","ex9.c",k,&k,NULL);CHKERRQ(ierr);
1951e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
196fbdc3dfeSToby Isaac   ierr = testOrthogonality(dim, deg);CHKERRQ(ierr);
197fbdc3dfeSToby Isaac   ierr = testDerivativesLegendre(dim, deg, k);CHKERRQ(ierr);
198fbdc3dfeSToby Isaac   ierr = PetscFinalize();
199fbdc3dfeSToby Isaac   return ierr;
200fbdc3dfeSToby Isaac }
201fbdc3dfeSToby Isaac 
202fbdc3dfeSToby Isaac /*TEST
203fbdc3dfeSToby Isaac 
204fbdc3dfeSToby Isaac   test:
205fbdc3dfeSToby Isaac     args: -dim {{1 2 3 4}}
206fbdc3dfeSToby Isaac 
207fbdc3dfeSToby Isaac TEST*/
208