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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 18*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 19*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M1,NULL)); 20*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&s,NULL)); 21*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ratio",&ratio,NULL)); 22*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL)); 23*9566063dSJacob Faibussowitsch PetscCall(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) { 38*9566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD,bx,M1,dof,s,NULL,&da_c)); 39*9566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD,bx,M2,dof,s,NULL,&da_f)); 40c4762a1bSJed Brown } else if (dim == 2) { 41*9566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M1,M1,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c)); 42*9566063dSJacob Faibussowitsch PetscCall(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) { 44*9566063dSJacob Faibussowitsch PetscCall(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)); 45*9566063dSJacob Faibussowitsch PetscCall(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"); 47*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da_c)); 48*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(da_c)); 49*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da_f)); 50*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(da_f)); 51c4762a1bSJed Brown 52*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da_c,&v_c)); 53*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da_f,&v_f)); 54c4762a1bSJed Brown 55*9566063dSJacob Faibussowitsch PetscCall(VecSet(v_c,one)); 56*9566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(da_c,da_f,&Interp,NULL)); 57*9566063dSJacob Faibussowitsch PetscCall(MatMult(Interp,v_c,v_f)); 58*9566063dSJacob Faibussowitsch PetscCall(VecView(v_f,PETSC_VIEWER_STDOUT_WORLD)); 59*9566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Interp,v_f,v_c)); 60*9566063dSJacob Faibussowitsch PetscCall(VecView(v_c,PETSC_VIEWER_STDOUT_WORLD)); 61c4762a1bSJed Brown 62*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp)); 63*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_c)); 64*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da_c)); 65*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v_f)); 66*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da_f)); 67*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 68b122ec5aSJacob Faibussowitsch return 0; 69c4762a1bSJed Brown } 70