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