1 static char help[] = "Spot test DMStag->DMDA routines in 3d\n\n";
2 #include <petscdm.h>
3 #include <petscdmstag.h>
4 
main(int argc,char ** argv)5 int main(int argc,char **argv)
6 {
7   PetscErrorCode  ierr;
8   DM              dm;
9   Vec             vec;
10 
11   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
12   ierr = DMStagCreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,4,4,4,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,3,3,3,3,DMSTAG_STENCIL_STAR,1,NULL,NULL,NULL,&dm);CHKERRQ(ierr);
13   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
14   ierr = DMSetUp(dm);CHKERRQ(ierr);
15   ierr = DMStagSetUniformCoordinatesProduct(dm,0.0,10.0,0.0,10.0,0.0,10.0);CHKERRQ(ierr);
16 
17   ierr = DMCreateGlobalVector(dm,&vec);CHKERRQ(ierr);
18   ierr = VecSet(vec,1.234);CHKERRQ(ierr);
19 
20   /* All element values */
21   {
22     DM  da;
23     Vec vecda;
24     ierr = DMStagVecSplitToDMDA(dm,vec,DMSTAG_ELEMENT,-3,&da,&vecda);CHKERRQ(ierr);
25     ierr = DMDestroy(&da);CHKERRQ(ierr);
26     ierr = VecDestroy(&vecda);CHKERRQ(ierr);
27   }
28 
29   /* Pad element values */
30   {
31     DM  da;
32     Vec vecda;
33     ierr = DMStagVecSplitToDMDA(dm,vec,DMSTAG_ELEMENT,-5,&da,&vecda);CHKERRQ(ierr);
34     ierr = DMDestroy(&da);CHKERRQ(ierr);
35     ierr = VecDestroy(&vecda);CHKERRQ(ierr);
36   }
37 
38   /* 2 element values */
39   {
40     DM  da;
41     Vec vecda;
42     ierr = DMStagVecSplitToDMDA(dm,vec,DMSTAG_ELEMENT,-2,&da,&vecda);CHKERRQ(ierr);
43     ierr = DMDestroy(&da);CHKERRQ(ierr);
44     ierr = VecDestroy(&vecda);CHKERRQ(ierr);
45   }
46 
47   /* One corner value */
48   {
49     DM  da;
50     Vec vecda;
51     ierr = DMStagVecSplitToDMDA(dm,vec,DMSTAG_FRONT_DOWN_LEFT,2,&da,&vecda);CHKERRQ(ierr);
52     ierr = DMDestroy(&da);CHKERRQ(ierr);
53     ierr = VecDestroy(&vecda);CHKERRQ(ierr);
54   }
55 
56   /* One edge value */
57   {
58     DM  da;
59     Vec vecda;
60     ierr = DMStagVecSplitToDMDA(dm,vec,DMSTAG_BACK_RIGHT,1,&da,&vecda);CHKERRQ(ierr);
61     ierr = DMDestroy(&da);CHKERRQ(ierr);
62     ierr = VecDestroy(&vecda);CHKERRQ(ierr);
63   }
64 
65   /* One face value */
66   {
67     DM  da;
68     Vec vecda;
69     ierr = DMStagVecSplitToDMDA(dm,vec,DMSTAG_DOWN,0,&da,&vecda);CHKERRQ(ierr);
70     ierr = DMDestroy(&da);CHKERRQ(ierr);
71     ierr = VecDestroy(&vecda);CHKERRQ(ierr);
72   }
73 
74   ierr = VecDestroy(&vec);CHKERRQ(ierr);
75   ierr = DMDestroy(&dm);CHKERRQ(ierr);
76   ierr = PetscFinalize();
77   return ierr;
78 }
79 
80 /*TEST
81 
82    test:
83       suffix: 1
84       nsize: 12
85       args: -stag_ranks_x 2 -stag_ranks_y 3 -stag_ranks_z 2
86 
87 TEST*/
88