Actual source code: ex32.c


  2: /*
  3: Laplacian in 2D. Modeled by the partial differential equation

  5:    div  grad u = f,  0 < x,y < 1,

  7: with forcing function

  9:    f = e^{-(1 - x)^2/\nu} e^{-(1 - y)^2/\nu}

 11: with pure Neumann boundary conditions

 13: The functions are cell-centered

 15: This uses multigrid to solve the linear system

 17:        Contributed by Andrei Draganescu <aidraga@sandia.gov>

 19: Note the nice multigrid convergence despite the fact it is only using
 20: peicewise constant interpolation/restriction. This is because cell-centered multigrid
 21: does not need the same rule:

 23:     polynomial degree(interpolation) + polynomial degree(restriction) + 2 > degree of PDE

 25: that vertex based multigrid needs.
 26: */

 28: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";

 30: #include <petscdm.h>
 31: #include <petscdmda.h>
 32: #include <petscksp.h>

 34: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
 35: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);

 37: typedef enum {DIRICHLET, NEUMANN} BCType;

 39: typedef struct {
 40:   PetscScalar nu;
 41:   BCType      bcType;
 42: } UserContext;

 44: int main(int argc,char **argv)
 45: {
 46:   KSP            ksp;
 47:   DM             da;
 48:   UserContext    user;
 49:   const char     *bcTypes[2] = {"dirichlet","neumann"};
 51:   PetscInt       bc;

 53:   PetscInitialize(&argc,&argv,(char*)0,help);
 54:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 55:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,12,12,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
 56:   DMSetFromOptions(da);
 57:   DMSetUp(da);
 58:   DMDASetInterpolationType(da, DMDA_Q0);

 60:   KSPSetDM(ksp,da);

 62:   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DM");
 63:   user.nu     = 0.1;
 64:   PetscOptionsScalar("-nu", "The width of the Gaussian source", "ex29.c", 0.1, &user.nu, NULL);
 65:   bc          = (PetscInt)NEUMANN;
 66:   PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);
 67:   user.bcType = (BCType)bc;
 68:   PetscOptionsEnd();

 70:   KSPSetComputeRHS(ksp,ComputeRHS,&user);
 71:   KSPSetComputeOperators(ksp,ComputeMatrix,&user);
 72:   KSPSetFromOptions(ksp);
 73:   KSPSolve(ksp,NULL,NULL);
 74:   KSPDestroy(&ksp);
 75:   DMDestroy(&da);
 76:   PetscFinalize();
 77:   return 0;
 78: }

 80: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
 81: {
 82:   UserContext    *user = (UserContext*)ctx;
 83:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
 84:   PetscScalar    Hx,Hy;
 85:   PetscScalar    **array;
 86:   DM             da;

 89:   KSPGetDM(ksp,&da);
 90:   DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
 91:   Hx   = 1.0 / (PetscReal)(mx);
 92:   Hy   = 1.0 / (PetscReal)(my);
 93:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
 94:   DMDAVecGetArray(da, b, &array);
 95:   for (j=ys; j<ys+ym; j++) {
 96:     for (i=xs; i<xs+xm; i++) {
 97:       array[j][i] = PetscExpScalar(-(((PetscReal)i+0.5)*Hx)*(((PetscReal)i+0.5)*Hx)/user->nu)*PetscExpScalar(-(((PetscReal)j+0.5)*Hy)*(((PetscReal)j+0.5)*Hy)/user->nu)*Hx*Hy;
 98:     }
 99:   }
100:   DMDAVecRestoreArray(da, b, &array);
101:   VecAssemblyBegin(b);
102:   VecAssemblyEnd(b);

104:   /* force right hand side to be consistent for singular matrix */
105:   /* note this is really a hack, normally the model would provide you with a consistent right handside */
106:   if (user->bcType == NEUMANN) {
107:     MatNullSpace nullspace;

109:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
110:     MatNullSpaceRemove(nullspace,b);
111:     MatNullSpaceDestroy(&nullspace);
112:   }
113:   return 0;
114: }

116: PetscErrorCode ComputeMatrix(KSP ksp, Mat J,Mat jac, void *ctx)
117: {
118:   UserContext    *user = (UserContext*)ctx;
119:   PetscInt       i,j,mx,my,xm,ym,xs,ys,num, numi, numj;
120:   PetscScalar    v[5],Hx,Hy,HydHx,HxdHy;
121:   MatStencil     row, col[5];
122:   DM             da;

125:   KSPGetDM(ksp,&da);
126:   DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
127:   Hx    = 1.0 / (PetscReal)(mx);
128:   Hy    = 1.0 / (PetscReal)(my);
129:   HxdHy = Hx/Hy;
130:   HydHx = Hy/Hx;
131:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
132:   for (j=ys; j<ys+ym; j++) {
133:     for (i=xs; i<xs+xm; i++) {
134:       row.i = i; row.j = j;
135:       if (i==0 || j==0 || i==mx-1 || j==my-1) {
136:         if (user->bcType == DIRICHLET) {
137:           v[0] = 2.0*(HxdHy + HydHx);
138:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
139:           SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Dirichlet boundary conditions not supported !");
140:         } else if (user->bcType == NEUMANN) {
141:           num = 0; numi=0; numj=0;
142:           if (j!=0) {
143:             v[num] = -HxdHy;
144:             col[num].i = i;
145:             col[num].j = j-1;
146:             num++; numj++;
147:           }
148:           if (i!=0) {
149:             v[num]     = -HydHx;
150:             col[num].i = i-1;
151:             col[num].j = j;
152:             num++; numi++;
153:           }
154:           if (i!=mx-1) {
155:             v[num]     = -HydHx;
156:             col[num].i = i+1;
157:             col[num].j = j;
158:             num++; numi++;
159:           }
160:           if (j!=my-1) {
161:             v[num]     = -HxdHy;
162:             col[num].i = i;
163:             col[num].j = j+1;
164:             num++; numj++;
165:           }
166:           v[num] = (PetscReal)(numj)*HxdHy + (PetscReal)(numi)*HydHx; col[num].i = i;   col[num].j = j;
167:           num++;
168:           MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
169:         }
170:       } else {
171:         v[0] = -HxdHy;              col[0].i = i;   col[0].j = j-1;
172:         v[1] = -HydHx;              col[1].i = i-1; col[1].j = j;
173:         v[2] = 2.0*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
174:         v[3] = -HydHx;              col[3].i = i+1; col[3].j = j;
175:         v[4] = -HxdHy;              col[4].i = i;   col[4].j = j+1;
176:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
177:       }
178:     }
179:   }
180:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
181:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
182:   if (user->bcType == NEUMANN) {
183:     MatNullSpace nullspace;

185:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
186:     MatSetNullSpace(J,nullspace);
187:     MatNullSpaceDestroy(&nullspace);
188:   }
189:   return 0;
190: }

192: /*TEST

194:    test:
195:       args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -pc_mg_levels 3 -mg_coarse_pc_factor_shift_type nonzero

197: TEST*/