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; 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]; 25*08401ef6SPierre Jolivet PetscCheck(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 } 289566063dSJacob Faibussowitsch PetscCall(PetscFree(p)); 299566063dSJacob Faibussowitsch PetscCall(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; 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)); 61fbdc3dfeSToby Isaac for (i = 0; i < Np; i++) for (j = 0; j < npoints; j++) proj[i * npoints + j] *= weights[j]; 62fbdc3dfeSToby Isaac 639566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 649566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1., 1.)); 65fbdc3dfeSToby Isaac 66fbdc3dfeSToby Isaac /* create a random coefficient vector */ 679566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Np, &lgndre_coeffs, Np, &pkd_coeffs)); 68fbdc3dfeSToby Isaac for (i = 0; i < Np; i++) { 699566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &lgndre_coeffs[i])); 70fbdc3dfeSToby Isaac } 71fbdc3dfeSToby Isaac 729566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dim, °tup, dim, &ktup)); 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(deg + 1, °rees)); 74fbdc3dfeSToby Isaac for (i = 0; i < deg + 1; i++) degrees[i] = i; 75fbdc3dfeSToby Isaac 76fbdc3dfeSToby Isaac /* project the lgndre_coeffs to pkd_coeffs */ 779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pkd_coeffs, Np)); 789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &points1d)); 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &B)); 80fbdc3dfeSToby Isaac for (d = 0; d < dim; d++) { 819566063dSJacob Faibussowitsch PetscCall(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]; 849566063dSJacob Faibussowitsch PetscCall(PetscDTLegendreEval(npoints, points1d, deg + 1, degrees, B[d], NULL, NULL)); 85fbdc3dfeSToby Isaac } 869566063dSJacob Faibussowitsch PetscCall(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 949566063dSJacob Faibussowitsch PetscCall(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++) { 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(B[i])); 106fbdc3dfeSToby Isaac } 1079566063dSJacob Faibussowitsch PetscCall(PetscFree(B)); 108fbdc3dfeSToby Isaac 109fbdc3dfeSToby Isaac /* create a random point in the biunit simplex */ 1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &point)); 111fbdc3dfeSToby Isaac for (i = 0; i < dim; i++) { 1129566063dSJacob Faibussowitsch PetscCall(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 1239566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Nk*Np, &pkd_jet_basis, Nk, &lgndre_jet, Nk, &pkd_jet)); 124fbdc3dfeSToby Isaac /* evaluate the jet at the point with PKD polynomials */ 1259566063dSJacob Faibussowitsch PetscCall(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 */ 1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim, &D)); 136fbdc3dfeSToby Isaac for (i = 0; i < dim; i++) { 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((deg + 1) * (k+1), &(D[i]))); 1389566063dSJacob Faibussowitsch PetscCall(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 1459566063dSJacob Faibussowitsch PetscCall(PetscDTIndexToGradedOrder(dim, i, degtup)); 146fbdc3dfeSToby Isaac for (j = 0; j < Nk; j++) { 147fbdc3dfeSToby Isaac PetscReal val = 1.; 148fbdc3dfeSToby Isaac 1499566063dSJacob Faibussowitsch PetscCall(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++) { 1579566063dSJacob Faibussowitsch PetscCall(PetscFree(D[i])); 158fbdc3dfeSToby Isaac } 1599566063dSJacob Faibussowitsch PetscCall(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 166*08401ef6SPierre 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); 167fbdc3dfeSToby Isaac } 168fbdc3dfeSToby Isaac 1699566063dSJacob Faibussowitsch PetscCall(PetscFree2(degtup,ktup)); 1709566063dSJacob Faibussowitsch PetscCall(PetscFree(degrees)); 1719566063dSJacob Faibussowitsch PetscCall(PetscFree3(pkd_jet_basis, lgndre_jet, pkd_jet)); 1729566063dSJacob Faibussowitsch PetscCall(PetscFree(point)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFree2(lgndre_coeffs, pkd_coeffs)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFree(proj)); 1759566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 1769566063dSJacob Faibussowitsch PetscCall(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 1859566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 186fbdc3dfeSToby Isaac dim = 3; 187fbdc3dfeSToby Isaac deg = 4; 188fbdc3dfeSToby Isaac k = 3; 1899566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PetscDTPKDEval() tests","none");PetscCall(ierr); 1909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dim", "Dimension of the simplex","ex9.c",dim,&dim,NULL)); 1919566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-degree", "The degree of the polynomial space","ex9.c",deg,°,NULL)); 1929566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-k", "The number of derivatives to use in the taylor test","ex9.c",k,&k,NULL)); 1939566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 1949566063dSJacob Faibussowitsch PetscCall(testOrthogonality(dim, deg)); 1959566063dSJacob Faibussowitsch PetscCall(testDerivativesLegendre(dim, deg, k)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 197b122ec5aSJacob 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