xref: /petsc/src/dm/dt/tests/ex15.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
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
testQuadrature(PetscInt dim,PetscInt degree,PetscDTSimplexQuadratureType type)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode testQuadrature(PetscInt dim, PetscInt degree, PetscDTSimplexQuadratureType type)
8d71ae5a4SJacob Faibussowitsch {
9d3c69ad0SToby Isaac   PetscInt         num_points;
10d3c69ad0SToby Isaac   const PetscReal *points;
11d3c69ad0SToby Isaac   const PetscReal *weights;
12d3c69ad0SToby Isaac   PetscInt         p_degree     = (degree + 1) / 2;
13d3c69ad0SToby Isaac   PetscInt         p_degree_min = degree - p_degree;
14d3c69ad0SToby Isaac   PetscInt         Nb, Nb_min;
15d3c69ad0SToby Isaac   PetscReal       *eval;
16d3c69ad0SToby Isaac   PetscQuadrature  quad;
17d3c69ad0SToby Isaac 
18d3c69ad0SToby Isaac   PetscFunctionBegin;
19d3c69ad0SToby Isaac   PetscCall(PetscDTSimplexQuadrature(dim, degree, type, &quad));
20d3c69ad0SToby Isaac   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &num_points, &points, &weights));
21d3c69ad0SToby Isaac   PetscCall(PetscDTBinomialInt(dim + p_degree, dim, &Nb));
22d3c69ad0SToby Isaac   PetscCall(PetscDTBinomialInt(dim + p_degree_min, dim, &Nb_min));
23d3c69ad0SToby Isaac   PetscCall(PetscMalloc1(num_points * Nb, &eval));
24d3c69ad0SToby Isaac   PetscCall(PetscDTPKDEvalJet(dim, num_points, points, p_degree, 0, eval));
25d3c69ad0SToby Isaac   for (PetscInt i = 0; i < Nb_min; i++) {
26d3c69ad0SToby Isaac     for (PetscInt j = i; j < Nb; j++) {
27d3c69ad0SToby Isaac       PetscReal I_exact = (i == j) ? 1.0 : 0.0;
28d3c69ad0SToby Isaac       PetscReal I_quad  = 0.0;
29d3c69ad0SToby Isaac       PetscReal err;
30d3c69ad0SToby Isaac 
31ad540459SPierre Jolivet       for (PetscInt q = 0; q < num_points; q++) I_quad += weights[q] * eval[i * num_points + q] * eval[j * num_points + q];
32d3c69ad0SToby Isaac       err = PetscAbsReal(I_exact - I_quad);
33300f1712SStefano Zampini       PetscCall(PetscInfo(quad, "Dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", method %s, error in <P_PKD(%" PetscInt_FMT "),P_PKD%" PetscInt_FMT "d)> = %g\n", dim, degree, PetscDTSimplexQuadratureTypes[type], i, j, (double)err));
34300f1712SStefano Zampini       PetscCheck(err < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", method %s, error in <P_PKD(%" PetscInt_FMT "),P_PKD(%" PetscInt_FMT ")> = %g", dim, degree, PetscDTSimplexQuadratureTypes[type], i, j, (double)err);
35d3c69ad0SToby Isaac     }
36d3c69ad0SToby Isaac   }
37d3c69ad0SToby Isaac   PetscCall(PetscFree(eval));
38d3c69ad0SToby Isaac   PetscCall(PetscQuadratureDestroy(&quad));
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40d3c69ad0SToby Isaac }
41d3c69ad0SToby Isaac 
main(int argc,char ** argv)42d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
43d71ae5a4SJacob Faibussowitsch {
44d3c69ad0SToby Isaac   const PetscInt dimdeg[] = {0, 20, 20, 20};
45d3c69ad0SToby Isaac 
46d3c69ad0SToby Isaac   PetscFunctionBeginUser;
47d3c69ad0SToby Isaac   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
48d3c69ad0SToby Isaac   for (PetscInt dim = 0; dim <= 3; dim++) {
49d3c69ad0SToby Isaac     for (PetscInt deg = 0; deg <= dimdeg[dim]; deg++) {
50d3c69ad0SToby Isaac       const PetscDTSimplexQuadratureType types[] = {PETSCDTSIMPLEXQUAD_DEFAULT, PETSCDTSIMPLEXQUAD_CONIC, PETSCDTSIMPLEXQUAD_MINSYM};
51d3c69ad0SToby Isaac 
5248a46eb9SPierre Jolivet       for (PetscInt t = 0; t < 3; t++) PetscCall(testQuadrature(dim, deg, types[t]));
53d3c69ad0SToby Isaac     }
54d3c69ad0SToby Isaac   }
55d3c69ad0SToby Isaac   PetscCall(PetscFinalize());
56d3c69ad0SToby Isaac   return 0;
57d3c69ad0SToby Isaac }
58d3c69ad0SToby Isaac 
59d3c69ad0SToby Isaac /*TEST
60d3c69ad0SToby Isaac 
61d3c69ad0SToby Isaac   test:
62d3c69ad0SToby Isaac     suffix: 0
63*3886731fSPierre Jolivet     output_file: output/empty.out
64d3c69ad0SToby Isaac 
65d3c69ad0SToby Isaac TEST*/
66