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