Actual source code: ex17.c
petsc-3.4.2 2013-07-02
2: static char help[] = "Tests DMDA interpolation for coarse DM on a subset of processors.\n\n";
4: #include <petscdmda.h>
8: int main(int argc,char **argv)
9: {
10: PetscInt M = 14,dof = 1,s = 1,ratio = 2,dim = 2;
12: DM da_c,da_f;
13: Vec v_c,v_f;
14: Mat I;
15: PetscScalar one = 1.0;
16: MPI_Comm comm_f, comm_c;
18: PetscInitialize(&argc,&argv,(char*)0,help);
20: PetscOptionsGetInt(NULL,"-dim",&dim,NULL);
21: PetscOptionsGetInt(NULL,"-M",&M,NULL);
22: PetscOptionsGetInt(NULL,"-sw",&s,NULL);
23: PetscOptionsGetInt(NULL,"-ratio",&ratio,NULL);
24: PetscOptionsGetInt(NULL,"-dof",&dof,NULL);
26: comm_f = PETSC_COMM_WORLD;
27: DMDASplitComm2d(comm_f,M,M,s,&comm_c);
29: /* Set up the array */
30: if (dim == 2) {
31: DMDACreate2d(comm_c,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c);
32: M = ratio*(M-1) + 1;
33: DMDACreate2d(comm_f,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_f);
34: }
36: DMCreateGlobalVector(da_c,&v_c);
37: DMCreateGlobalVector(da_f,&v_f);
39: VecSet(v_c,one);
40: DMCreateInterpolation(da_c,da_f,&I,NULL);
41: MatInterpolate(I,v_c,v_f);
42: VecView(v_f,PETSC_VIEWER_STDOUT_(comm_f));
43: MatRestrict(I,v_f,v_c);
44: VecView(v_c,PETSC_VIEWER_STDOUT_(comm_c));
46: MatDestroy(&I);
47: VecDestroy(&v_c);
48: DMDestroy(&da_c);
49: VecDestroy(&v_f);
50: DMDestroy(&da_f);
51: PetscFinalize();
52: return 0;
53: }