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.\n\n";
5: #include petscda.h
6: #include petscsys.h
10: int main(int argc,char **argv)
11: {
12: int M = 14,ierr,dof = 1,s = 1,ratio = 2,dim = 2;
13: DA 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(PETSC_NULL,"-dim",&dim,PETSC_NULL);
22: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
23: PetscOptionsGetInt(PETSC_NULL,"-sw",&s,PETSC_NULL);
24: PetscOptionsGetInt(PETSC_NULL,"-ratio",&ratio,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-dof",&dof,PETSC_NULL);
27: comm_f = PETSC_COMM_WORLD;
28: DASplitComm2d(comm_f,M,M,s,&comm_c);
29:
30: /* Set up the array */
31: if (dim == 2) {
32: DACreate2d(comm_c,DA_NONPERIODIC,DA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,&da_c);
33: M = ratio*(M-1) + 1;
34: DACreate2d(comm_f,DA_NONPERIODIC,DA_STENCIL_BOX,M,M,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,&da_f);
35: } else if (dim == 3) {
36: ;
37: }
39: DACreateGlobalVector(da_c,&v_c);
40: DACreateGlobalVector(da_f,&v_f);
42: VecSet(&one,v_c);
43: DAGetInterpolation(da_c,da_f,&I,PETSC_NULL);
44: MatInterpolate(I,v_c,v_f);
45: VecView(v_f,PETSC_VIEWER_STDOUT_(comm_f));
46: MatRestrict(I,v_f,v_c);
47: VecView(v_c,PETSC_VIEWER_STDOUT_(comm_c));
49: MatDestroy(I);
50: VecDestroy(v_c);
51: DADestroy(da_c);
52: VecDestroy(v_f);
53: DADestroy(da_f);
54: PetscFinalize();
55: return 0;
56: }
57: