xref: /petsc/src/dm/dt/tests/ex2.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Tests 1D cell-based discretization tools.\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   PetscInt       i,j,degrees[1000],ndegrees,nsrc_points,ntarget_points;
10c4762a1bSJed Brown   PetscReal      src_points[1000],target_points[1000],*R;
11c4762a1bSJed Brown   PetscBool      flg;
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
14c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Discretization tools test options",NULL);CHKERRQ(ierr);
15c4762a1bSJed Brown   {
16c4762a1bSJed Brown     ndegrees   = 1000;
17c4762a1bSJed Brown     degrees[0] = 1;
18c4762a1bSJed Brown     degrees[1] = 2;
19c4762a1bSJed Brown     degrees[2] = 3;
20*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsIntArray("-degrees","list of max degrees to evaluate","",degrees,&ndegrees,&flg));
21c4762a1bSJed Brown     if (!flg) ndegrees = 3;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     nsrc_points   = 1000;
24c4762a1bSJed Brown     src_points[0] = -1.;
25c4762a1bSJed Brown     src_points[1] = 0.;
26c4762a1bSJed Brown     src_points[2] = 1.;
27*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsRealArray("-src_points","list of points defining intervals on which to integrate","",src_points,&nsrc_points,&flg));
28c4762a1bSJed Brown     if (!flg) nsrc_points = 3;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown     ntarget_points   = 1000;
31c4762a1bSJed Brown     target_points[0] = -1.;
32c4762a1bSJed Brown     target_points[1] = 0.;
33c4762a1bSJed Brown     target_points[2] = 1.;
34*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsRealArray("-target_points","list of points defining intervals on which to integrate","",target_points,&ntarget_points,&flg));
35c4762a1bSJed Brown     if (!flg) ntarget_points = 3;
36c4762a1bSJed Brown   }
37c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
38c4762a1bSJed Brown 
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1((nsrc_points-1)*(ntarget_points-1),&R));
40c4762a1bSJed Brown   for (i=0; i<ndegrees; i++) {
41*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTReconstructPoly(degrees[i],nsrc_points-1,src_points,ntarget_points-1,target_points,R));
42c4762a1bSJed Brown     for (j=0; j<(ntarget_points-1)*(nsrc_points-1); j++) { /* Truncate to zero for nicer output */
43c4762a1bSJed Brown       if (PetscAbs(R[j]) < 10*PETSC_MACHINE_EPSILON) R[j] = 0;
44c4762a1bSJed Brown     }
45c4762a1bSJed Brown     for (j=0; j<ntarget_points-1; j++) {
46*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Degree %D target interval (%g,%g)\n",degrees[i],(double)target_points[j],(double)target_points[j+1]));
47*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRealView(nsrc_points-1,R+j*(nsrc_points-1),PETSC_VIEWER_STDOUT_WORLD));
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown   }
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(R));
51c4762a1bSJed Brown   ierr = PetscFinalize();
52c4762a1bSJed Brown   return ierr;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*TEST
56c4762a1bSJed Brown   test:
57c4762a1bSJed Brown     suffix: 1
58c4762a1bSJed Brown     args: -degrees 1,2,3 -target_points -0.3,0,.2 -src_points -1,-.3,0,.2,1
59c4762a1bSJed Brown TEST*/
60