1fbdc3dfeSToby Isaac const char help[] = "Tests PetscDTPKDEvalJet()"; 2fbdc3dfeSToby Isaac 3fbdc3dfeSToby Isaac #include <petscdt.h> 4fbdc3dfeSToby Isaac #include <petscblaslapack.h> 5fbdc3dfeSToby Isaac 69371c9d4SSatish Balay static PetscErrorCode testOrthogonality(PetscInt dim, PetscInt deg) { 7fbdc3dfeSToby Isaac PetscQuadrature q; 8fbdc3dfeSToby Isaac const PetscReal *points, *weights; 9fbdc3dfeSToby Isaac PetscInt Npoly, npoints, i, j, k; 10fbdc3dfeSToby Isaac PetscReal *p; 11fbdc3dfeSToby Isaac 12fbdc3dfeSToby Isaac PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q)); 149566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights)); 159566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + deg, dim, &Npoly)); 169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Npoly * npoints, &p)); 179566063dSJacob Faibussowitsch PetscCall(PetscDTPKDEvalJet(dim, npoints, points, deg, 0, p)); 18fbdc3dfeSToby Isaac for (i = 0; i < Npoly; i++) { 19fbdc3dfeSToby Isaac for (j = i; j < Npoly; j++) { 20fbdc3dfeSToby Isaac PetscReal integral = 0.; 21fbdc3dfeSToby Isaac PetscReal exact = (i == j) ? 1. : 0.; 22fbdc3dfeSToby Isaac 23fbdc3dfeSToby Isaac for (k = 0; k < npoints; k++) integral += weights[k] * p[i * npoints + k] * p[j * npoints + k]; 2463a3b9bcSJacob 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); 25fbdc3dfeSToby Isaac } 26fbdc3dfeSToby Isaac } 279566063dSJacob Faibussowitsch PetscCall(PetscFree(p)); 289566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 29fbdc3dfeSToby Isaac PetscFunctionReturn(0); 30fbdc3dfeSToby Isaac } 31fbdc3dfeSToby Isaac 329371c9d4SSatish Balay static PetscErrorCode testDerivativesLegendre(PetscInt dim, PetscInt deg, PetscInt k) { 33fbdc3dfeSToby Isaac PetscInt Np, Nk, i, j, l, d, npoints; 34fbdc3dfeSToby Isaac PetscRandom rand; 35fbdc3dfeSToby Isaac PetscReal *point; 36fbdc3dfeSToby Isaac PetscReal *lgndre_coeffs; 37fbdc3dfeSToby Isaac PetscReal *pkd_coeffs; 38fbdc3dfeSToby Isaac PetscReal *proj; 39fbdc3dfeSToby Isaac PetscReal **B; 40fbdc3dfeSToby Isaac PetscQuadrature q; 41fbdc3dfeSToby Isaac PetscReal *points1d; 42fbdc3dfeSToby Isaac PetscInt *degrees; 43fbdc3dfeSToby Isaac PetscInt *degtup, *ktup; 44fbdc3dfeSToby Isaac const PetscReal *points; 45fbdc3dfeSToby Isaac const PetscReal *weights; 46fbdc3dfeSToby Isaac PetscReal *lgndre_jet; 47fbdc3dfeSToby Isaac PetscReal **D; 48fbdc3dfeSToby Isaac PetscReal *pkd_jet, *pkd_jet_basis; 49fbdc3dfeSToby Isaac 50fbdc3dfeSToby Isaac PetscFunctionBegin; 519566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + deg, dim, &Np)); 529566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + k, dim, &Nk)); 53fbdc3dfeSToby Isaac 54fbdc3dfeSToby Isaac /* create the projector (because it is an orthonormal basis, the projector is the moment integrals) */ 559566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1., 1., &q)); 569566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights)); 579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints * Np, &proj)); 589566063dSJacob Faibussowitsch PetscCall(PetscDTPKDEvalJet(dim, npoints, points, deg, 0, proj)); 599371c9d4SSatish Balay for (i = 0; i < Np; i++) 609371c9d4SSatish Balay for (j = 0; j < npoints; j++) proj[i * npoints + j] *= weights[j]; 61fbdc3dfeSToby Isaac 629566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 639566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1., 1.)); 64fbdc3dfeSToby Isaac 65fbdc3dfeSToby Isaac /* create a random coefficient vector */ 669566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Np, &lgndre_coeffs, Np, &pkd_coeffs)); 67*48a46eb9SPierre Jolivet for (i = 0; i < Np; i++) PetscCall(PetscRandomGetValueReal(rand, &lgndre_coeffs[i])); 68fbdc3dfeSToby Isaac 699566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dim, °tup, dim, &ktup)); 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(deg + 1, °rees)); 71fbdc3dfeSToby Isaac for (i = 0; i < deg + 1; i++) degrees[i] = i; 72fbdc3dfeSToby Isaac 73fbdc3dfeSToby Isaac /* project the lgndre_coeffs to pkd_coeffs */ 749566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pkd_coeffs, Np)); 759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &points1d)); 769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &B)); 77fbdc3dfeSToby Isaac for (d = 0; d < dim; d++) { 789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((deg + 1) * npoints, &(B[d]))); 79fbdc3dfeSToby Isaac /* get this coordinate */ 80fbdc3dfeSToby Isaac for (i = 0; i < npoints; i++) points1d[i] = points[i * dim + d]; 819566063dSJacob Faibussowitsch PetscCall(PetscDTLegendreEval(npoints, points1d, deg + 1, degrees, B[d], NULL, NULL)); 82fbdc3dfeSToby Isaac } 839566063dSJacob Faibussowitsch PetscCall(PetscFree(points1d)); 84fbdc3dfeSToby Isaac for (i = 0; i < npoints; i++) { 85fbdc3dfeSToby Isaac PetscReal val = 0.; 86fbdc3dfeSToby Isaac 87fbdc3dfeSToby Isaac for (j = 0; j < Np; j++) { 88fbdc3dfeSToby Isaac PetscReal mul = lgndre_coeffs[j]; 89fbdc3dfeSToby Isaac PetscReal valj = 1.; 90fbdc3dfeSToby Isaac 919566063dSJacob Faibussowitsch PetscCall(PetscDTIndexToGradedOrder(dim, j, degtup)); 929371c9d4SSatish Balay for (l = 0; l < dim; l++) { valj *= B[l][i * (deg + 1) + degtup[l]]; } 93fbdc3dfeSToby Isaac val += mul * valj; 94fbdc3dfeSToby Isaac } 959371c9d4SSatish Balay for (j = 0; j < Np; j++) { pkd_coeffs[j] += proj[j * npoints + i] * val; } 96fbdc3dfeSToby Isaac } 97*48a46eb9SPierre Jolivet for (i = 0; i < dim; i++) PetscCall(PetscFree(B[i])); 989566063dSJacob Faibussowitsch PetscCall(PetscFree(B)); 99fbdc3dfeSToby Isaac 100fbdc3dfeSToby Isaac /* create a random point in the biunit simplex */ 1019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &point)); 102*48a46eb9SPierre Jolivet for (i = 0; i < dim; i++) PetscCall(PetscRandomGetValueReal(rand, &point[i])); 103fbdc3dfeSToby Isaac for (i = dim - 1; i > 0; i--) { 104fbdc3dfeSToby Isaac PetscReal val = point[i]; 105fbdc3dfeSToby Isaac PetscInt j; 106fbdc3dfeSToby Isaac 1079371c9d4SSatish Balay for (j = 0; j < i; j++) { point[j] = (point[j] + 1.) * (1. - val) * 0.5 - 1.; } 108fbdc3dfeSToby Isaac } 109fbdc3dfeSToby Isaac 1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Nk * Np, &pkd_jet_basis, Nk, &lgndre_jet, Nk, &pkd_jet)); 111fbdc3dfeSToby Isaac /* evaluate the jet at the point with PKD polynomials */ 1129566063dSJacob Faibussowitsch PetscCall(PetscDTPKDEvalJet(dim, 1, point, deg, k, pkd_jet_basis)); 113fbdc3dfeSToby Isaac for (i = 0; i < Nk; i++) { 114fbdc3dfeSToby Isaac PetscReal val = 0.; 1159371c9d4SSatish Balay for (j = 0; j < Np; j++) { val += pkd_coeffs[j] * pkd_jet_basis[j * Nk + i]; } 116fbdc3dfeSToby Isaac pkd_jet[i] = val; 117fbdc3dfeSToby Isaac } 118fbdc3dfeSToby Isaac 119fbdc3dfeSToby Isaac /* evaluate the 1D jets of the Legendre polynomials */ 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &D)); 121fbdc3dfeSToby Isaac for (i = 0; i < dim; i++) { 1229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((deg + 1) * (k + 1), &(D[i]))); 1239566063dSJacob Faibussowitsch PetscCall(PetscDTJacobiEvalJet(0., 0., 1, &(point[i]), deg, k, D[i])); 124fbdc3dfeSToby Isaac } 125fbdc3dfeSToby Isaac /* compile the 1D Legendre jets into the tensor Legendre jet */ 126fbdc3dfeSToby Isaac for (j = 0; j < Nk; j++) lgndre_jet[j] = 0.; 127fbdc3dfeSToby Isaac for (i = 0; i < Np; i++) { 128fbdc3dfeSToby Isaac PetscReal mul = lgndre_coeffs[i]; 129fbdc3dfeSToby Isaac 1309566063dSJacob Faibussowitsch PetscCall(PetscDTIndexToGradedOrder(dim, i, degtup)); 131fbdc3dfeSToby Isaac for (j = 0; j < Nk; j++) { 132fbdc3dfeSToby Isaac PetscReal val = 1.; 133fbdc3dfeSToby Isaac 1349566063dSJacob Faibussowitsch PetscCall(PetscDTIndexToGradedOrder(dim, j, ktup)); 1359371c9d4SSatish Balay for (l = 0; l < dim; l++) { val *= D[l][degtup[l] * (k + 1) + ktup[l]]; } 136fbdc3dfeSToby Isaac lgndre_jet[j] += mul * val; 137fbdc3dfeSToby Isaac } 138fbdc3dfeSToby Isaac } 139*48a46eb9SPierre Jolivet for (i = 0; i < dim; i++) PetscCall(PetscFree(D[i])); 1409566063dSJacob Faibussowitsch PetscCall(PetscFree(D)); 141fbdc3dfeSToby Isaac 142fbdc3dfeSToby Isaac for (i = 0; i < Nk; i++) { 143fbdc3dfeSToby Isaac PetscReal diff = lgndre_jet[i] - pkd_jet[i]; 144fbdc3dfeSToby Isaac PetscReal scale = 1. + PetscAbsReal(lgndre_jet[i]) + PetscAbsReal(pkd_jet[i]); 145fbdc3dfeSToby Isaac PetscReal tol = 10. * PETSC_SMALL * scale; 146fbdc3dfeSToby Isaac 14708401ef6SPierre 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); 148fbdc3dfeSToby Isaac } 149fbdc3dfeSToby Isaac 1509566063dSJacob Faibussowitsch PetscCall(PetscFree2(degtup, ktup)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFree(degrees)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree3(pkd_jet_basis, lgndre_jet, pkd_jet)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree(point)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFree2(lgndre_coeffs, pkd_coeffs)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFree(proj)); 1569566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 1579566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 158fbdc3dfeSToby Isaac PetscFunctionReturn(0); 159fbdc3dfeSToby Isaac } 160fbdc3dfeSToby Isaac 1619371c9d4SSatish Balay int main(int argc, char **argv) { 162fbdc3dfeSToby Isaac PetscInt dim, deg, k; 163fbdc3dfeSToby Isaac 164327415f7SBarry Smith PetscFunctionBeginUser; 1659566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 166fbdc3dfeSToby Isaac dim = 3; 167fbdc3dfeSToby Isaac deg = 4; 168fbdc3dfeSToby Isaac k = 3; 169d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PetscDTPKDEval() tests", "none"); 1709566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dim", "Dimension of the simplex", "ex9.c", dim, &dim, NULL)); 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-degree", "The degree of the polynomial space", "ex9.c", deg, °, NULL)); 1729566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-k", "The number of derivatives to use in the taylor test", "ex9.c", k, &k, NULL)); 173d0609cedSBarry Smith PetscOptionsEnd(); 1749566063dSJacob Faibussowitsch PetscCall(testOrthogonality(dim, deg)); 1759566063dSJacob Faibussowitsch PetscCall(testDerivativesLegendre(dim, deg, k)); 1769566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 177b122ec5aSJacob Faibussowitsch return 0; 178fbdc3dfeSToby Isaac } 179fbdc3dfeSToby Isaac 180fbdc3dfeSToby Isaac /*TEST 181fbdc3dfeSToby Isaac 182fbdc3dfeSToby Isaac test: 183fbdc3dfeSToby Isaac args: -dim {{1 2 3 4}} 184fbdc3dfeSToby Isaac 185fbdc3dfeSToby Isaac TEST*/ 186