xref: /petsc/src/dm/dt/tests/ex6.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Tests 1D Gauss-Lobatto-Legendre discretization on [-1, 1].\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdt.h>
4c4762a1bSJed Brown #include <petscviewer.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   PetscInt       n = 3,i;
10c4762a1bSJed Brown   PetscReal      *la_nodes,*la_weights,*n_nodes,*n_weights;
11c4762a1bSJed Brown 
12*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
14c4762a1bSJed Brown 
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&la_nodes));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&la_weights));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&n_nodes));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&n_weights));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTGaussLobattoLegendreQuadrature(n,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,la_nodes,la_weights));
20c4762a1bSJed Brown 
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"Gauss-Lobatto-Legendre nodes and weights computed via linear algebra: \n"));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(n,la_nodes,PETSC_VIEWER_STDOUT_SELF));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(n,la_weights,PETSC_VIEWER_STDOUT_SELF));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTGaussLobattoLegendreQuadrature(n,PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON,n_nodes,n_weights));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"Gauss-Lobatto-Legendre nodes and weights computed via Newton: \n"));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(n,n_nodes,PETSC_VIEWER_STDOUT_SELF));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(n,n_weights,PETSC_VIEWER_STDOUT_SELF));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   for (i=0; i<n; i++) {
30c4762a1bSJed Brown     la_nodes[i]   -= n_nodes[i];
31c4762a1bSJed Brown     la_weights[i] -= n_weights[i];
32c4762a1bSJed Brown   }
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"Difference: \n"));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(n,la_nodes,PETSC_VIEWER_STDOUT_SELF));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(n,la_weights,PETSC_VIEWER_STDOUT_SELF));
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(la_nodes));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(la_weights));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(n_nodes));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(n_weights));
41c4762a1bSJed Brown 
42*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
43*b122ec5aSJacob Faibussowitsch   return 0;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*TEST
47c4762a1bSJed Brown 
48c4762a1bSJed Brown    test:
49c4762a1bSJed Brown 
50c4762a1bSJed Brown TEST*/
51