Actual source code: ex15.c
petsc-3.3-p6 2013-02-11
2: static char help[] = "Tests DM interpolation.\n\n";
4: #include <petscdmda.h>
8: int main(int argc,char **argv)
9: {
10: PetscInt M1 = 3,M2,dof = 1,s = 1,ratio = 2,dim = 1;
12: DM da_c,da_f;
13: Vec v_c,v_f;
14: Mat I;
15: PetscScalar one = 1.0;
16: PetscBool pt;
17: DMDABoundaryType bx = DMDA_BOUNDARY_NONE,by = DMDA_BOUNDARY_NONE;
18:
19: PetscInitialize(&argc,&argv,(char*)0,help);
21: PetscOptionsGetInt(PETSC_NULL,"-dim",&dim,PETSC_NULL);
22: PetscOptionsGetInt(PETSC_NULL,"-M",&M1,PETSC_NULL);
23: PetscOptionsGetInt(PETSC_NULL,"-stencil_width",&s,PETSC_NULL);
24: PetscOptionsGetInt(PETSC_NULL,"-ratio",&ratio,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-dof",&dof,PETSC_NULL);
26: PetscOptionsGetBool(PETSC_NULL,"-periodic",(PetscBool*)&pt,PETSC_NULL);
28: if (pt) {
29: if (dim > 0) bx = DMDA_BOUNDARY_PERIODIC;
30: if (dim > 1) by = DMDA_BOUNDARY_PERIODIC;
31: if (dim > 2) bz = DMDA_BOUNDARY_PERIODIC;
32: }
33: if (bx == DMDA_BOUNDARY_NONE) {
34: M2 = ratio*(M1-1) + 1;
35: } else {
36: M2 = ratio*M1;
37: }
39: /* Set up the array */
40: if (dim == 1) {
41: DMDACreate1d(PETSC_COMM_WORLD,bx,M1,dof,s,PETSC_NULL,&da_c);
42: DMDACreate1d(PETSC_COMM_WORLD,bx,M2,dof,s,PETSC_NULL,&da_f);
43: } else if (dim == 2) {
44: DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M1,M1,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,&da_c);
45: DMDACreate2d(PETSC_COMM_WORLD,bx,by,DMDA_STENCIL_BOX,M2,M2,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,&da_f);
46: } else if (dim == 3) {
47: DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,DMDA_STENCIL_BOX,M1,M1,M1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,PETSC_NULL,&da_c);
48: DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,DMDA_STENCIL_BOX,M2,M2,M2,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,s,PETSC_NULL,PETSC_NULL,PETSC_NULL,&da_f);
49: }
51: DMCreateGlobalVector(da_c,&v_c);
52: DMCreateGlobalVector(da_f,&v_f);
54: VecSet(v_c,one);
55: DMCreateInterpolation(da_c,da_f,&I,PETSC_NULL);
56: MatMult(I,v_c,v_f);
57: VecView(v_f,PETSC_VIEWER_STDOUT_WORLD);
58: MatMultTranspose(I,v_f,v_c);
59: VecView(v_c,PETSC_VIEWER_STDOUT_WORLD);
61: MatDestroy(&I);
62: VecDestroy(&v_c);
63: DMDestroy(&da_c);
64: VecDestroy(&v_f);
65: DMDestroy(&da_f);
66: PetscFinalize();
67: return 0;
68: }
69: