Actual source code: ex32.c

petsc-master 2016-05-03
Report Typos and Errors
  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 <petscdm.h>
 36: #include <petscdmda.h>
 37: #include <petscksp.h>

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

 42: typedef enum {DIRICHLET, NEUMANN} BCType;

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

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

 60:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 61:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 62:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_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);
 79:   KSPSolve(ksp,NULL,NULL);
 80:   KSPDestroy(&ksp);
 81:   DMDestroy(&da);
 82:   PetscFinalize();
 83:   return ierr;
 84: }

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

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

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

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


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

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

198:     MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
199:     MatSetNullSpace(J,nullspace);
200:     MatNullSpaceDestroy(&nullspace);
201:   }
202:   return(0);
203: }