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