Actual source code: ex15.c

petsc-3.4.2 2013-07-02
  2: static char help[] = "Tests DMDA interpolation.\n\n";

  4: #include <petscdmda.h>

  8: int main(int argc,char **argv)
  9: {
 10:   PetscInt         M1 = 3,M2,dof = 1,s = 1,ratio = 2,dim = 1;
 11:   PetscErrorCode   ierr;
 12:   DM               da_c,da_f;
 13:   Vec              v_c,v_f;
 14:   Mat              Interp;
 15:   PetscScalar      one = 1.0;
 16:   PetscBool        pt;
 17:   DMDABoundaryType bx = DMDA_BOUNDARY_NONE,by = DMDA_BOUNDARY_NONE,bz = DMDA_BOUNDARY_NONE;

 19:   PetscInitialize(&argc,&argv,(char*)0,help);

 21:   PetscOptionsGetInt(NULL,"-dim",&dim,NULL);
 22:   PetscOptionsGetInt(NULL,"-M",&M1,NULL);
 23:   PetscOptionsGetInt(NULL,"-stencil_width",&s,NULL);
 24:   PetscOptionsGetInt(NULL,"-ratio",&ratio,NULL);
 25:   PetscOptionsGetInt(NULL,"-dof",&dof,NULL);
 26:   PetscOptionsGetBool(NULL,"-periodic",(PetscBool*)&pt,NULL);

 28:   if (pt) {
 29:     if (dim > 0) bx = DMDA_BOUNDARY_PERIODIC;
 30:     if (dim > 1) by = DMDA_BOUNDARY_PERIODIC;
 31:     if (dim > 2) bz = DMDA_BOUNDARY_PERIODIC;
 32:   }
 33:   if (bx == DMDA_BOUNDARY_NONE) {
 34:     M2 = ratio*(M1-1) + 1;
 35:   } else {
 36:     M2 = ratio*M1;
 37:   }

 39:   /* Set up the array */
 40:   if (dim == 1) {
 41:     DMDACreate1d(PETSC_COMM_WORLD,bx,M1,dof,s,NULL,&da_c);
 42:     DMDACreate1d(PETSC_COMM_WORLD,bx,M2,dof,s,NULL,&da_f);
 43:   } else if (dim == 2) {
 44:     DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M1,M1,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c);
 45:     DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M2,M2,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_f);
 46:   } else if (dim == 3) {
 47:     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);
 48:     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);
 49:   }

 51:   DMCreateGlobalVector(da_c,&v_c);
 52:   DMCreateGlobalVector(da_f,&v_f);

 54:   VecSet(v_c,one);
 55:   DMCreateInterpolation(da_c,da_f,&Interp,NULL);
 56:   MatMult(Interp,v_c,v_f);
 57:   VecView(v_f,PETSC_VIEWER_STDOUT_WORLD);
 58:   MatMultTranspose(Interp,v_f,v_c);
 59:   VecView(v_c,PETSC_VIEWER_STDOUT_WORLD);

 61:   MatDestroy(&Interp);
 62:   VecDestroy(&v_c);
 63:   DMDestroy(&da_c);
 64:   VecDestroy(&v_f);
 65:   DMDestroy(&da_f);
 66:   PetscFinalize();
 67:   return 0;
 68: }