Actual source code: ex5.c

  2: /* Program usage:  mpirun -np <procs> ex5 [-help] [all PETSc options] */

  4: static char help[] = "Bratu nonlinear PDE in 2d.\n\
  5: We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
  6: domain, using distributed arrays (DAs) to partition the parallel grid.\n\
  7: The command line options include:\n\
  8:   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
  9:      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n";

 11: /*T
 12:    Concepts: SNES^parallel Bratu example
 13:    Concepts: DA^using distributed arrays;
 14:    Processors: n
 15: T*/

 17: /* ------------------------------------------------------------------------

 19:     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
 20:     the partial differential equation
 21:   
 22:             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
 23:   
 24:     with boundary conditions
 25:    
 26:              u = 0  for  x = 0, x = 1, y = 0, y = 1.
 27:   
 28:     A finite difference approximation with the usual 5-point stencil
 29:     is used to discretize the boundary value problem to obtain a nonlinear 
 30:     system of equations.


 33:   ------------------------------------------------------------------------- */

 35: /* 
 36:    Include "petscda.h" so that we can use distributed arrays (DAs).
 37:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 38:    file automatically includes:
 39:      petsc.h       - base PETSc routines   petscvec.h - vectors
 40:      petscsys.h    - system routines       petscmat.h - matrices
 41:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 42:      petscviewer.h - viewers               petscpc.h  - preconditioners
 43:      petscksp.h   - linear solvers
 44: */
 45:  #include petscda.h
 46:  #include petscsnes.h

 48: /* 
 49:    User-defined application context - contains data needed by the 
 50:    application-provided call-back routines, FormJacobianLocal() and
 51:    FormFunctionLocal().
 52: */
 53: typedef struct {
 54:    DA          da;             /* distributed array data structure */
 55:    PassiveReal param;          /* test problem parameter */
 56: } AppCtx;

 58: /* 
 59:    User-defined routines
 60: */

 67: int main(int argc,char **argv)
 68: {
 69:   SNES                   snes;                 /* nonlinear solver */
 70:   Vec                    x,r;                  /* solution, residual vectors */
 71:   Mat                    A,J;                    /* Jacobian matrix */
 72:   AppCtx                 user;                 /* user-defined work context */
 73:   PetscInt               its;                  /* iterations for convergence */
 74:   PetscTruth             matlab_function = PETSC_FALSE;
 75:   PetscTruth             fd_jacobian = PETSC_FALSE,adic_jacobian=PETSC_FALSE;
 76:   PetscTruth             adicmf_jacobian = PETSC_FALSE;
 77:   PetscErrorCode         ierr;
 78:   PetscReal              bratu_lambda_max = 6.81,bratu_lambda_min = 0.;
 79:   MatFDColoring          matfdcoloring = 0;
 80:   ISColoring             iscoloring;

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:      Initialize program
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:      Initialize problem parameters
 90:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 91:   user.param = 6.0;
 92:   PetscOptionsGetReal(PETSC_NULL,"-par",&user.param,PETSC_NULL);
 93:   if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) {
 94:     SETERRQ(1,"Lambda is out of range");
 95:   }

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:      Create nonlinear solver context
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100:   SNESCreate(PETSC_COMM_WORLD,&snes);

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:      Create distributed array (DA) to manage parallel grid and vectors
104:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105:   DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,
106:                     1,1,PETSC_NULL,PETSC_NULL,&user.da);

108:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:      Extract global vectors from DA; then duplicate for remaining
110:      vectors that are the same types
111:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112:   DACreateGlobalVector(user.da,&x);
113:   VecDuplicate(x,&r);

115:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116:      Create matrix data structure; set Jacobian evaluation routine

118:      Set Jacobian matrix data structure and default Jacobian evaluation
119:      routine. User can override with:
120:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
121:                 (unless user explicitly sets preconditioner) 
122:      -snes_mf_operator : form preconditioning matrix as set by the user,
123:                          but use matrix-free approx for Jacobian-vector
124:                          products within Newton-Krylov method

126:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127:   /* J can be type of MATAIJ,MATBAIJ or MATSBAIJ */
128:   DAGetMatrix(user.da,MATAIJ,&J);
129: 
130:   A    = J;
131:   PetscOptionsGetTruth(PETSC_NULL,"-fd_jacobian",&fd_jacobian,0);
132:   PetscOptionsGetTruth(PETSC_NULL,"-adic_jacobian",&adic_jacobian,0);
133:   PetscOptionsGetTruth(PETSC_NULL,"-adicmf_jacobian",&adicmf_jacobian,0);
134: #if defined(PETSC_HAVE_ADIC)
135:   if (adicmf_jacobian) {
136:     DASetLocalAdicMFFunction(user.da,admf_FormFunctionLocal);
137:     MatRegisterDAAD();
138:     MatCreateDAAD(user.da,&A);
139:     MatDAADSetSNES(A,snes);
140:     MatDAADSetCtx(A,&user);
141:   }
142: #endif

144:   if (fd_jacobian) {
145:     DAGetColoring(user.da,IS_COLORING_LOCAL,&iscoloring);
146:     MatFDColoringCreate(J,iscoloring,&matfdcoloring);
147:     ISColoringDestroy(iscoloring);
148:     MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESDAFormFunction,&user);
149:     MatFDColoringSetFromOptions(matfdcoloring);
150:     SNESSetJacobian(snes,A,J,SNESDefaultComputeJacobianColor,matfdcoloring);
151: #if defined(PETSC_HAVE_ADIC)
152:   } else if (adic_jacobian) {
153:     DAGetColoring(user.da,IS_COLORING_GHOSTED,&iscoloring);
154:     MatSetColoring(J,iscoloring);
155:     ISColoringDestroy(iscoloring);
156:     SNESSetJacobian(snes,A,J,SNESDAComputeJacobianWithAdic,&user);
157: #endif
158:   } else {
159:     SNESSetJacobian(snes,A,J,SNESDAComputeJacobian,&user);
160:   }

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Set local function evaluation routine
164:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   DASetLocalFunction(user.da,(DALocalFunction1)FormFunctionLocal);
166:   DASetLocalJacobian(user.da,(DALocalFunction1)FormJacobianLocal);
167:   DASetLocalAdicFunction(user.da,ad_FormFunctionLocal);

169:   /* Decide which FormFunction to use */
170:   PetscOptionsGetTruth(PETSC_NULL,"-matlab_function",&matlab_function,0);

172:   SNESSetFunction(snes,r,SNESDAFormFunction,&user);
173: #if defined(PETSC_HAVE_MATLAB)
174:   if (matlab_function) {
175:     SNESSetFunction(snes,r,FormFunctionMatlab,&user);
176:   }
177: #endif

179:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:      Customize nonlinear solver; set runtime options
181:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182:   SNESSetFromOptions(snes);

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:      Evaluate initial guess
186:      Note: The user should initialize the vector, x, with the initial guess
187:      for the nonlinear solver prior to calling SNESSolve().  In particular,
188:      to employ an initial guess of zero, the user should explicitly set
189:      this vector to zero by calling VecSet().
190:   */
191:   FormInitialGuess(&user,x);

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:      Solve nonlinear system
195:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196:   SNESSolve(snes,PETSC_NULL,x);
197:   SNESGetIterationNumber(snes,&its);

199:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201:   PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D\n",its);

203:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:      Free work space.  All PETSc objects should be destroyed when they
205:      are no longer needed.
206:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

208:   if (A != J) {
209:     MatDestroy(A);
210:   }
211:   MatDestroy(J);
212:   if (matfdcoloring) {
213:     MatFDColoringDestroy(matfdcoloring);
214:   }
215:   VecDestroy(x);
216:   VecDestroy(r);
217:   SNESDestroy(snes);
218:   DADestroy(user.da);
219:   PetscFinalize();

221:   return(0);
222: }
223: /* ------------------------------------------------------------------- */
226: /* 
227:    FormInitialGuess - Forms initial approximation.

229:    Input Parameters:
230:    user - user-defined application context
231:    X - vector

233:    Output Parameter:
234:    X - vector
235:  */
236: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
237: {
238:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
240:   PetscReal      lambda,temp1,temp,hx,hy;
241:   PetscScalar    **x;

244:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
245:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

247:   lambda = user->param;
248:   hx     = 1.0/(PetscReal)(Mx-1);
249:   hy     = 1.0/(PetscReal)(My-1);
250:   temp1  = lambda/(lambda + 1.0);

252:   /*
253:      Get a pointer to vector data.
254:        - For default PETSc vectors, VecGetArray() returns a pointer to
255:          the data array.  Otherwise, the routine is implementation dependent.
256:        - You MUST call VecRestoreArray() when you no longer need access to
257:          the array.
258:   */
259:   DAVecGetArray(user->da,X,&x);

261:   /*
262:      Get local grid boundaries (for 2-dimensional DA):
263:        xs, ys   - starting grid indices (no ghost points)
264:        xm, ym   - widths of local grid (no ghost points)

266:   */
267:   DAGetCorners(user->da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

269:   /*
270:      Compute initial guess over the locally owned part of the grid
271:   */
272:   for (j=ys; j<ys+ym; j++) {
273:     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
274:     for (i=xs; i<xs+xm; i++) {

276:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
277:         /* boundary conditions are all zero Dirichlet */
278:         x[j][i] = 0.0;
279:       } else {
280:         x[j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
281:       }
282:     }
283:   }

285:   /*
286:      Restore vector
287:   */
288:   DAVecRestoreArray(user->da,X,&x);

290:   return(0);
291: }
292: /* ------------------------------------------------------------------- */
295: /* 
296:    FormFunctionLocal - Evaluates nonlinear function, F(x).

298:        Process adiC(36): FormFunctionLocal

300:  */
301: PetscErrorCode FormFunctionLocal(DALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
302: {
304:   PetscInt       i,j;
305:   PetscReal      two = 2.0,lambda,hx,hy,hxdhy,hydhx,sc;
306:   PetscScalar    u,uxx,uyy;


310:   lambda = user->param;
311:   hx     = 1.0/(PetscReal)(info->mx-1);
312:   hy     = 1.0/(PetscReal)(info->my-1);
313:   sc     = hx*hy*lambda;
314:   hxdhy  = hx/hy;
315:   hydhx  = hy/hx;
316:   /*
317:      Compute function over the locally owned part of the grid
318:   */
319:   for (j=info->ys; j<info->ys+info->ym; j++) {
320:     for (i=info->xs; i<info->xs+info->xm; i++) {
321:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
322:         f[j][i] = x[j][i];
323:       } else {
324:         u       = x[j][i];
325:         uxx     = (two*u - x[j][i-1] - x[j][i+1])*hydhx;
326:         uyy     = (two*u - x[j-1][i] - x[j+1][i])*hxdhy;
327:         f[j][i] = uxx + uyy - sc*PetscExpScalar(u);
328:       }
329:     }
330:   }

332:   PetscLogFlops(11*info->ym*info->xm);
333:   return(0);
334: }

338: /*
339:    FormJacobianLocal - Evaluates Jacobian matrix.
340: */
341: PetscErrorCode FormJacobianLocal(DALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user)
342: {
344:   PetscInt       i,j;
345:   MatStencil     col[5],row;
346:   PetscScalar    lambda,v[5],hx,hy,hxdhy,hydhx,sc;

349:   lambda = user->param;
350:   hx     = 1.0/(PetscReal)(info->mx-1);
351:   hy     = 1.0/(PetscReal)(info->my-1);
352:   sc     = hx*hy*lambda;
353:   hxdhy  = hx/hy;
354:   hydhx  = hy/hx;


357:   /* 
358:      Compute entries for the locally owned part of the Jacobian.
359:       - Currently, all PETSc parallel matrix formats are partitioned by
360:         contiguous chunks of rows across the processors. 
361:       - Each processor needs to insert only elements that it owns
362:         locally (but any non-local elements will be sent to the
363:         appropriate processor during matrix assembly). 
364:       - Here, we set all entries for a particular row at once.
365:       - We can set matrix entries either using either
366:         MatSetValuesLocal() or MatSetValues(), as discussed above.
367:   */
368:   for (j=info->ys; j<info->ys+info->ym; j++) {
369:     for (i=info->xs; i<info->xs+info->xm; i++) {
370:       row.j = j; row.i = i;
371:       /* boundary points */
372:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
373:         v[0] = 1.0;
374:         MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
375:       } else {
376:       /* interior grid points */
377:         v[0] = -hxdhy;                                           col[0].j = j - 1; col[0].i = i;
378:         v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
379:         v[2] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[2].j = row.j; col[2].i = row.i;
380:         v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
381:         v[4] = -hxdhy;                                           col[4].j = j + 1; col[4].i = i;
382:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
383:       }
384:     }
385:   }

387:   /* 
388:      Assemble matrix, using the 2-step process:
389:        MatAssemblyBegin(), MatAssemblyEnd().
390:   */
391:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
392:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
393:   /*
394:      Tell the matrix we will never add a new nonzero location to the
395:      matrix. If we do, it will generate an error.
396:   */
397:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
398:   return(0);
399: }


402: /*
403:       Variant of FormFunction() that computes the function in Matlab
404: */
405: #if defined(PETSC_HAVE_MATLAB)
406: PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr)
407: {
408:   AppCtx         *user = (AppCtx*)ptr;
410:   PetscInt       Mx,My;
411:   PetscReal      lambda,hx,hy;
412:   Vec            localX,localF;
413:   MPI_Comm       comm;

416:   DAGetLocalVector(user->da,&localX);
417:   DAGetLocalVector(user->da,&localF);
418:   PetscObjectSetName((PetscObject)localX,"localX");
419:   PetscObjectSetName((PetscObject)localF,"localF");
420:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
421:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

423:   lambda = user->param;
424:   hx     = 1.0/(PetscReal)(Mx-1);
425:   hy     = 1.0/(PetscReal)(My-1);

427:   PetscObjectGetComm((PetscObject)snes,&comm);
428:   /*
429:      Scatter ghost points to local vector,using the 2-step process
430:         DAGlobalToLocalBegin(),DAGlobalToLocalEnd().
431:      By placing code between these two statements, computations can be
432:      done while messages are in transition.
433:   */
434:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
435:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
436:   PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX);
437:   PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",hx,hy,lambda);
438:   PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF);

440:   /*
441:      Insert values into global vector
442:   */
443:   DALocalToGlobal(user->da,localF,INSERT_VALUES,F);
444:   DARestoreLocalVector(user->da,&localX);
445:   DARestoreLocalVector(user->da,&localF);
446:   return(0);
447: }
448: #endif