Actual source code: ex17.c
1: /*$Id: ex17.c,v 1.9 2001/08/07 03:04:42 balay Exp $*/
3: static char help[] = "Tests DA interpolation for coarse DA on a subset of processors.nn";
5: #include petscda.h
6: #include petscsys.h
8: int main(int argc,char **argv)
9: {
10: int M = 14,ierr,dof = 1,s = 1,ratio = 2,dim = 2;
11: DA 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);
19: PetscOptionsGetInt(PETSC_NULL,"-dim",&dim,PETSC_NULL);
20: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
21: PetscOptionsGetInt(PETSC_NULL,"-sw",&s,PETSC_NULL);
22: PetscOptionsGetInt(PETSC_NULL,"-ratio",&ratio,PETSC_NULL);
23: PetscOptionsGetInt(PETSC_NULL,"-dof",&dof,PETSC_NULL);
25: comm_f = PETSC_COMM_WORLD;
26: DASplitComm2d(comm_f,M,M,s,&comm_c);
27:
28: /* Set up the array */
29: if (dim == 2) {
30: DACreate2d(comm_c,DA_NONPERIODIC,DA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,&da_c);
31: M = ratio*(M-1) + 1;
32: DACreate2d(comm_f,DA_NONPERIODIC,DA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,&da_f);
33: } else if (dim == 3) {
34: ;
35: }
37: DACreateGlobalVector(da_c,&v_c);
38: DACreateGlobalVector(da_f,&v_f);
40: VecSet(&one,v_c);
41: DAGetInterpolation(da_c,da_f,&I,PETSC_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: DADestroy(da_c);
50: VecDestroy(v_f);
51: DADestroy(da_f);
52: PetscFinalize();
53: return 0;
54: }
55: