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
47: #ifdef HAVE_DA_HDF
48: PetscInt DAVecHDFOutput(DA,Vec,char*);
49: #endif
51: #define D_x(x,m,i,j) (p5 * (x[(j)][(i)+1].m - x[(j)][(i)-1].m) * dhx)
52: #define D_xm(x,m,i,j) ((x[(j)][(i)].m - x[(j)][(i)-1].m) * dhx)
53: #define D_xp(x,m,i,j) ((x[(j)][(i+1)].m - x[(j)][(i)].m) * dhx)
54: #define D_y(x,m,i,j) (p5 * (x[(j)+1][(i)].m - x[(j)-1][(i)].m) * dhy)
55: #define D_ym(x,m,i,j) ((x[(j)][(i)].m - x[(j)-1][(i)].m) * dhy)
56: #define D_yp(x,m,i,j) ((x[(j)+1][(i)].m - x[(j)][(i)].m) * dhy)
57: #define D_xx(x,m,i,j) ((x[(j)][(i)+1].m - two*x[(j)][(i)].m + x[(j)][(i)-1].m) * hydhx * dhxdhy)
58: #define D_yy(x,m,i,j) ((x[(j)+1][(i)].m - two*x[(j)][(i)].m + x[(j)-1][(i)].m) * hxdhy * dhxdhy)
59: #define Lapl(x,m,i,j) (D_xx(x,m,i,j) + D_yy(x,m,i,j))
60: #define lx (2.*PETSC_PI)
61: #define ly (4.*PETSC_PI)
62: #define sqr(a) ((a)*(a))
64: /*
65: User-defined routines and data structures
66: */
68: typedef struct {
69: PetscReal fnorm_ini,dt_ini;
70: PetscReal dt,dt_out;
71: PetscReal ptime;
72: PetscReal max_time;
73: PetscReal fnorm_ratio;
74: PetscInt ires,itstep;
75: PetscInt max_steps,print_freq;
76: PetscReal t;
77: PetscScalar fnorm;
79: PetscTruth ts_monitor; /* print information about each time step */
80: PetscReal dump_time; /* time to dump solution to a file */
81: PetscViewer socketviewer; /* socket to dump solution at each timestep for visualization */
82: } TstepCtx;
84: typedef struct {
85: PetscScalar phi,psi,U,F;
86: } Field;
88: typedef struct {
89: PassiveScalar phi,psi,U,F;
90: } PassiveField;
92: typedef struct {
93: PetscInt mglevels;
94: PetscInt cycles; /* number of time steps for integration */
95: PassiveReal nu,eta,d_e,rho_s; /* physical parameters */
96: PetscTruth draw_contours; /* flag - 1 indicates drawing contours */
97: PetscTruth second_order;
98: PetscTruth PreLoading;
99: } Parameter;
101: typedef struct {
102: Vec Xoldold,Xold;
103: TstepCtx *tsCtx;
104: Parameter *param;
105: } AppCtx;
118: int main(int argc,char **argv)
119: {
121: DMMG *dmmg; /* multilevel grid structure */
122: AppCtx *user; /* user-defined work context (one for each level) */
123: TstepCtx tsCtx; /* time-step parameters (one total) */
124: Parameter param; /* physical parameters (one total) */
125: PetscInt i,m,n,mx,my;
126: DA da;
127: PetscTruth flg;
128: PetscReal dt_ratio;
129: PetscInt dfill[16] = {1,0,1,0,
130: 0,1,0,1,
131: 0,0,1,1,
132: 0,1,1,1};
133: PetscInt ofill[16] = {1,0,0,0,
134: 0,1,0,0,
135: 1,1,1,1,
136: 1,1,1,1};
138: PetscInitialize(&argc,&argv,(char *)0,help);
141: PreLoadBegin(PETSC_TRUE,"SetUp");
143: param.PreLoading = PreLoading;
144: DMMGCreate(PETSC_COMM_WORLD,1,&user,&dmmg);
145: param.mglevels = DMMGGetLevels(dmmg);
147: /*
148: Create distributed array multigrid object (DMMG) to manage
149: parallel grid and vectors for principal unknowns (x) and
150: governing residuals (f)
151: */
152: DACreate2d(PETSC_COMM_WORLD, DA_XYPERIODIC, DA_STENCIL_STAR, -5, -5,
153: PETSC_DECIDE, PETSC_DECIDE, 4, 1, 0, 0, &da);
155: /* overwrite the default sparsity pattern toone specific for
156: this code's nonzero structure */
157: DASetBlockFills(da,dfill,ofill);
159: DMMGSetDM(dmmg,(DM)da);
160: DADestroy(da);
162: /* default physical parameters */
163: param.nu = 0;
164: param.eta = 1e-3;
165: param.d_e = 0.2;
166: param.rho_s = 0;
168: PetscOptionsGetReal(PETSC_NULL, "-viscosity", ¶m.nu,PETSC_NULL);
170: PetscOptionsGetReal(PETSC_NULL, "-resistivity", ¶m.eta,PETSC_NULL);
172: PetscOptionsGetReal(PETSC_NULL, "-skin_depth", ¶m.d_e,PETSC_NULL);
174: PetscOptionsGetReal(PETSC_NULL, "-larmor_radius", ¶m.rho_s,PETSC_NULL);
176: PetscOptionsHasName(PETSC_NULL, "-contours", ¶m.draw_contours);
178: PetscOptionsHasName(PETSC_NULL, "-second_order", ¶m.second_order);
180: DASetFieldName(DMMGGetDA(dmmg), 0, "phi");
182: DASetFieldName(DMMGGetDA(dmmg), 1, "psi");
184: DASetFieldName(DMMGGetDA(dmmg), 2, "U");
186: DASetFieldName(DMMGGetDA(dmmg), 3, "F");
188: /*======================================================================*/
189: /* Initialize stuff related to time stepping */
190: /*======================================================================*/
192: DAGetInfo(DMMGGetDA(dmmg),0,&mx,&my,0,0,0,0,0,0,0,0);
194: tsCtx.fnorm_ini = 0;
195: tsCtx.max_steps = 1000000;
196: tsCtx.max_time = 1.0e+12;
197: /* use for dt = min(dx,dy); multiplied by dt_ratio below */
198: tsCtx.dt = PetscMin(lx/mx,ly/my);
199: tsCtx.fnorm_ratio = 1.0e+10;
200: tsCtx.t = 0;
201: tsCtx.dt_out = 10;
202: tsCtx.print_freq = tsCtx.max_steps;
203: tsCtx.ts_monitor = PETSC_FALSE;
204: tsCtx.dump_time = -1.0;
206: PetscOptionsGetInt(PETSC_NULL, "-max_st", &tsCtx.max_steps,PETSC_NULL);
207: PetscOptionsGetReal(PETSC_NULL, "-ts_rtol", &tsCtx.fnorm_ratio,PETSC_NULL);
208: PetscOptionsGetInt(PETSC_NULL, "-print_freq", &tsCtx.print_freq,PETSC_NULL);
210: dt_ratio = 1.0;
211: PetscOptionsGetReal(PETSC_NULL, "-dt_ratio", &dt_ratio,PETSC_NULL);
212: tsCtx.dt *= dt_ratio;
214: PetscOptionsHasName(PETSC_NULL, "-ts_monitor", &tsCtx.ts_monitor);
215: PetscOptionsGetReal(PETSC_NULL, "-dump", &tsCtx.dump_time,PETSC_NULL);
217: tsCtx.socketviewer = 0;
218: PetscOptionsHasName(PETSC_NULL, "-socket_viewer", &flg);
219: if (flg && !PreLoading) {
220: tsCtx.ts_monitor = PETSC_TRUE;
221: PetscViewerSocketOpen(PETSC_COMM_WORLD,0,PETSC_DECIDE,&tsCtx.socketviewer);
222: }
223:
225: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226: Create user context, set problem data, create vector data structures.
227: Also, compute the initial guess.
228: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: /* create application context for each level */
230: PetscMalloc(param.mglevels*sizeof(AppCtx),&user);
231: for (i=0; i<param.mglevels; i++) {
232: /* create work vectors to hold previous time-step solution and
233: function value */
234: VecDuplicate(dmmg[i]->x, &user[i].Xoldold);
235: VecDuplicate(dmmg[i]->x, &user[i].Xold);
236: user[i].tsCtx = &tsCtx;
237: user[i].param = ¶m;
238: DMMGSetUser(dmmg,i,&user[i]);
239: }
240:
241: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242: Create nonlinear solver context
243:
244: Process adiC(20): AddTSTermLocal AddTSTermLocal2 FormFunctionLocal FormFunctionLocali
245: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246: DMMGSetSNESLocal(dmmg, FormFunctionLocal, 0,ad_FormFunctionLocal, admf_FormFunctionLocal);
247: DMMGSetSNESLocali(dmmg,FormFunctionLocali,ad_FormFunctionLocali,admf_FormFunctionLocali);
249: /* attach nullspace to each level of the preconditioner */
250: DMMGSetNullSpace(dmmg,PETSC_FALSE,1,CreateNullSpace);
252: if (PreLoading) {
253: PetscPrintf(PETSC_COMM_WORLD, "# viscosity = %g, resistivity = %g, "
254: "skin_depth # = %g, larmor_radius # = %g\n",
255: param.nu, param.eta, param.d_e, param.rho_s);
256: DAGetInfo(DMMGGetDA(dmmg),0,&m,&n,0,0,0,0,0,0,0,0);
257: PetscPrintf(PETSC_COMM_WORLD,"Problem size %D by %D\n",m,n);
258: 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);
259: }
261:
262:
263: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264: Solve the nonlinear system
265: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
266:
267: PreLoadStage("Solve");
269: if (param.draw_contours) {
270: VecView(((AppCtx*)DMMGGetUser(dmmg,param.mglevels-1))->Xold,PETSC_VIEWER_DRAW_WORLD);
271: }
273: Update(dmmg);
275: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276: Free work space. All PETSc objects should be destroyed when they
277: are no longer needed.
278: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
279:
280: for (i=0; i<param.mglevels; i++) {
281: VecDestroy(user[i].Xoldold);
282: VecDestroy(user[i].Xold);
283: }
284: PetscFree(user);
285: DMMGDestroy(dmmg);
287: PreLoadEnd();
288:
289: PetscFinalize();
290: return 0;
291: }
293: /* ------------------------------------------------------------------- */
296: /* ------------------------------------------------------------------- */
297: PetscErrorCode Gnuplot(DA da, Vec X, double mtime)
298: {
299: PetscInt i,j,xs,ys,xm,ym;
300: PetscInt xints,xinte,yints,yinte;
302: Field **x;
303: FILE *f;
304: char fname[PETSC_MAX_PATH_LEN];
305: PetscMPIInt rank;
307: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
309: sprintf(fname, "out-%g-%d.dat", mtime, rank);
311: f = fopen(fname, "w");
313: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
315: DAVecGetArray(da,X,&x);
317: xints = xs; xinte = xs+xm; yints = ys; yinte = ys+ym;
319: for (j=yints; j<yinte; j++) {
320: for (i=xints; i<xinte; i++) {
321: PetscFPrintf(PETSC_COMM_WORLD, f,
322: "%D %D %g %g %g %g %g %g\n",
323: i, j, 0.0, 0.0,
324: PetscAbsScalar(x[j][i].U), PetscAbsScalar(x[j][i].F),
325: PetscAbsScalar(x[j][i].phi), PetscAbsScalar(x[j][i].psi));
326: }
327: PetscFPrintf(PETSC_COMM_WORLD,f, "\n");
328: }
329: DAVecRestoreArray(da,X,&x);
330: fclose(f);
331: return 0;
332: }
334: /* ------------------------------------------------------------------- */
337: /* ------------------------------------------------------------------- */
338: PetscErrorCode Initialize(DMMG *dmmg)
339: {
340: AppCtx *appCtx = (AppCtx*)DMMGGetUser(dmmg,0);
341: Parameter *param = appCtx->param;
342: DA da;
344: PetscInt i,j,mx,my,xs,ys,xm,ym;
345: PetscReal two = 2.0,one = 1.0;
346: PetscReal hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
347: PetscReal d_e,rho_s,de2,xx,yy;
348: Field **x, **localx;
349: Vec localX;
350: PetscTruth flg;
353: PetscOptionsHasName(0,"-restart",&flg);
354: if (flg) {
355: PetscViewer viewer;
356: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"binaryoutput",PETSC_FILE_RDONLY,&viewer);
357: VecLoadIntoVector(viewer,dmmg[param->mglevels-1]->x);
358: PetscViewerDestroy(viewer);
359: return(0);
360: }
362: d_e = param->d_e;
363: rho_s = param->rho_s;
364: de2 = sqr(param->d_e);
366: da = (DA)(dmmg[param->mglevels-1]->dm);
367: DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);
369: dhx = mx/lx; dhy = my/ly;
370: hx = one/dhx; hy = one/dhy;
371: hxdhy = hx*dhy; hydhx = hy*dhx;
372: hxhy = hx*hy; dhxdhy = dhx*dhy;
374: /*
375: Get local grid boundaries (for 2-dimensional DA):
376: xs, ys - starting grid indices (no ghost points)
377: xm, ym - widths of local grid (no ghost points)
378: */
379: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
381: DAGetLocalVector(da,&localX);
382: /*
383: Get a pointer to vector data.
384: - For default PETSc vectors, VecGetArray() returns a pointer to
385: the data array. Otherwise, the routine is implementation dependent.
386: - You MUST call VecRestoreArray() when you no longer need access to
387: the array.
388: */
389: DAVecGetArray(da,dmmg[param->mglevels-1]->x,&x);
390: DAVecGetArray(da,localX,&localx);
392: /*
393: Compute initial solution over the locally owned part of the grid
394: */
395: {
396: PetscReal eps = lx/ly;
397: #if defined(PETSC_HAVE_ERF)
398: PetscReal pert = 1e-4;
399: #endif
400: PetscReal k = 1.*eps;
401: PetscReal gam;
403: if (d_e < rho_s) d_e = rho_s;
404: gam = k * d_e;
406: for (j=ys-1; j<ys+ym+1; j++) {
407: yy = j * hy;
408: for (i=xs-1; i<xs+xm+1; i++) {
409: xx = i * hx;
410: #if defined(PETSC_HAVE_ERF)
411: if (xx < -PETSC_PI/2) {
412: localx[j][i].phi = pert * gam / k * erf((xx + PETSC_PI) / (sqrt(2.0) * d_e)) * (-sin(k*yy));
413: } else if (xx < PETSC_PI/2) {
414: localx[j][i].phi = - pert * gam / k * erf(xx / (sqrt(2.0) * d_e)) * (-sin(k*yy));
415: } else if (xx < 3*PETSC_PI/2){
416: localx[j][i].phi = pert * gam / k * erf((xx - PETSC_PI) / (sqrt(2.0) * d_e)) * (-sin(k*yy));
417: } else {
418: localx[j][i].phi = - pert * gam / k * erf((xx - 2.*PETSC_PI) / (sqrt(2.0) * d_e)) * (-sin(k*yy));
419: }
420: #endif
421: #ifdef EQ
422: localx[j][i].psi = 0.;
423: #else
424: localx[j][i].psi = cos(xx);
425: #endif
426: }
427: }
428: for (j=ys; j<ys+ym; j++) {
429: for (i=xs; i<xs+xm; i++) {
430: x[j][i].U = Lapl(localx,phi,i,j);
431: x[j][i].F = localx[j][i].psi - de2 * Lapl(localx,psi,i,j);
432: x[j][i].phi = localx[j][i].phi;
433: x[j][i].psi = localx[j][i].psi;
434: }
435: }
436: }
437: /*
438: Restore vector
439: */
440: DAVecRestoreArray(da,dmmg[param->mglevels-1]->x,&x);
441:
442: DAVecRestoreArray(da,localX,&localx);
443:
444: DARestoreLocalVector(da,&localX);
445:
447: return(0);
448: }
452: PetscErrorCode ComputeMaxima(DA da, Vec X, PetscReal t)
453: {
455: PetscInt i,j,mx,my,xs,ys,xm,ym;
456: PetscInt xints,xinte,yints,yinte;
457: Field **x;
458: double norm[4] = {0,0,0,0};
459: double gnorm[4];
460: MPI_Comm comm;
464: DAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0);
466: DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
467:
469: xints = xs; xinte = xs+xm; yints = ys; yinte = ys+ym;
471: DAVecGetArray(da, X, &x);
473: for (j=yints; j<yinte; j++) {
474: for (i=xints; i<xinte; i++) {
475: norm[0] = PetscMax(norm[0],PetscAbsScalar(x[j][i].U));
476: norm[1] = PetscMax(norm[1],PetscAbsScalar(x[j][i].F));
477: norm[2] = PetscMax(norm[2],PetscAbsScalar(x[j][i].phi));
478: norm[3] = PetscMax(norm[3],PetscAbsScalar(x[j][i].psi));
479: }
480: }
482: DAVecRestoreArray(da,X,&x);
484: PetscObjectGetComm((PetscObject)da, &comm);
485: MPI_Allreduce(norm, gnorm, 4, MPI_DOUBLE, MPI_MAX, comm);
487: PetscFPrintf(PETSC_COMM_WORLD, stderr,"%g\t%g\t%g\t%g\t%g\n",t, gnorm[0], gnorm[1], gnorm[2], gnorm[3]);
489: return(0);
490: }
494: PetscErrorCode FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
495: {
496: AppCtx *user = (AppCtx*)ptr;
497: TstepCtx *tsCtx = user->tsCtx;
498: Parameter *param = user->param;
500: PetscInt xints,xinte,yints,yinte,i,j;
501: PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
502: PassiveReal de2,rhos2,nu,eta,dde2;
503: PassiveReal two = 2.0,one = 1.0,p5 = 0.5;
504: PassiveReal F_eq_x,By_eq;
505: PetscScalar xx;
506: PetscScalar vx,vy,avx,avy,vxp,vxm,vyp,vym;
507: PetscScalar Bx,By,aBx,aBy,Bxp,Bxm,Byp,Bym;
510: de2 = sqr(param->d_e);
511: rhos2 = sqr(param->rho_s);
512: nu = param->nu;
513: eta = param->eta;
514: dde2 = one/de2;
516: /*
517: Define mesh intervals ratios for uniform grid.
518: [Note: FD formulae below are normalized by multiplying through by
519: local volume element to obtain coefficients O(1) in two dimensions.]
520: */
521: dhx = info->mx/lx; dhy = info->my/ly;
522: hx = one/dhx; hy = one/dhy;
523: hxdhy = hx*dhy; hydhx = hy*dhx;
524: hxhy = hx*hy; dhxdhy = dhx*dhy;
526: xints = info->xs; xinte = info->xs+info->xm;
527: yints = info->ys; yinte = info->ys+info->ym;
529: /* Compute over the interior points */
530: for (j=yints; j<yinte; j++) {
531: for (i=xints; i<xinte; i++) {
532: #ifdef EQ
533: xx = i * hx;
534: F_eq_x = - (1. + de2) * sin(PetscAbsScalar(xx));
535: By_eq = sin(PetscAbsScalar(xx));
536: #else
537: F_eq_x = 0.;
538: By_eq = 0.;
539: #endif
541: /*
542: * convective coefficients for upwinding
543: */
545: vx = - D_y(x,phi,i,j);
546: vy = D_x(x,phi,i,j);
547: avx = PetscAbsScalar(vx); vxp = p5*(vx+avx); vxm = p5*(vx-avx);
548: avy = PetscAbsScalar(vy); vyp = p5*(vy+avy); vym = p5*(vy-avy);
549: #ifndef UPWINDING
550: vxp = vxm = p5*vx;
551: vyp = vym = p5*vy;
552: #endif
554: Bx = D_y(x,psi,i,j);
555: By = - D_x(x,psi,i,j) + By_eq;
556: aBx = PetscAbsScalar(Bx); Bxp = p5*(Bx+aBx); Bxm = p5*(Bx-aBx);
557: aBy = PetscAbsScalar(By); Byp = p5*(By+aBy); Bym = p5*(By-aBy);
558: #ifndef UPWINDING
559: Bxp = Bxm = p5*Bx;
560: Byp = Bym = p5*By;
561: #endif
563: /* Lap(phi) - U */
564: f[j][i].phi = (Lapl(x,phi,i,j) - x[j][i].U) * hxhy;
566: /* psi - d_e^2 * Lap(psi) - F */
567: f[j][i].psi = (x[j][i].psi - de2 * Lapl(x,psi,i,j) - x[j][i].F) * hxhy;
569: /* vx * U_x + vy * U_y - (B_x * F_x + B_y * F_y) / d_e^2
570: - nu Lap(U) */
571: f[j][i].U =
572: ((vxp * (D_xm(x,U,i,j)) + vxm * (D_xp(x,U,i,j)) +
573: vyp * (D_ym(x,U,i,j)) + vym * (D_yp(x,U,i,j))) -
574: (Bxp * (D_xm(x,F,i,j) + F_eq_x) + Bxm * (D_xp(x,F,i,j) + F_eq_x) +
575: Byp * (D_ym(x,F,i,j)) + Bym * (D_yp(x,F,i,j))) * dde2 -
576: nu * Lapl(x,U,i,j)) * hxhy;
577:
578: /* vx * F_x + vy * F_y - rho_s^2 * (B_x * U_x + B_y * U_y)
579: - eta * Lap(psi) */
580: f[j][i].F =
581: ((vxp * (D_xm(x,F,i,j) + F_eq_x) + vxm * (D_xp(x,F,i,j) + F_eq_x) +
582: vyp * (D_ym(x,F,i,j)) + vym * (D_yp(x,F,i,j))) -
583: (Bxp * (D_xm(x,U,i,j)) + Bxm * (D_xp(x,U,i,j)) +
584: Byp * (D_ym(x,U,i,j)) + Bym * (D_yp(x,U,i,j))) * rhos2 -
585: eta * Lapl(x,psi,i,j)) * hxhy;
586: }
587: }
589: /* Add time step contribution */
590: if (tsCtx->ires) {
591: if ((param->second_order) && (tsCtx->itstep > 0)){
592: AddTSTermLocal2(info,x,f,user);
593: } else {
594: AddTSTermLocal(info,x,f,user);
595: }
596: }
598: /*
599: Flop count (multiply-adds are counted as 2 operations)
600: */
601: /* PetscLogFlops(84*info->ym*info->xm); FIXME */
602: return(0);
603: }
605: /*---------------------------------------------------------------------*/
608: PetscErrorCode Update(DMMG *dmmg)
609: /*---------------------------------------------------------------------*/
610: {
611: AppCtx *user = (AppCtx *) DMMGGetUser(dmmg,0);
612: TstepCtx *tsCtx = user->tsCtx;
613: Parameter *param = user->param;
614: SNES snes;
615: PetscErrorCode ierr;
616: PetscInt its,lits,i;
617: PetscInt max_steps;
618: PetscInt nfailsCum = 0,nfails = 0;
619: static PetscInt ic_out;
620: PetscTruth ts_monitor = (tsCtx->ts_monitor && !param->PreLoading) ? PETSC_TRUE : PETSC_FALSE;
624: if (param->PreLoading)
625: max_steps = 1;
626: else
627: max_steps = tsCtx->max_steps;
628:
629: Initialize(dmmg);
631: snes = DMMGGetSNES(dmmg);
633: for (tsCtx->itstep = 0; tsCtx->itstep < max_steps; tsCtx->itstep++) {
635: if ((param->second_order) && (tsCtx->itstep > 0))
636: {
637: for (i=param->mglevels-1; i>=0 ;i--)
638: {
639: VecCopy(((AppCtx*)DMMGGetUser(dmmg,i))->Xold,((AppCtx*)DMMGGetUser(dmmg,i))->Xoldold);
640: }
641: }
643: for (i=param->mglevels-1; i>0 ;i--) {
644: MatRestrict(dmmg[i]->R, dmmg[i]->x, dmmg[i-1]->x);
645: VecPointwiseMult(dmmg[i]->Rscale,dmmg[i-1]->x,dmmg[i-1]->x);
646: VecCopy(dmmg[i]->x, ((AppCtx*)DMMGGetUser(dmmg,i))->Xold);
647: }
648: VecCopy(dmmg[0]->x, ((AppCtx*)DMMGGetUser(dmmg,0))->Xold);
650: DMMGSolve(dmmg);
654: if (tsCtx->itstep == 665000)
655: {
656: KSP ksp;
657: PC pc;
658: Mat mat, pmat;
659: MatStructure flag;
660: PetscViewer viewer;
661: char file[PETSC_MAX_PATH_LEN];
663: PetscStrcpy(file, "matrix");
665: PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, PETSC_FILE_CREATE, &viewer);
667: SNESGetKSP(snes, &ksp);
669: KSPGetPC(ksp, &pc);
671: PCGetOperators(pc, &mat, &pmat, &flag);
673: MatView(mat, viewer);
675: PetscViewerDestroy(viewer);
676: SETERRQ(1,"Done saving Jacobian");
677: }
680: tsCtx->t += tsCtx->dt;
682: /* save restart solution if requested at a particular time, then exit */
683: if (tsCtx->dump_time > 0.0 && tsCtx->t >= tsCtx->dump_time) {
684: Vec v = DMMGGetx(dmmg);
685: VecView(v,PETSC_VIEWER_BINARY_WORLD);
686: SETERRQ1(1,"Saved solution at time %g",tsCtx->t);
687: }
689: if (ts_monitor)
690: {
691: SNESGetIterationNumber(snes, &its);
692: SNESGetNumberLinearIterations(snes, &lits);
693: SNESGetNumberUnsuccessfulSteps(snes, &nfails);
694: SNESGetFunctionNorm(snes, &tsCtx->fnorm);
696: nfailsCum += nfails;
697: if (nfailsCum >= 2)
698: SETERRQ(1, "unable to find a newton step");
700: PetscPrintf(PETSC_COMM_WORLD,
701: "time step = %D, time = %g, number of nonlinear steps = %D, "
702: "number of linear steps = %D, norm of the function = %g\n",
703: tsCtx->itstep + 1, tsCtx->t, its, lits, PetscAbsScalar(tsCtx->fnorm));
705: /* send solution over to Matlab, to be visualized (using ex29.m) */
706: if (!param->PreLoading && tsCtx->socketviewer)
707: {
708: Vec v;
709: SNESGetSolution(snes, &v);
710: VecView(v, tsCtx->socketviewer);
711: }
712: }
714: if (!param->PreLoading) {
715: if (param->draw_contours) {
716: VecView(DMMGGetx(dmmg),PETSC_VIEWER_DRAW_WORLD);
717: }
719: if (1 && ts_monitor) {
720: /* compute maxima */
721: ComputeMaxima((DA) dmmg[param->mglevels-1]->dm,dmmg[param->mglevels-1]->x,tsCtx->t);
722: /* output */
723: if (ic_out++ == (int)(tsCtx->dt_out / tsCtx->dt + .5)) {
724: #ifdef HAVE_DA_HDF
725: char fname[PETSC_MAX_PATH_LEN];
727: sprintf(fname, "out-%g.hdf", tsCtx->t);
728: DAVecHDFOutput(DMMGGetDA(dmmg), DMMGGetx(dmmg), fname);
729: #else
730: /*
731: Gnuplot(DMMGGetDA(dmmg), DMMGGetx(dmmg), tsCtx->t);
732:
733: */
734: #endif
735: ic_out = 1;
736: }
737: }
738: }
739: } /* End of time step loop */
740:
741: if (!param->PreLoading){
742: SNESGetFunctionNorm(snes,&tsCtx->fnorm);
743: PetscPrintf(PETSC_COMM_WORLD, "timesteps %D fnorm = %g\n",
744: tsCtx->itstep, PetscAbsScalar(tsCtx->fnorm));
745: }
747: if (param->PreLoading) {
748: tsCtx->fnorm_ini = 0.0;
749: }
750:
751: return(0);
752: }
754: /*---------------------------------------------------------------------*/
757: PetscErrorCode AddTSTermLocal(DALocalInfo* info,Field **x,Field **f,AppCtx *user)
758: /*---------------------------------------------------------------------*/
759: {
760: TstepCtx *tsCtx = user->tsCtx;
761: DA da = info->da;
763: PetscInt i,j;
764: PetscInt xints,xinte,yints,yinte;
765: PassiveReal hx,hy,dhx,dhy,hxhy;
766: PassiveReal one = 1.0;
767: PassiveScalar dtinv;
768: PassiveField **xold;
772: xints = info->xs; xinte = info->xs+info->xm;
773: yints = info->ys; yinte = info->ys+info->ym;
775: dhx = info->mx/lx; dhy = info->my/ly;
776: hx = one/dhx; hy = one/dhy;
777: hxhy = hx*hy;
779: dtinv = hxhy/(tsCtx->dt);
781: DAVecGetArray(da,user->Xold,&xold);
782: for (j=yints; j<yinte; j++) {
783: for (i=xints; i<xinte; i++) {
784: f[j][i].U += dtinv*(x[j][i].U-xold[j][i].U);
785: f[j][i].F += dtinv*(x[j][i].F-xold[j][i].F);
786: }
787: }
788: DAVecRestoreArray(da,user->Xold,&xold);
790: return(0);
791: }
793: /*---------------------------------------------------------------------*/
796: PetscErrorCode AddTSTermLocal2(DALocalInfo* info,Field **x,Field **f,AppCtx *user)
797: /*---------------------------------------------------------------------*/
798: {
799: TstepCtx *tsCtx = user->tsCtx;
800: DA da = info->da;
802: PetscInt i,j,xints,xinte,yints,yinte;
803: PassiveReal hx,hy,dhx,dhy,hxhy;
804: PassiveReal one = 1.0, onep5 = 1.5, two = 2.0, p5 = 0.5;
805: PassiveScalar dtinv;
806: PassiveField **xoldold,**xold;
810: xints = info->xs; xinte = info->xs+info->xm;
811: yints = info->ys; yinte = info->ys+info->ym;
813: dhx = info->mx/lx; dhy = info->my/ly;
814: hx = one/dhx; hy = one/dhy;
815: hxhy = hx*hy;
817: dtinv = hxhy/(tsCtx->dt);
819: DAVecGetArray(da,user->Xoldold,&xoldold);
820: DAVecGetArray(da,user->Xold,&xold);
821: for (j=yints; j<yinte; j++) {
822: for (i=xints; i<xinte; i++) {
823: f[j][i].U += dtinv * (onep5 * x[j][i].U - two * xold[j][i].U +
824: p5 * xoldold[j][i].U);
825: f[j][i].F += dtinv * (onep5 * x[j][i].F - two * xold[j][i].F +
826: p5 * xoldold[j][i].F);
827: }
828: }
829: DAVecRestoreArray(da,user->Xoldold,&xoldold);
830: DAVecRestoreArray(da,user->Xold,&xold);
832: return(0);
833: }
835: /* Creates the null space of the Jacobian for a particular level */
838: PetscErrorCode CreateNullSpace(DMMG dmmg,Vec *nulls)
839: {
841: PetscInt i,N,rstart,rend;
842: PetscScalar scale,*vx;
845: VecGetSize(nulls[0],&N);
846: scale = 2.0/sqrt((PetscReal)N);
847: VecGetArray(nulls[0],&vx);
848: VecGetOwnershipRange(nulls[0],&rstart,&rend);
849: for (i=rstart; i<rend; i++) {
850: if (!(i % 4)) vx[i-rstart] = scale;
851: else vx[i-rstart] = 0.0;
852: }
853: VecRestoreArray(nulls[0],&vx);
854: return(0);
855: }
857: /*
858: This is an experimental function and can be safely ignored.
859: */
860: PetscErrorCode FormFunctionLocali(DALocalInfo *info,MatStencil *st,Field **x,PetscScalar *f,void *ptr)
861: {
862: AppCtx *user = (AppCtx*)ptr;
863: TstepCtx *tsCtx = user->tsCtx;
864: Parameter *param = user->param;
866: PetscInt i,j,c;
867: PetscInt xints,xinte,yints,yinte;
868: PassiveReal hx,hy,dhx,dhy,hxdhy,hydhx,hxhy,dhxdhy;
869: PassiveReal de2,rhos2,nu,eta,dde2;
870: PassiveReal two = 2.0,one = 1.0,p5 = 0.5;
871: PassiveReal F_eq_x,By_eq,dtinv;
872: PetscScalar xx;
873: PetscScalar vx,vy,avx,avy,vxp,vxm,vyp,vym;
874: PetscScalar Bx,By,aBx,aBy,Bxp,Bxm,Byp,Bym;
875: PassiveField **xold;
878: de2 = sqr(param->d_e);
879: rhos2 = sqr(param->rho_s);
880: nu = param->nu;
881: eta = param->eta;
882: dde2 = one/de2;
884: /*
885: Define mesh intervals ratios for uniform grid.
886: [Note: FD formulae below are normalized by multiplying through by
887: local volume element to obtain coefficients O(1) in two dimensions.]
888: */
889: dhx = info->mx/lx; dhy = info->my/ly;
890: hx = one/dhx; hy = one/dhy;
891: hxdhy = hx*dhy; hydhx = hy*dhx;
892: hxhy = hx*hy; dhxdhy = dhx*dhy;
894: xints = info->xs; xinte = info->xs+info->xm;
895: yints = info->ys; yinte = info->ys+info->ym;
898: i = st->i; j = st->j; c = st->c;
900: #ifdef EQ
901: xx = i * hx;
902: F_eq_x = - (1. + de2) * sin(PetscAbsScalar(xx));
903: By_eq = sin(PetscAbsScalar(xx));
904: #else
905: F_eq_x = 0.;
906: By_eq = 0.;
907: #endif
909: /*
910: * convective coefficients for upwinding
911: */
913: vx = - D_y(x,phi,i,j);
914: vy = D_x(x,phi,i,j);
915: avx = PetscAbsScalar(vx); vxp = p5*(vx+avx); vxm = p5*(vx-avx);
916: avy = PetscAbsScalar(vy); vyp = p5*(vy+avy); vym = p5*(vy-avy);
917: #ifndef UPWINDING
918: vxp = vxm = p5*vx;
919: vyp = vym = p5*vy;
920: #endif
922: Bx = D_y(x,psi,i,j);
923: By = - D_x(x,psi,i,j) + By_eq;
924: aBx = PetscAbsScalar(Bx); Bxp = p5*(Bx+aBx); Bxm = p5*(Bx-aBx);
925: aBy = PetscAbsScalar(By); Byp = p5*(By+aBy); Bym = p5*(By-aBy);
926: #ifndef UPWINDING
927: Bxp = Bxm = p5*Bx;
928: Byp = Bym = p5*By;
929: #endif
931: DAVecGetArray(info->da,user->Xold,&xold);
932: dtinv = hxhy/(tsCtx->dt);
933: switch(c) {
935: case 0:
936: /* Lap(phi) - U */
937: *f = (Lapl(x,phi,i,j) - x[j][i].U) * hxhy;
938: break;
940: case 1:
941: /* psi - d_e^2 * Lap(psi) - F */
942: *f = (x[j][i].psi - de2 * Lapl(x,psi,i,j) - x[j][i].F) * hxhy;
943: break;
945: case 2:
946: /* vx * U_x + vy * U_y - (B_x * F_x + B_y * F_y) / d_e^2
947: - nu Lap(U) */
948: *f =
949: ((vxp * (D_xm(x,U,i,j)) + vxm * (D_xp(x,U,i,j)) +
950: vyp * (D_ym(x,U,i,j)) + vym * (D_yp(x,U,i,j))) -
951: (Bxp * (D_xm(x,F,i,j) + F_eq_x) + Bxm * (D_xp(x,F,i,j) + F_eq_x) +
952: Byp * (D_ym(x,F,i,j)) + Bym * (D_yp(x,F,i,j))) * dde2 -
953: nu * Lapl(x,U,i,j)) * hxhy;
954: *f += dtinv*(x[j][i].U-xold[j][i].U);
955: break;
957: case 3:
958: /* vx * F_x + vy * F_y - rho_s^2 * (B_x * U_x + B_y * U_y)
959: - eta * Lap(psi) */
960: *f =
961: ((vxp * (D_xm(x,F,i,j) + F_eq_x) + vxm * (D_xp(x,F,i,j) + F_eq_x) +
962: vyp * (D_ym(x,F,i,j)) + vym * (D_yp(x,F,i,j))) -
963: (Bxp * (D_xm(x,U,i,j)) + Bxm * (D_xp(x,U,i,j)) +
964: Byp * (D_ym(x,U,i,j)) + Bym * (D_yp(x,U,i,j))) * rhos2 -
965: eta * Lapl(x,psi,i,j)) * hxhy;
966: *f += dtinv*(x[j][i].F-xold[j][i].F);
967: break;
968: }
969: DAVecRestoreArray(info->da,user->Xold,&xold);
972: /*
973: Flop count (multiply-adds are counted as 2 operations)
974: */
975: /* PetscLogFlops(84*info->ym*info->xm); FIXME */
976: return(0);
977: }