Actual source code: ex17.c
petsc-3.9.4 2018-09-11
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>
7: int main(int argc,char **argv)
8: {
9: PetscInt M = 14,dof = 1,s = 1,ratio = 2,dim = 2;
11: DM da_c,da_f;
12: Vec v_c,v_f;
13: Mat I;
14: PetscScalar one = 1.0;
15: MPI_Comm comm_f, comm_c;
17: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
18: PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);
19: PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
20: PetscOptionsGetInt(NULL,NULL,"-sw",&s,NULL);
21: PetscOptionsGetInt(NULL,NULL,"-ratio",&ratio,NULL);
22: PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);
24: comm_f = PETSC_COMM_WORLD;
25: DMDASplitComm2d(comm_f,M,M,s,&comm_c);
27: /* Set up the array */
28: if (dim == 2) {
29: DMDACreate2d(comm_c,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_c);
30: M = ratio*(M-1) + 1;
31: DMDACreate2d(comm_f,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,NULL,NULL,&da_f);
32: }
33: DMSetFromOptions(da_c);
34: DMSetUp(da_c);CHKERRQ(ierr)
35: DMSetFromOptions(da_f);
36: DMSetUp(da_f);
38: DMCreateGlobalVector(da_c,&v_c);
39: DMCreateGlobalVector(da_f,&v_f);
41: VecSet(v_c,one);
42: DMCreateInterpolation(da_c,da_f,&I,NULL);
43: MatInterpolate(I,v_c,v_f);
44: VecView(v_f,PETSC_VIEWER_STDOUT_(comm_f));
45: MatRestrict(I,v_f,v_c);
46: VecView(v_c,PETSC_VIEWER_STDOUT_(comm_c));
48: MatDestroy(&I);
49: VecDestroy(&v_c);
50: DMDestroy(&da_c);
51: VecDestroy(&v_f);
52: DMDestroy(&da_f);
53: PetscFinalize();
54: return ierr;
55: }