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