Actual source code: ex29.c

  1: /*T
  2:    Concepts: KSP^solving a system of linear equations
  3:    Concepts: KSP^Laplacian, 2d
  4:    Processors: n
  5: T*/

  7: /*
  8: Added at the request of Marc Garbey.

 10: Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation

 12:    div \rho grad u = f,  0 < x,y < 1,

 14: with forcing function

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

 18: with Dirichlet boundary conditions

 20:    u = f(x,y) for x = 0, x = 1, y = 0, y = 1

 22: or pure Neumman boundary conditions

 24: This uses multigrid to solve the linear system
 25: */

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

 29:  #include petscda.h
 30:  #include petscksp.h
 31:  #include petscmg.h
 32:  #include petscdmmg.h


 38: typedef enum {DIRICHLET, NEUMANN} BCType;

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

 47: int main(int argc,char **argv)
 48: {
 49:   DMMG           *dmmg;
 50:   DA             da;
 51:   UserContext    user;
 52:   PetscReal      norm;
 53:   PetscScalar    mone = -1.0;
 54:   const char     *bcTypes[2] = {"dirichlet","neumann"};
 56:   PetscInt       l,bc;

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

 60:   DMMGCreate(PETSC_COMM_WORLD,3,PETSC_NULL,&dmmg);
 61:   DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-3,-3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);
 62:   DMMGSetDM(dmmg,(DM)da);
 63:   DADestroy(da);
 64:   for (l = 0; l < DMMGGetLevels(dmmg); l++) {
 65:     DMMGSetUser(dmmg,l,&user);
 66:   }

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

 76:   DMMGSetKSP(dmmg,ComputeRHS,ComputeJacobian);
 77:   if (user.bcType == NEUMANN) {
 78:     DMMGSetNullSpace(dmmg,PETSC_TRUE,0,PETSC_NULL);
 79:   }

 81:   DMMGSolve(dmmg);

 83:   MatMult(DMMGGetJ(dmmg),DMMGGetx(dmmg),DMMGGetr(dmmg));
 84:   VecAXPY(DMMGGetr(dmmg),mone,DMMGGetRHS(dmmg));
 85:   VecNorm(DMMGGetr(dmmg),NORM_2,&norm);
 86:   /* PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",norm); */
 87:   VecAssemblyBegin(DMMGGetx(dmmg));
 88:   VecAssemblyEnd(DMMGGetx(dmmg));
 89:   VecView_VTK(DMMGGetRHS(dmmg), "rhs.vtk", bcTypes[user.bcType]);
 90:   VecView_VTK(DMMGGetx(dmmg), "solution.vtk", bcTypes[user.bcType]);

 92:   DMMGDestroy(dmmg);
 93:   PetscFinalize();

 95:   return 0;
 96: }

100: PetscErrorCode ComputeRHS(DMMG dmmg, Vec b)
101: {
102:   DA             da = (DA)dmmg->dm;
103:   UserContext    *user = (UserContext *) dmmg->user;
105:   PetscInt       i,j,mx,my,xm,ym,xs,ys;
106:   PetscScalar    Hx,Hy;
107:   PetscScalar    **array;

110:   DAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0);
111:   Hx   = 1.0 / (PetscReal)(mx-1);
112:   Hy   = 1.0 / (PetscReal)(my-1);
113:   DAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
114:   DAVecGetArray(da, b, &array);
115:   for (j=ys; j<ys+ym; j++){
116:     for(i=xs; i<xs+xm; i++){
117:       array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
118:     }
119:   }
120:   DAVecRestoreArray(da, b, &array);
121:   VecAssemblyBegin(b);
122:   VecAssemblyEnd(b);

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

129:     KSPGetNullSpace(dmmg->ksp,&nullspace);
130:     MatNullSpaceRemove(nullspace,b,PETSC_NULL);
131:   }
132:   return(0);
133: }

135: 
138: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscScalar *rho)
139: {
141:   if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
142:     *rho = 100.0;
143:   } else {
144:     *rho = 1.0;
145:   }
146:   return(0);
147: }

151: PetscErrorCode ComputeJacobian(DMMG dmmg, Mat jac)
152: {
153:   DA             da = (DA) dmmg->dm;
154:   UserContext    *user = (UserContext *) dmmg->user;
156:   PetscInt       i,j,mx,my,xm,ym,xs,ys,num;
157:   PetscScalar    v[5],Hx,Hy,HydHx,HxdHy,rho;
158:   MatStencil     row, col[5];

161:   DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);
162:   Hx    = 1.0 / (PetscReal)(mx-1);
163:   Hy    = 1.0 / (PetscReal)(my-1);
164:   HxdHy = Hx/Hy;
165:   HydHx = Hy/Hx;
166:   DAGetCorners(da,&xs,&ys,0,&xm,&ym,0);
167:   for (j=ys; j<ys+ym; j++){
168:     for(i=xs; i<xs+xm; i++){
169:       row.i = i; row.j = j;
170:       if (i==0 || j==0 || i==mx-1 || j==my-1) {
171:        ComputeRho(i, j, mx, my, &rho);
172:        if (user->bcType == DIRICHLET) {
173:            v[0] = 2.0*rho*(HxdHy + HydHx);
174:           MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
175:         } else if (user->bcType == NEUMANN) {
176:           num = 0;
177:           if (j!=0) {
178:             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j-1;
179:             num++;
180:           }
181:           if (i!=0) {
182:             v[num] = -rho*HydHx;              col[num].i = i-1; col[num].j = j;
183:             num++;
184:           }
185:           if (i!=mx-1) {
186:             v[num] = -rho*HydHx;              col[num].i = i+1; col[num].j = j;
187:             num++;
188:           }
189:           if (j!=my-1) {
190:             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j+1;
191:             num++;
192:           }
193:           v[num]   = (num/2.0)*rho*(HxdHy + HydHx); col[num].i = i;   col[num].j = j;
194:           num++;
195:           MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
196:         }
197:       } else {
198:         v[0] = -rho*HxdHy;              col[0].i = i;   col[0].j = j-1;
199:         v[1] = -rho*HydHx;              col[1].i = i-1; col[1].j = j;
200:         v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
201:         v[3] = -rho*HydHx;              col[3].i = i+1; col[3].j = j;
202:         v[4] = -rho*HxdHy;              col[4].i = i;   col[4].j = j+1;
203:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
204:       }
205:     }
206:   }
207:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
208:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
209:   return(0);
210: }

214: PetscErrorCode VecView_VTK(Vec x, const char filename[], const char bcName[])
215: {
216:   MPI_Comm           comm;
217:   DA                 da;
218:   Vec                coords;
219:   PetscViewer        viewer;
220:   PetscScalar       *array, *values;
221:   PetscInt           n, N, maxn, mx, my, dof;
222:   PetscInt           i, p;
223:   MPI_Status         status;
224:   PetscMPIInt        rank, size, tag;
225:   PetscErrorCode     ierr;

228:   PetscObjectGetComm((PetscObject) x, &comm);
229:   PetscViewerASCIIOpen(comm, filename, &viewer);

231:   VecGetSize(x, &N);
232:   VecGetLocalSize(x, &n);
233:   PetscObjectQuery((PetscObject) x, "DA", (PetscObject *) &da);
234:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");

236:   DAGetInfo(da, 0, &mx, &my, 0,0,0,0, &dof,0,0,0);

238:   PetscViewerASCIIPrintf(viewer, "# vtk DataFile Version 2.0\n");
239:   PetscViewerASCIIPrintf(viewer, "Inhomogeneous Poisson Equation with %s boundary conditions\n", bcName);
240:   PetscViewerASCIIPrintf(viewer, "ASCII\n");
241:   /* get coordinates of nodes */
242:   DAGetCoordinates(da, &coords);
243:   if (!coords) {
244:     DASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0);
245:     DAGetCoordinates(da, &coords);
246:   }
247:   PetscViewerASCIIPrintf(viewer, "DATASET RECTILINEAR_GRID\n");
248:   PetscViewerASCIIPrintf(viewer, "DIMENSIONS %d %d %d\n", mx, my, 1);
249:   VecGetArray(coords, &array);
250:   PetscViewerASCIIPrintf(viewer, "X_COORDINATES %d double\n", mx);
251:   for(i = 0; i < mx; i++) {
252:     PetscViewerASCIIPrintf(viewer, "%g ", PetscRealPart(array[i*2]));
253:   }
254:   PetscViewerASCIIPrintf(viewer, "\n");
255:   PetscViewerASCIIPrintf(viewer, "Y_COORDINATES %d double\n", my);
256:   for(i = 0; i < my; i++) {
257:     PetscViewerASCIIPrintf(viewer, "%g ", PetscRealPart(array[i*mx*2+1]));
258:   }
259:   PetscViewerASCIIPrintf(viewer, "\n");
260:   PetscViewerASCIIPrintf(viewer, "Z_COORDINATES %d double\n", 1);
261:   PetscViewerASCIIPrintf(viewer, "%g\n", 0.0);
262:   VecRestoreArray(coords, &array);
263:   PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", N);
264:   PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", dof);
265:   PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");
266:   VecGetArray(x, &array);
267:   /* Determine maximum message to arrive */
268:   MPI_Comm_rank(comm, &rank);
269:   MPI_Comm_size(comm, &size);
270:   MPI_Reduce(&n, &maxn, 1, MPIU_INT, MPI_MAX, 0, comm);
271:   tag  = ((PetscObject) viewer)->tag;
272:   if (!rank) {
273:     PetscMalloc((maxn+1) * sizeof(PetscScalar), &values);
274:     for(i = 0; i < n; i++) {
275:       PetscViewerASCIIPrintf(viewer, "%g\n", PetscRealPart(array[i]));
276:     }
277:     for(p = 1; p < size; p++) {
278:       MPI_Recv(values, (PetscMPIInt) n, MPIU_SCALAR, p, tag, comm, &status);
279:       MPI_Get_count(&status, MPIU_SCALAR, &n);
280:       for(i = 0; i < n; i++) {
281:         PetscViewerASCIIPrintf(viewer, "%g\n", PetscRealPart(array[i]));
282:       }
283:     }
284:     PetscFree(values);
285:   } else {
286:     MPI_Send(array, n, MPIU_SCALAR, 0, tag, comm);
287:   }
288:   VecRestoreArray(x, &array);
289:   PetscViewerFlush(viewer);
290:   PetscViewerDestroy(viewer);
291:   return(0);
292: }