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 DM da_c,da_f; 11c4762a1bSJed Brown Vec v_c,v_f; 12c4762a1bSJed Brown Mat Interp; 13c4762a1bSJed Brown PetscScalar one = 1.0; 14c4762a1bSJed Brown PetscBool pt; 15c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE,bz = DM_BOUNDARY_NONE; 16c4762a1bSJed Brown 17*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M1,NULL)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&s,NULL)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ratio",&ratio,NULL)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-periodic",(PetscBool*)&pt,NULL)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown if (pt) { 26c4762a1bSJed Brown if (dim > 0) bx = DM_BOUNDARY_PERIODIC; 27c4762a1bSJed Brown if (dim > 1) by = DM_BOUNDARY_PERIODIC; 28c4762a1bSJed Brown if (dim > 2) bz = DM_BOUNDARY_PERIODIC; 29c4762a1bSJed Brown } 30c4762a1bSJed Brown if (bx == DM_BOUNDARY_NONE) { 31c4762a1bSJed Brown M2 = ratio*(M1-1) + 1; 32c4762a1bSJed Brown } else { 33c4762a1bSJed Brown M2 = ratio*M1; 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Set up the array */ 37c4762a1bSJed Brown if (dim == 1) { 385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,bx,M1,dof,s,NULL,&da_c)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,bx,M2,dof,s,NULL,&da_f)); 40c4762a1bSJed Brown } else if (dim == 2) { 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M1,M1,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M2,M2,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_f)); 43c4762a1bSJed Brown } else if (dim == 3) { 445f80ce2aSJacob 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)); 455f80ce2aSJacob 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)); 46c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"dim must be 1,2, or 3"); 475f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da_c)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da_c)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da_f)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da_f)); 51c4762a1bSJed Brown 525f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da_c,&v_c)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da_f,&v_f)); 54c4762a1bSJed Brown 555f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(v_c,one)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(da_c,da_f,&Interp,NULL)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Interp,v_c,v_f)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v_f,PETSC_VIEWER_STDOUT_WORLD)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(Interp,v_f,v_c)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v_c,PETSC_VIEWER_STDOUT_WORLD)); 61c4762a1bSJed Brown 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Interp)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_c)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da_c)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v_f)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da_f)); 67*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 68*b122ec5aSJacob Faibussowitsch return 0; 69c4762a1bSJed Brown } 70