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