1d3c69ad0SToby Isaac const char help[] = "Test PetscDTSimplexQuadrature()"; 2d3c69ad0SToby Isaac 3d3c69ad0SToby Isaac #include <petscdt.h> 4d3c69ad0SToby Isaac 5d3c69ad0SToby Isaac // if we trust the PKD polynomials (tested in ex9), then we can see how well the quadrature integrates 6d3c69ad0SToby Isaac // the mass matrix, which should be the identity 79371c9d4SSatish Balay static PetscErrorCode testQuadrature(PetscInt dim, PetscInt degree, PetscDTSimplexQuadratureType type) { 8d3c69ad0SToby Isaac PetscInt num_points; 9d3c69ad0SToby Isaac const PetscReal *points; 10d3c69ad0SToby Isaac const PetscReal *weights; 11d3c69ad0SToby Isaac PetscInt p_degree = (degree + 1) / 2; 12d3c69ad0SToby Isaac PetscInt p_degree_min = degree - p_degree; 13d3c69ad0SToby Isaac PetscInt Nb, Nb_min; 14d3c69ad0SToby Isaac PetscReal *eval; 15d3c69ad0SToby Isaac PetscQuadrature quad; 16d3c69ad0SToby Isaac 17d3c69ad0SToby Isaac PetscFunctionBegin; 18d3c69ad0SToby Isaac PetscCall(PetscDTSimplexQuadrature(dim, degree, type, &quad)); 19d3c69ad0SToby Isaac PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &num_points, &points, &weights)); 20d3c69ad0SToby Isaac PetscCall(PetscDTBinomialInt(dim + p_degree, dim, &Nb)); 21d3c69ad0SToby Isaac PetscCall(PetscDTBinomialInt(dim + p_degree_min, dim, &Nb_min)); 22d3c69ad0SToby Isaac PetscCall(PetscMalloc1(num_points * Nb, &eval)); 23d3c69ad0SToby Isaac PetscCall(PetscDTPKDEvalJet(dim, num_points, points, p_degree, 0, eval)); 24d3c69ad0SToby Isaac for (PetscInt i = 0; i < Nb_min; i++) { 25d3c69ad0SToby Isaac for (PetscInt j = i; j < Nb; j++) { 26d3c69ad0SToby Isaac PetscReal I_exact = (i == j) ? 1.0 : 0.0; 27d3c69ad0SToby Isaac PetscReal I_quad = 0.0; 28d3c69ad0SToby Isaac PetscReal err; 29d3c69ad0SToby Isaac 309371c9d4SSatish Balay for (PetscInt q = 0; q < num_points; q++) { I_quad += weights[q] * eval[i * num_points + q] * eval[j * num_points + q]; } 31d3c69ad0SToby Isaac err = PetscAbsReal(I_exact - I_quad); 32d3c69ad0SToby Isaac PetscCall(PetscInfo(quad, "Dimension %d, degree %d, method %s, error in <P_PKD(%d),P_PKD(%d)> = %g\n", (int)dim, (int)degree, PetscDTSimplexQuadratureTypes[type], (int)i, (int)j, (double)err)); 33d3c69ad0SToby Isaac PetscCheck(err < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dimension %d, degree %d, method %s, error in <P_PKD(%d),P_PKD(%d)> = %g", (int)dim, (int)degree, PetscDTSimplexQuadratureTypes[type], (int)i, (int)j, (double)err); 34d3c69ad0SToby Isaac } 35d3c69ad0SToby Isaac } 36d3c69ad0SToby Isaac PetscCall(PetscFree(eval)); 37d3c69ad0SToby Isaac PetscCall(PetscQuadratureDestroy(&quad)); 38d3c69ad0SToby Isaac PetscFunctionReturn(0); 39d3c69ad0SToby Isaac } 40d3c69ad0SToby Isaac 419371c9d4SSatish Balay int main(int argc, char **argv) { 42d3c69ad0SToby Isaac const PetscInt dimdeg[] = {0, 20, 20, 20}; 43d3c69ad0SToby Isaac 44d3c69ad0SToby Isaac PetscFunctionBeginUser; 45d3c69ad0SToby Isaac PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 46d3c69ad0SToby Isaac for (PetscInt dim = 0; dim <= 3; dim++) { 47d3c69ad0SToby Isaac for (PetscInt deg = 0; deg <= dimdeg[dim]; deg++) { 48d3c69ad0SToby Isaac const PetscDTSimplexQuadratureType types[] = {PETSCDTSIMPLEXQUAD_DEFAULT, PETSCDTSIMPLEXQUAD_CONIC, PETSCDTSIMPLEXQUAD_MINSYM}; 49d3c69ad0SToby Isaac 50*48a46eb9SPierre Jolivet for (PetscInt t = 0; t < 3; t++) PetscCall(testQuadrature(dim, deg, types[t])); 51d3c69ad0SToby Isaac } 52d3c69ad0SToby Isaac } 53d3c69ad0SToby Isaac PetscCall(PetscFinalize()); 54d3c69ad0SToby Isaac return 0; 55d3c69ad0SToby Isaac } 56d3c69ad0SToby Isaac 57d3c69ad0SToby Isaac /*TEST 58d3c69ad0SToby Isaac 59d3c69ad0SToby Isaac test: 60d3c69ad0SToby Isaac suffix: 0 61d3c69ad0SToby Isaac 62d3c69ad0SToby Isaac TEST*/ 63