xref: /petsc/src/dm/dt/tests/ex6.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Tests 1D Gauss-Lobatto-Legendre discretization on [-1, 1].\n\n";
2 
3 #include <petscdt.h>
4 #include <petscviewer.h>
5 
6 int main(int argc, char **argv) {
7   PetscInt   n = 3, i;
8   PetscReal *la_nodes, *la_weights, *n_nodes, *n_weights;
9 
10   PetscFunctionBeginUser;
11   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
12   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
13 
14   PetscCall(PetscMalloc1(n, &la_nodes));
15   PetscCall(PetscMalloc1(n, &la_weights));
16   PetscCall(PetscMalloc1(n, &n_nodes));
17   PetscCall(PetscMalloc1(n, &n_weights));
18   PetscCall(PetscDTGaussLobattoLegendreQuadrature(n, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, la_nodes, la_weights));
19 
20   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Gauss-Lobatto-Legendre nodes and weights computed via linear algebra: \n"));
21   PetscCall(PetscRealView(n, la_nodes, PETSC_VIEWER_STDOUT_SELF));
22   PetscCall(PetscRealView(n, la_weights, PETSC_VIEWER_STDOUT_SELF));
23   PetscCall(PetscDTGaussLobattoLegendreQuadrature(n, PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON, n_nodes, n_weights));
24   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Gauss-Lobatto-Legendre nodes and weights computed via Newton: \n"));
25   PetscCall(PetscRealView(n, n_nodes, PETSC_VIEWER_STDOUT_SELF));
26   PetscCall(PetscRealView(n, n_weights, PETSC_VIEWER_STDOUT_SELF));
27 
28   for (i = 0; i < n; i++) {
29     la_nodes[i] -= n_nodes[i];
30     la_weights[i] -= n_weights[i];
31   }
32   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Difference: \n"));
33   PetscCall(PetscRealView(n, la_nodes, PETSC_VIEWER_STDOUT_SELF));
34   PetscCall(PetscRealView(n, la_weights, PETSC_VIEWER_STDOUT_SELF));
35 
36   PetscCall(PetscFree(la_nodes));
37   PetscCall(PetscFree(la_weights));
38   PetscCall(PetscFree(n_nodes));
39   PetscCall(PetscFree(n_weights));
40 
41   PetscCall(PetscFinalize());
42   return 0;
43 }
44 
45 /*TEST
46 
47    test:
48 
49 TEST*/
50