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