xref: /petsc/src/dm/tests/ex15.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests DMDA interpolation.\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 
7 int main(int argc,char **argv)
8 {
9   PetscInt         M1 = 3,M2,dof = 1,s = 1,ratio = 2,dim = 1;
10   PetscErrorCode   ierr;
11   DM               da_c,da_f;
12   Vec              v_c,v_f;
13   Mat              Interp;
14   PetscScalar      one = 1.0;
15   PetscBool        pt;
16   DMBoundaryType   bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE,bz = DM_BOUNDARY_NONE;
17 
18   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
20   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M1,NULL));
21   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&s,NULL));
22   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ratio",&ratio,NULL));
23   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL));
24   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-periodic",(PetscBool*)&pt,NULL));
25 
26   if (pt) {
27     if (dim > 0) bx = DM_BOUNDARY_PERIODIC;
28     if (dim > 1) by = DM_BOUNDARY_PERIODIC;
29     if (dim > 2) bz = DM_BOUNDARY_PERIODIC;
30   }
31   if (bx == DM_BOUNDARY_NONE) {
32     M2 = ratio*(M1-1) + 1;
33   } else {
34     M2 = ratio*M1;
35   }
36 
37   /* Set up the array */
38   if (dim == 1) {
39     CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,bx,M1,dof,s,NULL,&da_c));
40     CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,bx,M2,dof,s,NULL,&da_f));
41   } else if (dim == 2) {
42     CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M1,M1,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c));
43     CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M2,M2,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_f));
44   } else if (dim == 3) {
45     CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,DMDA_STENCIL_BOX,M1,M1,M1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,NULL,&da_c));
46     CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,DMDA_STENCIL_BOX,M2,M2,M2,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,NULL,&da_f));
47   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"dim must be 1,2, or 3");
48   CHKERRQ(DMSetFromOptions(da_c));
49   CHKERRQ(DMSetUp(da_c));
50   CHKERRQ(DMSetFromOptions(da_f));
51   CHKERRQ(DMSetUp(da_f));
52 
53   CHKERRQ(DMCreateGlobalVector(da_c,&v_c));
54   CHKERRQ(DMCreateGlobalVector(da_f,&v_f));
55 
56   CHKERRQ(VecSet(v_c,one));
57   CHKERRQ(DMCreateInterpolation(da_c,da_f,&Interp,NULL));
58   CHKERRQ(MatMult(Interp,v_c,v_f));
59   CHKERRQ(VecView(v_f,PETSC_VIEWER_STDOUT_WORLD));
60   CHKERRQ(MatMultTranspose(Interp,v_f,v_c));
61   CHKERRQ(VecView(v_c,PETSC_VIEWER_STDOUT_WORLD));
62 
63   CHKERRQ(MatDestroy(&Interp));
64   CHKERRQ(VecDestroy(&v_c));
65   CHKERRQ(DMDestroy(&da_c));
66   CHKERRQ(VecDestroy(&v_f));
67   CHKERRQ(DMDestroy(&da_f));
68   ierr = PetscFinalize();
69   return ierr;
70 }
71