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