Actual source code: ex32.c

petsc-3.4.5 2014-06-29
  1: /*T
  2:    Concepts: KSP^solving a system of linear equations
  3:    Concepts: KSP^Laplacian, 2d
  4:    Processors: n
  5: T*/

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

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

 12: with forcing function

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

 16: with pure Neumann boundary conditions

 18: The functions are cell-centered

 20: This uses multigrid to solve the linear system

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

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

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

 30: that vertex based multigrid needs.
 31: */

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

 35: #include <petscdmda.h>
 36: #include <petscksp.h>

 38: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,MatStructure*,void*);
 39: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);

 41: typedef enum {DIRICHLET, NEUMANN} BCType;

 43: typedef struct {
 44:   PetscScalar nu;
 45:   BCType      bcType;
 46: } UserContext;

 50: int main(int argc,char **argv)
 51: {
 52:   KSP            ksp;
 53:   DM             da;
 54:   UserContext    user;
 55:   const char     *bcTypes[2] = {"dirichlet","neumann"};
 57:   PetscInt       bc;

 59:   PetscInitialize(&argc,&argv,(char*)0,help);

 61:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 62:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,12,12,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
 63:   DMDASetInterpolationType(da, DMDA_Q0);

 65:   KSPSetDM(ksp,da);


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

 76:   KSPSetComputeRHS(ksp,ComputeRHS,&user);
 77:   KSPSetComputeOperators(ksp,ComputeMatrix,&user);
 78:   KSPSetFromOptions(ksp);

 80:   KSPSolve(ksp,NULL,NULL);

 82:   KSPDestroy(&ksp);
 83:   DMDestroy(&da);
 84:   PetscFinalize();
 85:   return 0;
 86: }

 90: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
 91: {
 92:   UserContext    *user = (UserContext*)ctx;
 94:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
 95:   PetscScalar    Hx,Hy;
 96:   PetscScalar    **array;
 97:   DM             da;

100:   KSPGetDM(ksp,&da);
101:   DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);
102:   Hx   = 1.0 / (PetscReal)(mx);
103:   Hy   = 1.0 / (PetscReal)(my);
104:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
105:   DMDAVecGetArray(da, b, &array);
106:   for (j=ys; j<ys+ym; j++) {
107:     for (i=xs; i<xs+xm; i++) {
108:       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;
109:     }
110:   }
111:   DMDAVecRestoreArray(da, b, &array);
112:   VecAssemblyBegin(b);
113:   VecAssemblyEnd(b);

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

120:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
121:     MatNullSpaceRemove(nullspace,b,NULL);
122:     MatNullSpaceDestroy(&nullspace);
123:   }
124:   return(0);
125: }


130: PetscErrorCode ComputeMatrix(KSP ksp, Mat J,Mat jac,MatStructure *str, void *ctx)
131: {
132:   UserContext    *user = (UserContext*)ctx;
134:   PetscInt       i,j,mx,my,xm,ym,xs,ys,num, numi, numj;
135:   PetscScalar    v[5],Hx,Hy,HydHx,HxdHy;
136:   MatStencil     row, col[5];
137:   DM             da;

140:   KSPGetDM(ksp,&da);
141:   DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);
142:   Hx    = 1.0 / (PetscReal)(mx);
143:   Hy    = 1.0 / (PetscReal)(my);
144:   HxdHy = Hx/Hy;
145:   HydHx = Hy/Hx;
146:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
147:   for (j=ys; j<ys+ym; j++) {
148:     for (i=xs; i<xs+xm; i++) {
149:       row.i = i; row.j = j;
150:       if (i==0 || j==0 || i==mx-1 || j==my-1) {
151:         if (user->bcType == DIRICHLET) {
152:           v[0] = 2.0*(HxdHy + HydHx);
153:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
154:           SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Dirichlet boundary conditions not supported !\n");
155:         } else if (user->bcType == NEUMANN) {
156:           num = 0; numi=0; numj=0;
157:           if (j!=0) {
158:             v[num] = -HxdHy;
159:             col[num].i = i;
160:             col[num].j = j-1;
161:             num++; numj++;
162:           }
163:           if (i!=0) {
164:             v[num]     = -HydHx;
165:             col[num].i = i-1;
166:             col[num].j = j;
167:             num++; numi++;
168:           }
169:           if (i!=mx-1) {
170:             v[num]     = -HydHx;
171:             col[num].i = i+1;
172:             col[num].j = j;
173:             num++; numi++;
174:           }
175:           if (j!=my-1) {
176:             v[num]     = -HxdHy;
177:             col[num].i = i;
178:             col[num].j = j+1;
179:             num++; numj++;
180:           }
181:           v[num] = (PetscReal)(numj)*HxdHy + (PetscReal)(numi)*HydHx; col[num].i = i;   col[num].j = j;
182:           num++;
183:           MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
184:         }
185:       } else {
186:         v[0] = -HxdHy;              col[0].i = i;   col[0].j = j-1;
187:         v[1] = -HydHx;              col[1].i = i-1; col[1].j = j;
188:         v[2] = 2.0*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
189:         v[3] = -HydHx;              col[3].i = i+1; col[3].j = j;
190:         v[4] = -HxdHy;              col[4].i = i;   col[4].j = j+1;
191:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
192:       }
193:     }
194:   }
195:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
196:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
197:   if (user->bcType == NEUMANN) {
198:     MatNullSpace nullspace;

200:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
201:     MatSetNullSpace(jac,nullspace);
202:     MatNullSpaceDestroy(&nullspace);
203:   }
204:   return(0);
205: }