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