Actual source code: ex26.c

  2: static char help[] = "Grad-Shafranov solver for one dimensional CHI equilibrium.\n\
  3: The command line options include:\n\
  4:   -n <n> number of grid points\n\
  5:   -psi_axis <axis> \n\
  6:   -r_min <min> \n\
  7:   -param <param> \n\n";

  9: /*T
 10:    Concepts: SNES^parallel CHI equilibrium
 11:    Concepts: DA^using distributed arrays;
 12:    Processors: n
 13: T*/

 15: /* ------------------------------------------------------------------------

 17:    Grad-Shafranov solver for one dimensional CHI equilibrium
 18:   
 19:     A finite difference approximation with the usual 3-point stencil
 20:     is used to discretize the boundary value problem to obtain a nonlinear 
 21:     system of equations.

 23:     Contributed by Xianzhu Tang, LANL

 25:     An interesting feature of this example is that as you refine the grid
 26:     (with a larger -n <n> you cannot force the residual norm as small. This
 27:     appears to be due to "NOISE" in the function, the FormFunctionLocal() cannot
 28:     be computed as accurately with a finer grid.
 29:   ------------------------------------------------------------------------- */

 31:  #include petscda.h
 32:  #include petscsnes.h

 34: /* 
 35:    User-defined application context - contains data needed by the 
 36:    application-provided call-back routines, FormJacobian() and
 37:    FormFunction().
 38: */
 39: typedef struct {
 40:   DA            da;               /* distributed array data structure */
 41:   Vec           psi,r;            /* solution, residual vectors */
 42:   Mat           A,J;              /* Jacobian matrix */
 43:   Vec           coordinates;      /* grid coordinates */
 44:   PassiveReal   psi_axis,psi_bdy;
 45:   PassiveReal   r_min;
 46:   PassiveReal   param;            /* test problem parameter */
 47: } AppCtx;

 49: #define GdGdPsi(r,u)      (((r) < 0.05) ? 0.0 : (user->param*((r)-0.05)*(1.0-(u)*(u))*(1.0-(u)*(u))))
 50: #define CurrentWire(r)    (((r) < .05) ? -3.E2 : 0.0)
 51: #define GdGdPsiPrime(r,u) (((r) < 0.05) ? 0.0 : -4.*(user->param*((r)-0.05)*(1.0-(u)*(u)))*u)

 53: /* 
 54:    User-defined routines
 55: */

 62: int main(int argc,char **argv)
 63: {
 64:   SNES                   snes;                 /* nonlinear solver */
 65:   AppCtx                 user;                 /* user-defined work context */
 66:   PetscInt               its;                  /* iterations for convergence */
 67:   PetscTruth             fd_jacobian = PETSC_FALSE;
 68:   PetscTruth             adicmf_jacobian = PETSC_FALSE;
 69:   PetscInt               grids = 100, dof = 1, stencil_width = 1;
 70:   PetscErrorCode         ierr;
 71:   PetscReal              fnorm;
 72:   MatFDColoring          matfdcoloring = 0;
 73:   ISColoring             iscoloring;

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:      Initialize program
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 79:   PetscInitialize(&argc,&argv,(char *)0,help);
 80: 
 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:      Initialize problem parameters
 83:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   user.psi_axis=0.0;
 86:   PetscOptionsGetReal(PETSC_NULL,"-psi_axis",&user.psi_axis,PETSC_NULL);
 87:   user.psi_bdy=1.0;
 88:   PetscOptionsGetReal(PETSC_NULL,"-psi_bdy",&user.psi_bdy,PETSC_NULL);
 89:   user.r_min=0.0;
 90:   PetscOptionsGetReal(PETSC_NULL,"-r_min",&user.r_min,PETSC_NULL);
 91:   user.param=-1.E1;
 92:   PetscOptionsGetReal(PETSC_NULL,"-param",&user.param,PETSC_NULL);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Create nonlinear solver context
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   SNESCreate(PETSC_COMM_WORLD,&snes);
 98: 
 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100:      Create distributed array (DA) to manage parallel grid and vectors
101:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102:   PetscOptionsGetInt(PETSC_NULL,"-n",&grids,PETSC_NULL);
103:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,grids,dof,stencil_width,PETSC_NULL,&user.da);
104: 
105:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:      Extract global vectors from DA; then duplicate for remaining
107:      vectors that are the same types
108:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109:   DACreateGlobalVector(user.da,&user.psi);
110:   VecDuplicate(user.psi,&user.r);

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:      Set local function evaluation routine
114:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

116:   DASetLocalFunction(user.da,(DALocalFunction1)FormFunctionLocal);

118:   SNESSetFunction(snes,user.r,SNESDAFormFunction,&user);

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:      Create matrix data structure; set Jacobian evaluation routine

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

131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132:   PetscOptionsGetLogical(PETSC_NULL,"-fd_jacobian",&fd_jacobian,0);
133:   /*
134:        Note that fd_jacobian DOES NOT compute the finite difference approximation to 
135:     the ENTIRE Jacobian. Rather it removes the global coupling from the Jacobian and
136:     computes the finite difference approximation only for the "local" coupling.

138:        Thus running with fd_jacobian and not -snes_mf_operator or -adicmf_jacobian
139:     won't converge.
140:   */
141:   if (!fd_jacobian) {
142:     MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,grids,grids,&user.J);
143:     MatSetType(user.J,MATAIJ);
144:     MatSetFromOptions(user.J);
145:     MatSeqAIJSetPreallocation(user.J,5,PETSC_NULL);
146:     MatMPIAIJSetPreallocation(user.J,5,PETSC_NULL,3,PETSC_NULL);
147:     user.A    = user.J;
148:   } else {
149:     DAGetMatrix(user.da,MATAIJ,&user.J);
150:     user.A    = user.J;
151:   }

153:   PetscOptionsGetLogical(PETSC_NULL,"-adicmf_jacobian",&adicmf_jacobian,0);
154: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
155:   if (adicmf_jacobian) {
156:     DASetLocalAdicMFFunction(user.da,admf_FormFunctionLocal);
157:     MatRegisterDAAD();
158:     MatCreateDAAD(user.da,&user.A);
159:     MatDAADSetSNES(user.A,snes);
160:     MatDAADSetCtx(user.A,&user);
161:   }
162: #endif

164:   if (fd_jacobian) {
165:     DAGetColoring(user.da,IS_COLORING_LOCAL,&iscoloring);
166:     MatFDColoringCreate(user.J,iscoloring,&matfdcoloring);
167:     ISColoringDestroy(iscoloring);
168:     MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESDAFormFunction,&user);
169:     MatFDColoringSetFromOptions(matfdcoloring);
170:     SNESSetJacobian(snes,user.A,user.J,SNESDefaultComputeJacobianColor,matfdcoloring);
171:   } else {
172:     SNESSetJacobian(snes,user.A,user.J,FormJacobian,&user);
173:   }

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176:      Customize nonlinear solver; set runtime options
177:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178:   SNESSetFromOptions(snes);

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

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190:      Solve nonlinear system
191:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192:   SNESSolve(snes,user.psi);
193:   SNESGetIterationNumber(snes,&its);

195:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196:      Explicitly check norm of the residual of the solution
197:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198:   SNESDAFormFunction(snes,user.psi,user.r,(void*)&user);
199:   VecNorm(user.r,NORM_MAX,&fnorm);
200:   PetscPrintf(PETSC_COMM_WORLD,"Number of Newton iterations = %D fnorm %g\n",its,fnorm);

202:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
203:      Output the solution vector
204:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205:   {
206:     PetscViewer view_out;
207:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"psi.binary",PETSC_FILE_CREATE,&view_out);
208:     VecView(user.psi,view_out);
209:     PetscViewerDestroy(view_out);
210:     PetscViewerASCIIOpen(PETSC_COMM_WORLD,"psi.out",&view_out);
211:     VecView(user.psi,view_out);
212:     PetscViewerDestroy(view_out);
213:   }

215:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216:      Free work space.  All PETSc objects should be destroyed when they
217:      are no longer needed.
218:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

220:   if (user.A != user.J) {
221:     MatDestroy(user.A);
222:   }
223:   MatDestroy(user.J);
224:   if (matfdcoloring) {
225:     MatFDColoringDestroy(matfdcoloring);
226:   }
227:   VecDestroy(user.psi);
228:   VecDestroy(user.r);
229:   SNESDestroy(snes);
230:   DADestroy(user.da);
231:   PetscFinalize();

233:   return(0);
234: }
235: /* ------------------------------------------------------------------- */
238: /* 
239:    FormInitialGuess - Forms initial approximation.

241:    Input Parameters:
242:    user - user-defined application context
243:    X - vector

245:    Output Parameter:
246:    X - vector
247:  */
248: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
249: {
251:   PetscInt       i,Mx,xs,xm;
252:   PetscScalar    *x;

255:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
256:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

258:   /*
259:      Get a pointer to vector data.
260:        - For default PETSc vectors, VecGetArray() returns a pointer to
261:          the data array.  Otherwise, the routine is implementation dependent.
262:        - You MUST call VecRestoreArray() when you no longer need access to
263:          the array.
264:   */
265:   DAVecGetArray(user->da,X,&x);

267:   /*
268:      Get local grid boundaries (for 2-dimensional DA):
269:        xs, ys   - starting grid indices (no ghost points)
270:        xm, ym   - widths of local grid (no ghost points)

272:   */
273:   DAGetCorners(user->da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);

275:   /*
276:      Compute initial guess over the locally owned part of the grid
277:   */
278:   for (i=xs; i<xs+xm; i++) {
279:     x[i] = user->psi_axis + i*(user->psi_bdy - user->psi_axis)/(PetscReal)(Mx-1);
280:   }

282:   /*
283:      Restore vector
284:   */
285:   DAVecRestoreArray(user->da,X,&x);

287:   /* 
288:      Check to see if we can import an initial guess from disk
289:   */
290:   {
291:     char         filename[PETSC_MAX_PATH_LEN];
292:     PetscTruth   flg;
293:     PetscViewer  view_in;
294:     PetscReal    fnorm;
295:     Vec          Y;
296:     PetscOptionsGetString(PETSC_NULL,"-initialGuess",filename,256,&flg);
297:     if (flg) {
298:       PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,PETSC_FILE_RDONLY,&view_in);
299:       VecLoad(view_in,PETSC_NULL,&Y);
300:       PetscViewerDestroy(view_in);
301:       VecMax(Y,PETSC_NULL,&user->psi_bdy);
302:       SNESDAFormFunction(PETSC_NULL,Y,user->r,(void*)user);
303:       VecNorm(user->r,NORM_2,&fnorm);
304:       PetscPrintf(PETSC_COMM_WORLD,"In initial guess: psi_bdy = %f, fnorm = %g.\n",user->psi_bdy,fnorm);
305:       VecCopy(Y,X);
306:       VecDestroy(Y);
307:     }
308:   }

310:   return(0);
311: }
312: /* ------------------------------------------------------------------- */
315: /* 
316:    FormFunctionLocal - Evaluates nonlinear function, F(x).
317:    
318:    Process adiC(36): FormFunctionLocal
319:    
320: */
321: PetscErrorCode FormFunctionLocal(DALocalInfo *info,PetscScalar *x,PetscScalar *f,AppCtx *user)
322: {
324:   PetscInt       i,xint,xend;
325:   PetscReal      hx,dhx,r;
326:   PetscScalar    u,uxx,min = 100000.0,max = -10000.0;
327:   PetscScalar    psi_0=0.0, psi_a=1.0;
328: 
330:   for (i=info->xs; i<info->xs + info->xm; i++) {
331:     if (x[i] > max) max = x[i];
332:     if (x[i] < min) min = x[i];
333:   }
334:   PetscGlobalMax(&max,&psi_a,PETSC_COMM_WORLD);
335:   PetscGlobalMin(&min,&psi_0,PETSC_COMM_WORLD);

337:   hx     = 1.0/(PetscReal)(info->mx-1);  dhx    = 1.0/hx;
338: 
339:   /*
340:     Compute function over the locally owned part of the grid
341:   */
342:   if (info->xs == 0) {
343:     xint = info->xs + 1; f[0] = (4.*x[1] - x[2] - 3.*x[0])*dhx; /* f[0] = x[0] - user->psi_axis; */
344:   }
345:   else {
346:     xint = info->xs;
347:   }
348:   if ((info->xs+info->xm) == info->mx) {
349:     xend = info->mx - 1; f[info->mx-1] = -(x[info->mx-1] - user->psi_bdy)*dhx;
350:   }
351:   else {
352:     xend = info->xs + info->xm;
353:   }
354: 
355:   for (i=xint; i<xend; i++) {
356:     r       = i*hx + user->r_min;   /* r coordinate */
357:     u       = (x[i]-psi_0)/(psi_a - psi_0);
358:     uxx     = ((r+0.5*hx)*(x[i+1]-x[i]) - (r-0.5*hx)*(x[i]-x[i-1]))*dhx; /* r*nabla^2\psi */
359:     f[i] = uxx  + r*GdGdPsi(r,u)*hx  + r*CurrentWire(r)*hx ;
360:   }

362:   PetscLogFlops(11*info->ym*info->xm);
363:   return(0);
364: }
365: /* ------------------------------------------------------------------- */
368: /*
369:    FormJacobian - Evaluates Jacobian matrix.

371:    Input Parameters:
372: .  snes - the SNES context
373: .  x - input vector
374: .  ptr - optional user-defined context, as set by SNESSetJacobian()

376:    Output Parameters:
377: .  A - Jacobian matrix
378: .  B - optionally different preconditioning matrix
379: .  flag - flag indicating matrix structure

381: */
382: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
383: {
384:   AppCtx         *user = (AppCtx*)ptr;  /* user-defined application context */
385:   Mat            jac = *B;                /* Jacobian matrix */
386:   Vec            localX;
388:   PetscInt       col[6],row,i,xs,xm,Mx,xint,xend,imin, imax;
389:   PetscScalar    v[6],hx,dhx,*x;
390:   PetscReal      r, u, psi_0=0.0, psi_a=1.0;
391:   PetscTruth     assembled;

394:   MatAssembled(*B,&assembled);
395:   if (assembled) {
396:     MatZeroEntries(*B);
397:   }

399:   DAGetLocalVector(user->da,&localX);
400:   DAGetInfo(user->da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
401:                    PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

403:   hx     = 1.0/(PetscReal)(Mx-1);  dhx    = 1.0/hx;

405:   imin = 0; imax = Mx-1;
406:   VecMin(X,&imin,&psi_0);
407:   VecMax(X,&imax,&psi_a);
408:   PetscPrintf(PETSC_COMM_WORLD,"psi_0(%D)=%g, psi_a(%D)=%g.\n",imin,psi_0,imax,psi_a);

410:   /*
411:      Scatter ghost points to local vector, using the 2-step process
412:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
413:      By placing code between these two statements, computations can be
414:      done while messages are in transition.
415:   */
416:   DAGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
417:   DAGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);

419:   /*
420:      Get pointer to vector data
421:   */
422:   DAVecGetArray(user->da,localX,&x);

424:   /*
425:      Get local grid boundaries
426:   */
427:   DAGetCorners(user->da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);

429:   /* 
430:      Compute entries for the locally owned part of the Jacobian.
431:       - Currently, all PETSc parallel matrix formats are partitioned by
432:         contiguous chunks of rows across the processors. 
433:       - Each processor needs to insert only elements that it owns
434:         locally (but any non-local elements will be sent to the
435:         appropriate processor during matrix assembly). 
436:       - Here, we set all entries for a particular row at once.
437:       - We can set matrix entries either using either
438:         MatSetValuesLocal() or MatSetValues(), as discussed above.
439:   */
440:   if (xs == 0) {
441:     xint = xs + 1; /* f[0] = 4.*x[1] - x[2] - 3.*x[0]; */
442:     row  = 0;     /* first row */
443:     v[0] = -3.0*dhx;                                              col[0]=row;
444:     v[1] = 4.0*dhx;                                               col[1]=row+1;
445:     v[2] = -1.0*dhx;                                              col[2]=row+2;
446:     MatSetValues(jac,1,&row,3,col,v,ADD_VALUES);
447:   }
448:   else {
449:     xint = xs;
450:   }
451:   if ((xs+xm) == Mx) {
452:     xend = Mx - 1;   /* f[Mx-1] = x[Mx-1] - user->psi_bdy; */
453:     row  = Mx - 1;  /* last row */
454:     v[0] = -1.0*dhx;
455:     MatSetValue(jac,row,row,v[0],ADD_VALUES);
456:   }
457:   else {
458:     xend = xs + xm;
459:   }

461:   for (i=xint; i<xend; i++) {
462:     r       = i*hx + user->r_min;   /* r coordinate */
463:     u       = (x[i]-psi_0)/(psi_a - psi_0);
464:     /* uxx     = ((r+0.5*hx)*(x[i+1]-x[i]) - (r-0.5*hx)*(x[i]-x[i-1]))*dhx*dhx; */ /* r*nabla^2\psi */
465:     row  = i;
466:     v[0] = (r-0.5*hx)*dhx;                                                              col[0] = i-1;
467:     v[1] = -2.*r*dhx + hx*r*GdGdPsiPrime(r,u)/(psi_a - psi_0);                          col[1] = i;
468:     v[2] = (r+0.5*hx)*dhx;                                                              col[2] = i+1;
469:     v[3] = hx*r*GdGdPsiPrime(r,u)*(x[i] - psi_a)/((psi_a - psi_0)*(psi_a - psi_0));     col[3] = imin;
470:     v[4] = hx*r*GdGdPsiPrime(r,u)*(psi_0 - x[i])/((psi_a - psi_0)*(psi_a - psi_0));     col[4] = imax;
471:     MatSetValues(jac,1,&row,5,col,v,ADD_VALUES);
472:   }
473: 
474:   DAVecRestoreArray(user->da,localX,&x);
475:   DARestoreLocalVector(user->da,&localX);

477:   /* 
478:      Assemble matrix, using the 2-step process:
479:        MatAssemblyBegin(), MatAssemblyEnd().
480:   */
481:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
482:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);

484:   /*
485:      Normally since the matrix has already been assembled above; this
486:      would do nothing. But in the matrix free mode -snes_mf_operator
487:      this tells the "matrix-free" matrix that a new linear system solve
488:      is about to be done.
489:   */
490:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
491:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

493:   /*
494:      Set flag to indicate that the Jacobian matrix retains an identical
495:      nonzero structure throughout all nonlinear iterations (although the
496:      values of the entries change). Thus, we can save some work in setting
497:      up the preconditioner (e.g., no need to redo symbolic factorization for
498:      ILU/ICC preconditioners).
499:       - If the nonzero structure of the matrix is different during
500:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
501:         must be used instead.  If you are unsure whether the matrix
502:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
503:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
504:         believes your assertion and does not check the structure
505:         of the matrix.  If you erroneously claim that the structure
506:         is the same when it actually is not, the new preconditioner
507:         will not function correctly.  Thus, use this optimization
508:         feature with caution!
509:   */
510:   *flag = SAME_NONZERO_PATTERN;

512:   /*
513:      Tell the matrix we will never add a new nonzero location to the
514:      matrix. If we do, it will generate an error.
515:   */

517:   return(0);
518: }