Actual source code: ex29.c

  2: /*#define HAVE_DA_HDF*/

  4: /* solve the equations for the perturbed magnetic field only */
  5: #define EQ 

  7: /* turning this on causes instability?!? */
  8: /* #define UPWINDING  */

 10: static char help[] = "Hall MHD with in two dimensions with time stepping and multigrid.\n\n\
 11: -options_file ex29.options\n\
 12: other PETSc options\n\
 13: -resistivity <eta>\n\
 14: -viscosity <nu>\n\
 15: -skin_depth <d_e>\n\
 16: -larmor_radius <rho_s>\n\
 17: -contours : draw contour plots of solution\n\n";

 19: /*T
 20:    Concepts: SNES^solving a system of nonlinear equations (parallel multicomponent example);
 21:    Concepts: DA^using distributed arrays;
 22:    Concepts: multicomponent
 23:    Processors: n
 24: T*/

 26: /* ------------------------------------------------------------------------

 28:     We thank Kai Germaschewski for contributing the FormFunctionLocal()

 30:   ------------------------------------------------------------------------- */

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

 48: #ifdef HAVE_DA_HDF
 49: PetscInt DAVecHDFOutput(DA,Vec,char*);
 50: #endif

 52: #define D_x(x,m,i,j)  (p5 * (x[(j)][(i)+1].m - x[(j)][(i)-1].m) * dhx)
 53: #define D_xm(x,m,i,j) ((x[(j)][(i)].m - x[(j)][(i)-1].m) * dhx)
 54: #define D_xp(x,m,i,j) ((x[(j)][(i+1)].m - x[(j)][(i)].m) * dhx)
 55: #define D_y(x,m,i,j)  (p5 * (x[(j)+1][(i)].m - x[(j)-1][(i)].m) * dhy)
 56: #define D_ym(x,m,i,j) ((x[(j)][(i)].m - x[(j)-1][(i)].m) * dhy)
 57: #define D_yp(x,m,i,j) ((x[(j)+1][(i)].m - x[(j)][(i)].m) * dhy)
 58: #define D_xx(x,m,i,j) ((x[(j)][(i)+1].m - two*x[(j)][(i)].m + x[(j)][(i)-1].m) * hydhx * dhxdhy)
 59: #define D_yy(x,m,i,j) ((x[(j)+1][(i)].m - two*x[(j)][(i)].m + x[(j)-1][(i)].m) * hxdhy * dhxdhy)
 60: #define Lapl(x,m,i,j) (D_xx(x,m,i,j) + D_yy(x,m,i,j))
 61: #define lx            (2.*3.1415926)
 62: #define ly            (4.*3.1415926)
 63: #define sqr(a)        ((a)*(a))

 65: /* 
 66:    User-defined routines and data structures
 67: */

 69: typedef struct {
 70:   PetscReal      fnorm_ini,dt_ini;
 71:   PetscReal      dt,dt_out;
 72:   PetscReal      ptime;
 73:   PetscReal      max_time;
 74:   PetscReal      fnorm_ratio;
 75:   PetscInt       ires,itstep;
 76:   PetscInt       max_steps,print_freq;
 77:   PetscReal      t;
 78:   PetscReal      fnorm;

 80:   PetscTruth     ts_monitor;           /* print information about each time step */
 81:   PetscReal      dump_time;            /* time to dump solution to a file */
 82:   PetscViewer    socketviewer;         /* socket to dump solution at each timestep for visualization */
 83: } TstepCtx;

 85: typedef struct {
 86:   PetscScalar phi,psi,U,F;
 87: } Field;

 89: typedef struct {
 90:   PassiveScalar phi,psi,U,F;
 91: } PassiveField;

 93: typedef struct {
 94:   PetscInt     mglevels;
 95:   PetscInt     cycles;           /* number of time steps for integration */
 96:   PassiveReal  nu,eta,d_e,rho_s; /* physical parameters */
 97:   PetscTruth   draw_contours;    /* flag - 1 indicates drawing contours */
 98:   PetscTruth   second_order;
 99:   PetscTruth   PreLoading;
100: } Parameters;

102: typedef struct {
103:   Vec          Xoldold,Xold;
104:   TstepCtx     *tsCtx;
105:   Parameters    *param;
106: } AppCtx;


120: int main(int argc,char **argv)
121: {
123:   DMMG           *dmmg;       /* multilevel grid structure */
124:   AppCtx         *user;       /* user-defined work context (one for each level) */
125:   TstepCtx       tsCtx;       /* time-step parameters (one total) */
126:   Parameters      param;       /* physical parameters (one total) */
127:   PetscInt       i,m,n,mx,my;
128:   DA             da;
129:   PetscReal      dt_ratio;
130:   PetscInt       dfill[16] = {1,0,1,0,
131:                               0,1,0,1,
132:                               0,0,1,1,
133:                               0,1,1,1};
134:   PetscInt       ofill[16] = {1,0,0,0,
135:                               0,1,0,0,
136:                               1,1,1,1,
137:                               1,1,1,1};

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


142:   PreLoadBegin(PETSC_TRUE,"SetUp");

144:   param.PreLoading = PreLoading;
145:     DMMGCreate(PETSC_COMM_WORLD,1,&user,&dmmg);
146:     param.mglevels = DMMGGetLevels(dmmg);

148:     /*
149:       Create distributed array multigrid object (DMMG) to manage
150:       parallel grid and vectors for principal unknowns (x) and
151:       governing residuals (f)
152:     */
153:     DACreate2d(PETSC_COMM_WORLD, DA_XYPERIODIC, DA_STENCIL_STAR, -5, -5,
154:                       PETSC_DECIDE, PETSC_DECIDE, 4, 1, 0, 0, &da);

156:     /* overwrite the default sparsity pattern toone specific for
157:        this code's nonzero structure */
158:     DASetBlockFills(da,dfill,ofill);

160:     DMMGSetDM(dmmg,(DM)da);
161:     DADestroy(da);

163:     /* default physical parameters */
164:     param.nu    = 0;
165:     param.eta   = 1e-3;
166:     param.d_e   = 0.2;
167:     param.rho_s = 0;
168:     param.second_order = PETSC_FALSE;

170:     PetscOptionsGetReal(PETSC_NULL, "-viscosity", &param.nu,PETSC_NULL);

172:     PetscOptionsGetReal(PETSC_NULL, "-resistivity", &param.eta,PETSC_NULL);

174:     PetscOptionsGetReal(PETSC_NULL, "-skin_depth", &param.d_e,PETSC_NULL);

176:     PetscOptionsGetReal(PETSC_NULL, "-larmor_radius", &param.rho_s,PETSC_NULL);

178:     PetscOptionsHasName(PETSC_NULL, "-contours", &param.draw_contours);

180:     PetscOptionsHasName(PETSC_NULL, "-second_order", &param.second_order);

182:     DASetFieldName(DMMGGetDA(dmmg), 0, "phi");

184:     DASetFieldName(DMMGGetDA(dmmg), 1, "psi");

186:     DASetFieldName(DMMGGetDA(dmmg), 2, "U");

188:     DASetFieldName(DMMGGetDA(dmmg), 3, "F");

190:     /*======================================================================*/
191:     /* Initialize stuff related to time stepping */
192:     /*======================================================================*/

194:     DAGetInfo(DMMGGetDA(dmmg),0,&mx,&my,0,0,0,0,0,0,0,0);

196:     tsCtx.fnorm_ini   = 0;
197:     tsCtx.max_steps   = 1000000;
198:     tsCtx.max_time    = 1.0e+12;
199:     /* use for dt = min(dx,dy); multiplied by dt_ratio below */
200:     tsCtx.dt          = PetscMin(lx/mx,ly/my);
201:     tsCtx.fnorm_ratio = 1.0e+10;
202:     tsCtx.t           = 0;
203:     tsCtx.dt_out      = 10;
204:     tsCtx.print_freq  = tsCtx.max_steps;
205:     tsCtx.ts_monitor  = PETSC_FALSE;
206:     tsCtx.dump_time   = -1.0;

208:     PetscOptionsGetInt(PETSC_NULL, "-max_st", &tsCtx.max_steps,PETSC_NULL);
209:     PetscOptionsGetReal(PETSC_NULL, "-ts_rtol", &tsCtx.fnorm_ratio,PETSC_NULL);
210:     PetscOptionsGetInt(PETSC_NULL, "-print_freq", &tsCtx.print_freq,PETSC_NULL);

212:     dt_ratio = 1.0;
213:     PetscOptionsGetReal(PETSC_NULL, "-dt_ratio", &dt_ratio,PETSC_NULL);
214:     tsCtx.dt *= dt_ratio;

216:     PetscOptionsHasName(PETSC_NULL, "-ts_monitor", &tsCtx.ts_monitor);
217:     PetscOptionsGetReal(PETSC_NULL, "-dump", &tsCtx.dump_time,PETSC_NULL);


220:     tsCtx.socketviewer = 0;
221: #if defined(PETSC_USE_SOCKET_VIEWER)
222:     {
223:       PetscTruth flg;
224:       PetscOptionsHasName(PETSC_NULL, "-socket_viewer", &flg);
225:       if (flg && !PreLoading) {
226:         tsCtx.ts_monitor = PETSC_TRUE;
227:         PetscViewerSocketOpen(PETSC_COMM_WORLD,0,PETSC_DECIDE,&tsCtx.socketviewer);
228:       }
229:     }
230: #endif

232:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233:        Create user context, set problem data, create vector data structures.
234:        Also, compute the initial guess.
235:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236:     /* create application context for each level */
237:     PetscMalloc(param.mglevels*sizeof(AppCtx),&user);
238:     for (i=0; i<param.mglevels; i++) {
239:       /* create work vectors to hold previous time-step solution and
240:          function value */
241:       VecDuplicate(dmmg[i]->x, &user[i].Xoldold);
242:       VecDuplicate(dmmg[i]->x, &user[i].Xold);
243:       user[i].tsCtx = &tsCtx;
244:       user[i].param = &param;
245:       DMMGSetUser(dmmg,i,&user[i]);
246:     }
247: 
248:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249:        Create nonlinear solver context
250:        
251:        Process adiC(20):  AddTSTermLocal AddTSTermLocal2 FormFunctionLocal FormFunctionLocali
252:        Process blockadiC(4):  FormFunctionLocali4
253:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
254:     DMMGSetSNESLocal(dmmg, FormFunctionLocal, 0,ad_FormFunctionLocal, admf_FormFunctionLocal);
255:     DMMGSetFromOptions(dmmg);
256:     DMMGSetSNESLocali(dmmg,FormFunctionLocali,0,admf_FormFunctionLocali);
257:     DMMGSetSNESLocalib(dmmg,FormFunctionLocali4,0,admfb_FormFunctionLocali4);
258: 
259:    /* attach nullspace to each level of the preconditioner */
260:     DMMGSetNullSpace(dmmg,PETSC_FALSE,1,CreateNullSpace);
261: 
262:     PetscPrintf(PETSC_COMM_WORLD, "finish setupNull!");

264:     if (PreLoading) {
265:       PetscPrintf(PETSC_COMM_WORLD, "# viscosity = %G, resistivity = %G, "
266:                          "skin_depth # = %G, larmor_radius # = %G\n",
267:                          param.nu, param.eta, param.d_e, param.rho_s);
268:       DAGetInfo(DMMGGetDA(dmmg),0,&m,&n,0,0,0,0,0,0,0,0);
269:       PetscPrintf(PETSC_COMM_WORLD,"Problem size %D by %D\n",m,n);
270:       PetscPrintf(PETSC_COMM_WORLD,"dx %G dy %G dt %G ratio dt/min(dx,dy) %G\n",lx/mx,ly/my,tsCtx.dt,dt_ratio);
271:     }

273: 
274: 
275:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276:        Solve the nonlinear system
277:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278: 
279:     PreLoadStage("Solve");

281:     if (param.draw_contours) {
282:       VecView(((AppCtx*)DMMGGetUser(dmmg,param.mglevels-1))->Xold,PETSC_VIEWER_DRAW_WORLD);
283:     }

285:     Update(dmmg);

287:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288:        Free work space.  All PETSc objects should be destroyed when they
289:        are no longer needed.
290:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
291: 
292:     for (i=0; i<param.mglevels; i++) {
293:       VecDestroy(user[i].Xoldold);
294:       VecDestroy(user[i].Xold);
295:     }
296:     PetscFree(user);
297:     DMMGDestroy(dmmg);

299:     PreLoadEnd();
300: 
301:   PetscFinalize();
302:   return 0;
303: }

305: /* ------------------------------------------------------------------- */
308: /* ------------------------------------------------------------------- */
309: PetscErrorCode Gnuplot(DA da, Vec X, double mtime)
310: {
311:   PetscInt       i,j,xs,ys,xm,ym;
312:   PetscInt       xints,xinte,yints,yinte;
314:   Field          **x;
315:   FILE           *f;
316:   char           fname[PETSC_MAX_PATH_LEN];
317:   PetscMPIInt    rank;

319:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

321:   sprintf(fname, "out-%g-%d.dat", mtime, rank);

323:   f = fopen(fname, "w");

325:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

327:   DAVecGetArray(da,X,&x);

329:   xints = xs; xinte = xs+xm; yints = ys; yinte = ys+ym;

331:   for (j=yints; j<yinte; j++) {
332:     for (i=xints; i<xinte; i++) {
333:       PetscFPrintf(PETSC_COMM_WORLD, f,
334:                           "%D %D %G %G %G %G %G %G\n",
335:                           i, j, 0.0, 0.0,
336:                           PetscAbsScalar(x[j][i].U), PetscAbsScalar(x[j][i].F),
337:                           PetscAbsScalar(x[j][i].phi), PetscAbsScalar(x[j][i].psi));
338:     }
339:     PetscFPrintf(PETSC_COMM_WORLD,f, "\n");
340:   }
341:   DAVecRestoreArray(da,X,&x);
342:   fclose(f);
343:   return 0;
344: }

346: /* ------------------------------------------------------------------- */
349: /* ------------------------------------------------------------------- */
350: PetscErrorCode Initialize(DMMG *dmmg)
351: {
352:   AppCtx         *appCtx = (AppCtx*)DMMGGetUser(dmmg,0);
353:   Parameters      *param = appCtx->param;
354:   DA             da;
356:   PetscInt       i,j,mx,my,xs,ys,xm,ym;
357:   PetscReal      two = 2.0,one = 1.0;
358:   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
359:   PetscReal      d_e,rho_s,de2,xx,yy;
360:   Field          **x, **localx;
361:   Vec            localX;
362:   PetscTruth     flg;

365:   PetscOptionsHasName(0,"-restart",&flg);
366:   if (flg) {
367:     PetscViewer viewer;
368:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"binaryoutput",FILE_MODE_READ,&viewer);
369:     VecLoadIntoVector(viewer,dmmg[param->mglevels-1]->x);
370:     PetscViewerDestroy(viewer);
371:     return(0);
372:   }

374:   d_e   = param->d_e;
375:   rho_s = param->rho_s;
376:   de2   = sqr(param->d_e);

378:   da   = (DA)(dmmg[param->mglevels-1]->dm);
379:   DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);

381:   dhx   = mx/lx;              dhy = my/ly;
382:   hx    = one/dhx;             hy = one/dhy;
383:   hxdhy = hx*dhy;           hydhx = hy*dhx;
384:   hxhy  = hx*hy;           dhxdhy = dhx*dhy;

386:   /*
387:      Get local grid boundaries (for 2-dimensional DA):
388:        xs, ys   - starting grid indices (no ghost points)
389:        xm, ym   - widths of local grid (no ghost points)
390:   */
391:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);

393:   DAGetLocalVector(da,&localX);
394:   /*
395:      Get a pointer to vector data.
396:        - For default PETSc vectors, VecGetArray() returns a pointer to
397:          the data array.  Otherwise, the routine is implementation dependent.
398:        - You MUST call VecRestoreArray() when you no longer need access to
399:          the array.
400:   */
401:   DAVecGetArray(da,dmmg[param->mglevels-1]->x,&x);
402:   DAVecGetArray(da,localX,&localx);

404:   /*
405:      Compute initial solution over the locally owned part of the grid
406:   */
407: #if defined (PETSC_HAVE_ERF)
408:   {
409:     PetscReal eps = lx/ly;
410:     PetscReal pert = 1e-4;
411:     PetscReal k = 1.*eps;
412:     PetscReal gam;

414:     if (d_e < rho_s) d_e = rho_s;
415:     gam = k * d_e;
416:     for (j=ys-1; j<ys+ym+1; j++) {
417:       yy = j * hy;
418:       for (i=xs-1; i<xs+xm+1; i++) {
419:         xx = i * hx;
420:         if (xx < -3.1416926/2) {
421:           localx[j][i].phi = pert * gam / k * erf((xx + 3.1415926) / (sqrt(2.0) * d_e)) * (-sin(k*yy));
422:         } else if (xx < 3.1415926/2) {
423:           localx[j][i].phi = - pert * gam / k * erf(xx / (sqrt(2.0) * d_e)) * (-sin(k*yy));
424:         } else if (xx < 3*3.1415926/2){
425:           localx[j][i].phi = pert * gam / k * erf((xx - 3.1415926) / (sqrt(2.0) * d_e)) * (-sin(k*yy));
426:         } else {
427:           localx[j][i].phi = - pert * gam / k * erf((xx - 2.*3.1415926) / (sqrt(2.0) * d_e)) * (-sin(k*yy));
428:         }
429: #ifdef EQ
430:         localx[j][i].psi = 0.;
431: #else
432:         localx[j][i].psi = cos(xx);
433: #endif
434:       }
435:     }

437:     for (j=ys; j<ys+ym; j++) {
438:       for (i=xs; i<xs+xm; i++) {
439:         x[j][i].U   = Lapl(localx,phi,i,j);
440:         x[j][i].F   = localx[j][i].psi - de2 * Lapl(localx,psi,i,j);
441:         x[j][i].phi = localx[j][i].phi;
442:         x[j][i].psi = localx[j][i].psi;
443:       }
444:     }
445:   }
446: #else
447:   SETERRQ(1,"erf() not available on this machine");
448: #endif

450:   /*
451:      Restore vector
452:   */
453:   DAVecRestoreArray(da,dmmg[param->mglevels-1]->x,&x);
454: 
455:   DAVecRestoreArray(da,localX,&localx);
456: 
457:   DARestoreLocalVector(da,&localX);
458: 

460:   return(0);
461: }

465: PetscErrorCode ComputeMaxima(DA da, Vec X, PetscReal t)
466: {
468:   PetscInt       i,j,mx,my,xs,ys,xm,ym;
469:   PetscInt       xints,xinte,yints,yinte;
470:   Field          **x;
471:   double         norm[4] = {0,0,0,0};
472:   double         gnorm[4];
473:   MPI_Comm       comm;


477:   DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);

479:   DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
480: 

482:   xints = xs; xinte = xs+xm; yints = ys; yinte = ys+ym;

484:   DAVecGetArray(da, X, &x);

486:   for (j=yints; j<yinte; j++) {
487:     for (i=xints; i<xinte; i++) {
488:       norm[0] = PetscMax(norm[0],PetscAbsScalar(x[j][i].U));
489:       norm[1] = PetscMax(norm[1],PetscAbsScalar(x[j][i].F));
490:       norm[2] = PetscMax(norm[2],PetscAbsScalar(x[j][i].phi));
491:       norm[3] = PetscMax(norm[3],PetscAbsScalar(x[j][i].psi));
492:     }
493:   }

495:   DAVecRestoreArray(da,X,&x);

497:   PetscObjectGetComm((PetscObject)da, &comm);
498:   MPI_Allreduce(norm, gnorm, 4, MPI_DOUBLE, MPI_MAX, comm);

500:   PetscFPrintf(PETSC_COMM_WORLD, stderr,"%G\t%G\t%G\t%G\t%G\n",t, gnorm[0], gnorm[1], gnorm[2], gnorm[3]);

502:   return(0);
503: }

507: PetscErrorCode FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
508: {
509:   AppCtx         *user = (AppCtx*)ptr;
510:   TstepCtx       *tsCtx = user->tsCtx;
511:   Parameters      *param = user->param;
513:   PetscInt       xints,xinte,yints,yinte,i,j;
514:   PassiveReal    hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
515:   PassiveReal    de2,rhos2,nu,eta,dde2;
516:   PassiveReal    two = 2.0,one = 1.0,p5 = 0.5;
517:   PassiveReal    F_eq_x,By_eq;
518:   PetscScalar    xx;
519:   PetscScalar    vx,vy,avx,avy,vxp,vxm,vyp,vym;
520:   PetscScalar    Bx,By,aBx,aBy,Bxp,Bxm,Byp,Bym;

523:   de2     = sqr(param->d_e);
524:   rhos2   = sqr(param->rho_s);
525:   nu      = param->nu;
526:   eta     = param->eta;
527:   dde2    = one/de2;

529:   /* 
530:      Define mesh intervals ratios for uniform grid.
531:      [Note: FD formulae below are normalized by multiplying through by
532:      local volume element to obtain coefficients O(1) in two dimensions.]
533:   */
534:   dhx   = info->mx/lx;        dhy   = info->my/ly;
535:   hx    = one/dhx;             hy   = one/dhy;
536:   hxdhy = hx*dhy;           hydhx   = hy*dhx;
537:   hxhy  = hx*hy;             dhxdhy = dhx*dhy;

539:   xints = info->xs; xinte = info->xs+info->xm;
540:   yints = info->ys; yinte = info->ys+info->ym;

542:   /* Compute over the interior points */
543:   for (j=yints; j<yinte; j++) {
544:     for (i=xints; i<xinte; i++) {
545: #ifdef EQ
546:       xx = i * hx;
547:       F_eq_x = - (1. + de2) * sin(PetscAbsScalar(xx));
548:       By_eq = sin(PetscAbsScalar(xx));
549: #else
550:       F_eq_x = 0.;
551:       By_eq = 0.;
552: #endif

554:       /*
555:        * convective coefficients for upwinding
556:        */

558:       vx  = - D_y(x,phi,i,j);
559:       vy  =   D_x(x,phi,i,j);
560:       avx = PetscAbsScalar(vx); vxp = p5*(vx+avx); vxm = p5*(vx-avx);
561:       avy = PetscAbsScalar(vy); vyp = p5*(vy+avy); vym = p5*(vy-avy);
562: #ifndef UPWINDING
563:       vxp = vxm = p5*vx;
564:       vyp = vym = p5*vy;
565: #endif

567:       Bx  =   D_y(x,psi,i,j);
568:       By  = - D_x(x,psi,i,j) + By_eq;
569:       aBx = PetscAbsScalar(Bx); Bxp = p5*(Bx+aBx); Bxm = p5*(Bx-aBx);
570:       aBy = PetscAbsScalar(By); Byp = p5*(By+aBy); Bym = p5*(By-aBy);
571: #ifndef UPWINDING
572:       Bxp = Bxm = p5*Bx;
573:       Byp = Bym = p5*By;
574: #endif

576:       /* Lap(phi) - U */
577:       f[j][i].phi = (Lapl(x,phi,i,j) - x[j][i].U) * hxhy;

579:       /* psi - d_e^2 * Lap(psi) - F */
580:       f[j][i].psi = (x[j][i].psi - de2 * Lapl(x,psi,i,j) - x[j][i].F) * hxhy;

582:       /* vx * U_x + vy * U_y - (B_x * F_x + B_y * F_y) / d_e^2
583:          - nu Lap(U) */
584:       f[j][i].U  =
585:         ((vxp * (D_xm(x,U,i,j)) + vxm * (D_xp(x,U,i,j)) +
586:           vyp * (D_ym(x,U,i,j)) + vym * (D_yp(x,U,i,j))) -
587:          (Bxp * (D_xm(x,F,i,j) + F_eq_x) + Bxm * (D_xp(x,F,i,j) + F_eq_x) +
588:           Byp * (D_ym(x,F,i,j)) + Bym * (D_yp(x,F,i,j))) * dde2 -
589:          nu * Lapl(x,U,i,j)) * hxhy;
590: 
591:       /* vx * F_x + vy * F_y - rho_s^2 * (B_x * U_x + B_y * U_y)
592:          - eta * Lap(psi) */
593:       f[j][i].F  =
594:         ((vxp * (D_xm(x,F,i,j) + F_eq_x) + vxm * (D_xp(x,F,i,j) + F_eq_x) +
595:           vyp * (D_ym(x,F,i,j)) + vym * (D_yp(x,F,i,j))) -
596:          (Bxp * (D_xm(x,U,i,j)) + Bxm * (D_xp(x,U,i,j)) +
597:           Byp * (D_ym(x,U,i,j)) + Bym * (D_yp(x,U,i,j))) * rhos2 -
598:          eta * Lapl(x,psi,i,j)) * hxhy;
599:     }
600:   }

602:   /* Add time step contribution */
603:   if (tsCtx->ires) {
604:     if ((param->second_order) && (tsCtx->itstep > 0)){
605:       AddTSTermLocal2(info,x,f,user);
606:     } else {
607:       AddTSTermLocal(info,x,f,user);
608:     }
609:   }

611:   /*
612:      Flop count (multiply-adds are counted as 2 operations)
613:   */
614:   /*  PetscLogFlops(84*info->ym*info->xm); FIXME */
615:   return(0);
616: }

618: /*---------------------------------------------------------------------*/
621: PetscErrorCode Update(DMMG *dmmg)
622: /*---------------------------------------------------------------------*/
623: {
624:   AppCtx          *user = (AppCtx *) DMMGGetUser(dmmg,0);
625:   TstepCtx        *tsCtx = user->tsCtx;
626:   Parameters       *param = user->param;
627:   SNES            snes;
628:   PetscErrorCode  ierr;
629:   PetscInt        its,lits,i;
630:   PetscInt        max_steps;
631:   PetscInt        nfailsCum = 0,nfails = 0;
632:   static PetscInt ic_out;
633:   PetscTruth      ts_monitor = (tsCtx->ts_monitor && !param->PreLoading) ? PETSC_TRUE : PETSC_FALSE;


637:   if (param->PreLoading)
638:     max_steps = 1;
639:   else
640:     max_steps = tsCtx->max_steps;
641: 
642:   Initialize(dmmg);

644:   snes = DMMGGetSNES(dmmg);

646:   for (tsCtx->itstep = 0; tsCtx->itstep < max_steps; tsCtx->itstep++) {
647: 
648:     PetscPrintf(PETSC_COMM_WORLD, "time step=%d!\n",tsCtx->itstep);
649:     if ((param->second_order) && (tsCtx->itstep > 0))
650:     {
651:       for (i=param->mglevels-1; i>=0 ;i--)
652:       {
653:         VecCopy(((AppCtx*)DMMGGetUser(dmmg,i))->Xold,((AppCtx*)DMMGGetUser(dmmg,i))->Xoldold);
654:       }
655:     }

657:     for (i=param->mglevels-1; i>0 ;i--) {
658:       MatRestrict(dmmg[i]->R, dmmg[i]->x, dmmg[i-1]->x);
659:       VecPointwiseMult(dmmg[i-1]->x,dmmg[i-1]->x,dmmg[i]->Rscale);
660:       VecCopy(dmmg[i]->x, ((AppCtx*)DMMGGetUser(dmmg,i))->Xold);
661:     }
662:     VecCopy(dmmg[0]->x, ((AppCtx*)DMMGGetUser(dmmg,0))->Xold);

664:     DMMGSolve(dmmg);



668:     if (tsCtx->itstep == 665000)
669:     {
670:       KSP          ksp;
671:       PC           pc;
672:       Mat          mat, pmat;
673:       MatStructure flag;
674:       PetscViewer  viewer;
675:       char         file[PETSC_MAX_PATH_LEN];

677:       PetscStrcpy(file, "matrix");

679:       PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_WRITE, &viewer);

681:       SNESGetKSP(snes, &ksp);

683:       KSPGetPC(ksp, &pc);

685:       PCGetOperators(pc, &mat, &pmat, &flag);

687:       MatView(mat, viewer);

689:       PetscViewerDestroy(viewer);
690:       SETERRQ(1,"Done saving Jacobian");
691:     }


694:     tsCtx->t += tsCtx->dt;

696:     /* save restart solution if requested at a particular time, then exit */
697:     if (tsCtx->dump_time > 0.0 && tsCtx->t >= tsCtx->dump_time) {
698:       Vec v = DMMGGetx(dmmg);
699:       VecView(v,PETSC_VIEWER_BINARY_WORLD);
700:       SETERRQ1(1,"Saved solution at time %G",tsCtx->t);
701:     }

703:     if (ts_monitor)
704:     {
705:       SNESGetIterationNumber(snes, &its);
706:       SNESGetLinearSolveIterations(snes, &lits);
707:       SNESGetNonlinearStepFailures(snes, &nfails);
708:       SNESGetFunctionNorm(snes, &tsCtx->fnorm);

710:       nfailsCum += nfails;
711:       if (nfailsCum >= 2)
712:         SETERRQ(1, "unable to find a newton step");

714:       PetscPrintf(PETSC_COMM_WORLD,
715:                          "time step = %D, time = %G, number of nonlinear steps = %D, "
716:                          "number of linear steps = %D, norm of the function = %G\n",
717:                          tsCtx->itstep + 1, tsCtx->t, its, lits, PetscAbsScalar(tsCtx->fnorm));

719:       /* send solution over to Matlab, to be visualized (using ex29.m) */
720:       if (!param->PreLoading && tsCtx->socketviewer)
721:       {
722:         Vec v;
723:         SNESGetSolution(snes, &v);
724: #if defined(PETSC_USE_SOCKET_VIEWER)
725:         VecView(v, tsCtx->socketviewer);
726: #endif
727:       }
728:     }

730:     if (!param->PreLoading) {
731:       if (param->draw_contours) {
732:         VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
733:       }

735:       if (1 && ts_monitor) {
736:         /* compute maxima */
737:         ComputeMaxima((DA) dmmg[param->mglevels-1]->dm,dmmg[param->mglevels-1]->x,tsCtx->t);
738:         /* output */
739:         if (ic_out++ == (int)(tsCtx->dt_out / tsCtx->dt + .5)) {
740: #ifdef HAVE_DA_HDF
741:           char fname[PETSC_MAX_PATH_LEN];

743:           sprintf(fname, "out-%G.hdf", tsCtx->t);
744:           DAVecHDFOutput(DMMGGetDA(dmmg), DMMGGetx(dmmg), fname);
745: #else
746: /*
747:           Gnuplot(DMMGGetDA(dmmg), DMMGGetx(dmmg), tsCtx->t);
748:           
749: */
750: #endif
751:           ic_out = 1;
752:         }
753:       }
754:     }
755:   } /* End of time step loop */
756: 
757:   if (!param->PreLoading){
758:     SNESGetFunctionNorm(snes,&tsCtx->fnorm);
759:     PetscPrintf(PETSC_COMM_WORLD, "timesteps %D fnorm = %G\n",
760:                        tsCtx->itstep, PetscAbsScalar(tsCtx->fnorm));
761:   }

763:   if (param->PreLoading) {
764:     tsCtx->fnorm_ini = 0.0;
765:   }
766: 
767:   return(0);
768: }

770: /*---------------------------------------------------------------------*/
773: PetscErrorCode AddTSTermLocal(DALocalInfo* info,Field **x,Field **f,AppCtx *user)
774: /*---------------------------------------------------------------------*/
775: {
776:   TstepCtx       *tsCtx = user->tsCtx;
777:   DA             da = info->da;
779:   PetscInt       i,j;
780:   PetscInt       xints,xinte,yints,yinte;
781:   PassiveReal    hx,hy,dhx,dhy,hxhy;
782:   PassiveReal    one = 1.0;
783:   PassiveScalar  dtinv;
784:   PassiveField   **xold;


788:   xints = info->xs; xinte = info->xs+info->xm;
789:   yints = info->ys; yinte = info->ys+info->ym;

791:   dhx  = info->mx/lx;            dhy = info->my/ly;
792:   hx   = one/dhx;                 hy = one/dhy;
793:   hxhy = hx*hy;

795:   dtinv = hxhy/(tsCtx->dt);

797:   DAVecGetArray(da,user->Xold,&xold);
798:   for (j=yints; j<yinte; j++) {
799:     for (i=xints; i<xinte; i++) {
800:       f[j][i].U += dtinv*(x[j][i].U-xold[j][i].U);
801:       f[j][i].F += dtinv*(x[j][i].F-xold[j][i].F);
802:     }
803:   }
804:   DAVecRestoreArray(da,user->Xold,&xold);

806:   return(0);
807: }

809: /*---------------------------------------------------------------------*/
812: PetscErrorCode AddTSTermLocal2(DALocalInfo* info,Field **x,Field **f,AppCtx *user)
813: /*---------------------------------------------------------------------*/
814: {
815:   TstepCtx       *tsCtx = user->tsCtx;
816:   DA             da = info->da;
818:   PetscInt       i,j,xints,xinte,yints,yinte;
819:   PassiveReal    hx,hy,dhx,dhy,hxhy;
820:   PassiveReal    one = 1.0, onep5 = 1.5, two = 2.0, p5 = 0.5;
821:   PassiveScalar  dtinv;
822:   PassiveField   **xoldold,**xold;


826:   xints = info->xs; xinte = info->xs+info->xm;
827:   yints = info->ys; yinte = info->ys+info->ym;

829:   dhx  = info->mx/lx;            dhy = info->my/ly;
830:   hx   = one/dhx;                 hy = one/dhy;
831:   hxhy = hx*hy;

833:   dtinv = hxhy/(tsCtx->dt);

835:   DAVecGetArray(da,user->Xoldold,&xoldold);
836:   DAVecGetArray(da,user->Xold,&xold);
837:   for (j=yints; j<yinte; j++) {
838:     for (i=xints; i<xinte; i++) {
839:       f[j][i].U += dtinv * (onep5 * x[j][i].U - two * xold[j][i].U +
840:                             p5 * xoldold[j][i].U);
841:       f[j][i].F += dtinv * (onep5 * x[j][i].F - two * xold[j][i].F +
842:                             p5 * xoldold[j][i].F);
843:     }
844:   }
845:   DAVecRestoreArray(da,user->Xoldold,&xoldold);
846:   DAVecRestoreArray(da,user->Xold,&xold);

848:   return(0);
849: }

851: /* Creates the null space of the Jacobian for a particular level */
854: PetscErrorCode CreateNullSpace(DMMG dmmg,Vec *nulls)
855: {
857:   PetscInt       i,N,rstart,rend;
858:   PetscScalar    scale,*vx;

861:   VecGetSize(nulls[0],&N);
862:   scale = 2.0/sqrt((PetscReal)N);
863:   VecGetArray(nulls[0],&vx);
864:   VecGetOwnershipRange(nulls[0],&rstart,&rend);
865:   for (i=rstart; i<rend; i++) {
866:     if (!(i % 4)) vx[i-rstart] = scale;
867:     else          vx[i-rstart] = 0.0;
868:   }
869:   VecRestoreArray(nulls[0],&vx);
870:   return(0);
871: }

873: /*
874:     This is an experimental function and can be safely ignored.
875: */
876: PetscErrorCode FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
877: {
878:   AppCtx         *user = (AppCtx*)ptr;
879:   TstepCtx       *tsCtx = user->tsCtx;
880:   Parameters      *param = user->param;
882:   PetscInt       i,j,c;
883:   PetscInt       xints,xinte,yints,yinte;
884:   PassiveReal    hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
885:   PassiveReal    de2,rhos2,nu,eta,dde2;
886:   PassiveReal    two = 2.0,one = 1.0,p5 = 0.5;
887:   PassiveReal    F_eq_x,By_eq,dtinv;
888:   PetscScalar    xx;
889:   PetscScalar    vx,vy,avx,avy,vxp,vxm,vyp,vym;
890:   PetscScalar    Bx,By,aBx,aBy,Bxp,Bxm,Byp,Bym;
891:   PassiveField   **xold;

894:   de2     = sqr(param->d_e);
895:   rhos2   = sqr(param->rho_s);
896:   nu      = param->nu;
897:   eta     = param->eta;
898:   dde2    = one/de2;

900:   /* 
901:      Define mesh intervals ratios for uniform grid.
902:      [Note: FD formulae below are normalized by multiplying through by
903:      local volume element to obtain coefficients O(1) in two dimensions.]
904:   */
905:   dhx   = info->mx/lx;        dhy   = info->my/ly;
906:   hx    = one/dhx;             hy   = one/dhy;
907:   hxdhy = hx*dhy;           hydhx   = hy*dhx;
908:   hxhy  = hx*hy;             dhxdhy = dhx*dhy;

910:   xints = info->xs; xinte = info->xs+info->xm;
911:   yints = info->ys; yinte = info->ys+info->ym;


914:   i = st->i; j = st->j; c = st->c;

916: #ifdef EQ
917:       xx = i * hx;
918:       F_eq_x = - (1. + de2) * sin(PetscAbsScalar(xx));
919:       By_eq = sin(PetscAbsScalar(xx));
920: #else
921:       F_eq_x = 0.;
922:       By_eq = 0.;
923: #endif

925:       /*
926:        * convective coefficients for upwinding
927:        */

929:       vx  = - D_y(x,phi,i,j);
930:       vy  =   D_x(x,phi,i,j);
931:       avx = PetscAbsScalar(vx); vxp = p5*(vx+avx); vxm = p5*(vx-avx);
932:       avy = PetscAbsScalar(vy); vyp = p5*(vy+avy); vym = p5*(vy-avy);
933: #ifndef UPWINDING
934:       vxp = vxm = p5*vx;
935:       vyp = vym = p5*vy;
936: #endif

938:       Bx  =   D_y(x,psi,i,j);
939:       By  = - D_x(x,psi,i,j) + By_eq;
940:       aBx = PetscAbsScalar(Bx); Bxp = p5*(Bx+aBx); Bxm = p5*(Bx-aBx);
941:       aBy = PetscAbsScalar(By); Byp = p5*(By+aBy); Bym = p5*(By-aBy);
942: #ifndef UPWINDING
943:       Bxp = Bxm = p5*Bx;
944:       Byp = Bym = p5*By;
945: #endif

947:       DAVecGetArray(info->da,user->Xold,&xold);
948:       dtinv = hxhy/(tsCtx->dt);
949:       switch(c) {

951:         case 0:
952:           /* Lap(phi) - U */
953:           *f = (Lapl(x,phi,i,j) - x[j][i].U) * hxhy;
954:           break;

956:         case 1:
957:           /* psi - d_e^2 * Lap(psi) - F */
958:           *f = (x[j][i].psi - de2 * Lapl(x,psi,i,j) - x[j][i].F) * hxhy;
959:           break;

961:         case 2:
962:           /* vx * U_x + vy * U_y - (B_x * F_x + B_y * F_y) / d_e^2
963:             - nu Lap(U) */
964:           *f  =
965:             ((vxp * (D_xm(x,U,i,j)) + vxm * (D_xp(x,U,i,j)) +
966:               vyp * (D_ym(x,U,i,j)) + vym * (D_yp(x,U,i,j))) -
967:              (Bxp * (D_xm(x,F,i,j) + F_eq_x) + Bxm * (D_xp(x,F,i,j) + F_eq_x) +
968:               Byp * (D_ym(x,F,i,j)) + Bym * (D_yp(x,F,i,j))) * dde2 -
969:              nu * Lapl(x,U,i,j)) * hxhy;
970:           *f += dtinv*(x[j][i].U-xold[j][i].U);
971:           break;

973:         case 3:
974:           /* vx * F_x + vy * F_y - rho_s^2 * (B_x * U_x + B_y * U_y)
975:            - eta * Lap(psi) */
976:           *f  =
977:             ((vxp * (D_xm(x,F,i,j) + F_eq_x) + vxm * (D_xp(x,F,i,j) + F_eq_x) +
978:              vyp * (D_ym(x,F,i,j)) + vym * (D_yp(x,F,i,j))) -
979:             (Bxp * (D_xm(x,U,i,j)) + Bxm * (D_xp(x,U,i,j)) +
980:              Byp * (D_ym(x,U,i,j)) + Bym * (D_yp(x,U,i,j))) * rhos2 -
981:             eta * Lapl(x,psi,i,j)) * hxhy;
982:           *f += dtinv*(x[j][i].F-xold[j][i].F);
983:           break;
984:       }
985:       DAVecRestoreArray(info->da,user->Xold,&xold);


988:   /*
989:      Flop count (multiply-adds are counted as 2 operations)
990:   */
991:   /*  PetscLogFlops(84*info->ym*info->xm); FIXME */
992:   return(0);
993: }
994: /*
995:     This is an experimental function and can be safely ignored.
996: */
997: PetscErrorCode FormFunctionLocali4(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *ff,void *ptr)
998: {
999:   Field          *f  = (Field *)ff;
1000:   AppCtx         *user = (AppCtx*)ptr;
1001:   TstepCtx       *tsCtx = user->tsCtx;
1002:   Parameters     *param = user->param;
1004:   PetscInt       i,j;
1005:   PetscInt       xints,xinte,yints,yinte;
1006:   PassiveReal    hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
1007:   PassiveReal    de2,rhos2,nu,eta,dde2;
1008:   PassiveReal    two = 2.0,one = 1.0,p5 = 0.5;
1009:   PassiveReal    F_eq_x,By_eq,dtinv;
1010:   PetscScalar    xx;
1011:   PetscScalar    vx,vy,avx,avy,vxp,vxm,vyp,vym;
1012:   PetscScalar    Bx,By,aBx,aBy,Bxp,Bxm,Byp,Bym;
1013:   PassiveField   **xold;

1016:   de2     = sqr(param->d_e);
1017:   rhos2   = sqr(param->rho_s);
1018:   nu      = param->nu;
1019:   eta     = param->eta;
1020:   dde2    = one/de2;

1022:   /* 
1023:      Define mesh intervals ratios for uniform grid.
1024:      [Note: FD formulae below are normalized by multiplying through by
1025:      local volume element to obtain coefficients O(1) in two dimensions.]
1026:   */
1027:   dhx   = info->mx/lx;        dhy   = info->my/ly;
1028:   hx    = one/dhx;             hy   = one/dhy;
1029:   hxdhy = hx*dhy;           hydhx   = hy*dhx;
1030:   hxhy  = hx*hy;             dhxdhy = dhx*dhy;

1032:   xints = info->xs; xinte = info->xs+info->xm;
1033:   yints = info->ys; yinte = info->ys+info->ym;


1036:   i = st->i; j = st->j;

1038: #ifdef EQ
1039:       xx = i * hx;
1040:       F_eq_x = - (1. + de2) * sin(PetscAbsScalar(xx));
1041:       By_eq = sin(PetscAbsScalar(xx));
1042: #else
1043:       F_eq_x = 0.;
1044:       By_eq = 0.;
1045: #endif

1047:       /*
1048:        * convective coefficients for upwinding
1049:        */

1051:       vx  = - D_y(x,phi,i,j);
1052:       vy  =   D_x(x,phi,i,j);
1053:       avx = PetscAbsScalar(vx); vxp = p5*(vx+avx); vxm = p5*(vx-avx);
1054:       avy = PetscAbsScalar(vy); vyp = p5*(vy+avy); vym = p5*(vy-avy);
1055: #ifndef UPWINDING
1056:       vxp = vxm = p5*vx;
1057:       vyp = vym = p5*vy;
1058: #endif

1060:       Bx  =   D_y(x,psi,i,j);
1061:       By  = - D_x(x,psi,i,j) + By_eq;
1062:       aBx = PetscAbsScalar(Bx); Bxp = p5*(Bx+aBx); Bxm = p5*(Bx-aBx);
1063:       aBy = PetscAbsScalar(By); Byp = p5*(By+aBy); Bym = p5*(By-aBy);
1064: #ifndef UPWINDING
1065:       Bxp = Bxm = p5*Bx;
1066:       Byp = Bym = p5*By;
1067: #endif

1069:       DAVecGetArray(info->da,user->Xold,&xold);
1070:       dtinv = hxhy/(tsCtx->dt);
1071: 

1073: 
1074:           /* Lap(phi) - U */
1075:           f->phi = (Lapl(x,phi,i,j) - x[j][i].U) * hxhy;
1076: 

1078: 
1079:           /* psi - d_e^2 * Lap(psi) - F */
1080:           f->psi = (x[j][i].psi - de2 * Lapl(x,psi,i,j) - x[j][i].F) * hxhy;
1081: 

1083: 
1084:           /* vx * U_x + vy * U_y - (B_x * F_x + B_y * F_y) / d_e^2
1085:             - nu Lap(U) */
1086:           f->U  =
1087:             ((vxp * (D_xm(x,U,i,j)) + vxm * (D_xp(x,U,i,j)) +
1088:               vyp * (D_ym(x,U,i,j)) + vym * (D_yp(x,U,i,j))) -
1089:              (Bxp * (D_xm(x,F,i,j) + F_eq_x) + Bxm * (D_xp(x,F,i,j) + F_eq_x) +
1090:               Byp * (D_ym(x,F,i,j)) + Bym * (D_yp(x,F,i,j))) * dde2 -
1091:              nu * Lapl(x,U,i,j)) * hxhy;
1092:           f->U += dtinv*(x[j][i].U-xold[j][i].U);
1093: 

1095: 
1096:           /* vx * F_x + vy * F_y - rho_s^2 * (B_x * U_x + B_y * U_y)
1097:            - eta * Lap(psi) */
1098:           f->F  =
1099:             ((vxp * (D_xm(x,F,i,j) + F_eq_x) + vxm * (D_xp(x,F,i,j) + F_eq_x) +
1100:              vyp * (D_ym(x,F,i,j)) + vym * (D_yp(x,F,i,j))) -
1101:             (Bxp * (D_xm(x,U,i,j)) + Bxm * (D_xp(x,U,i,j)) +
1102:              Byp * (D_ym(x,U,i,j)) + Bym * (D_yp(x,U,i,j))) * rhos2 -
1103:             eta * Lapl(x,psi,i,j)) * hxhy;
1104:           f->F += dtinv*(x[j][i].F-xold[j][i].F);
1105: 
1106: 
1107:       DAVecRestoreArray(info->da,user->Xold,&xold);


1110:   /*
1111:      Flop count (multiply-adds are counted as 2 operations)
1112:   */
1113:   /*  PetscLogFlops(84*info->ym*info->xm); FIXME */
1114:   return(0);
1115: }