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 DM da_c,da_f; 11 Vec v_c,v_f; 12 Mat Interp; 13 PetscScalar one = 1.0; 14 PetscBool pt; 15 DMBoundaryType bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE,bz = DM_BOUNDARY_NONE; 16 17 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 18 PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 19 PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M1,NULL)); 20 PetscCall(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&s,NULL)); 21 PetscCall(PetscOptionsGetInt(NULL,NULL,"-ratio",&ratio,NULL)); 22 PetscCall(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL)); 23 PetscCall(PetscOptionsGetBool(NULL,NULL,"-periodic",(PetscBool*)&pt,NULL)); 24 25 if (pt) { 26 if (dim > 0) bx = DM_BOUNDARY_PERIODIC; 27 if (dim > 1) by = DM_BOUNDARY_PERIODIC; 28 if (dim > 2) bz = DM_BOUNDARY_PERIODIC; 29 } 30 if (bx == DM_BOUNDARY_NONE) { 31 M2 = ratio*(M1-1) + 1; 32 } else { 33 M2 = ratio*M1; 34 } 35 36 /* Set up the array */ 37 if (dim == 1) { 38 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,bx,M1,dof,s,NULL,&da_c)); 39 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,bx,M2,dof,s,NULL,&da_f)); 40 } else if (dim == 2) { 41 PetscCall(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M1,M1,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c)); 42 PetscCall(DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M2,M2,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_f)); 43 } else if (dim == 3) { 44 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 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)); 46 } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"dim must be 1,2, or 3"); 47 PetscCall(DMSetFromOptions(da_c)); 48 PetscCall(DMSetUp(da_c)); 49 PetscCall(DMSetFromOptions(da_f)); 50 PetscCall(DMSetUp(da_f)); 51 52 PetscCall(DMCreateGlobalVector(da_c,&v_c)); 53 PetscCall(DMCreateGlobalVector(da_f,&v_f)); 54 55 PetscCall(VecSet(v_c,one)); 56 PetscCall(DMCreateInterpolation(da_c,da_f,&Interp,NULL)); 57 PetscCall(MatMult(Interp,v_c,v_f)); 58 PetscCall(VecView(v_f,PETSC_VIEWER_STDOUT_WORLD)); 59 PetscCall(MatMultTranspose(Interp,v_f,v_c)); 60 PetscCall(VecView(v_c,PETSC_VIEWER_STDOUT_WORLD)); 61 62 PetscCall(MatDestroy(&Interp)); 63 PetscCall(VecDestroy(&v_c)); 64 PetscCall(DMDestroy(&da_c)); 65 PetscCall(VecDestroy(&v_f)); 66 PetscCall(DMDestroy(&da_f)); 67 PetscCall(PetscFinalize()); 68 return 0; 69 } 70