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