Actual source code: ex65.c

petsc-dev 2014-08-28
Report Typos and Errors
  1: static char help[] = "2D coupled Allen-Cahn and Cahn-Hilliard equation for constant mobility. Only c_v and eta are considered.\n\
  2: Runtime options include:\n\
  3: -xmin <xmin>\n\
  4: -xmax <xmax>\n\
  5: -ymin <ymin>\n\
  6: -T <T>, where <T> is the end time for the time domain simulation\n\
  7: -dt <dt>,where <dt> is the step size for the numerical integration\n\
  8: -gamma <gamma>\n\
  9: -theta_c <theta_c>\n\n";

 11: /*
 12: ./ex65 -ksp_type fgmres  -snes_atol 1.e-13  -da_refine 6  -VG 10   -pc_type mg -pc_mg_galerkin -log_summary -dt .000001 -mg_coarse_pc_type redundant -mg_coarse_redundant_pc_type svd  -ksp_rtol 1.e-13 -snes_linesearch_type basic -T .0020  -voidgrowth -da_refine 5 -draw_fileds 0,1,2 -dt .0001 -T 1 -da_grid_x 4 -da_grid_y 4 -periodic -snes_rtol 1.e-13 -ksp_atol 1.e-13  -snes_vi_ignore_function_sign -domain 1
 13:  */

 15: #include <petscsnes.h>
 16: #include <petscdmda.h>

 18: typedef struct {
 19:   PetscReal   dt,T; /* Time step and end time */
 20:   PetscReal   dtevent;  /* time scale of radiation events, roughly one event per dtevent */
 21:   PetscInt    maxevents; /* once this number of events is reached no more events are generated */
 22:   PetscReal   initv;
 23:   PetscReal   initeta;
 24:   DM          da1,da1_clone,da2;
 25:   Mat         M;    /* Jacobian matrix */
 26:   Mat         M_0;
 27:   Vec         q,wv,cv,eta,DPsiv,DPsieta,logcv,logcv2,Pv,Pi,Piv;
 28:   Vec         phi1,phi2,Phi2D_V,Sv,Si; /* for twodomain modeling */
 29:   Vec         work1,work2;
 30:   PetscScalar Mv,L,kaeta,kav,Evf,A,B,cv0,Sv_scalar,VG;
 31:   PetscScalar Svr,Sir,cv_eq,ci_eq; /* for twodomain modeling */
 32:   Vec         work3,work4; /* for twodomain modeling*/
 33:   PetscReal   xmin,xmax,ymin,ymax;
 34:   PetscInt    nx;
 35:   PetscBool   graphics;
 36:   PetscBool   periodic;
 37:   PetscBool   lumpedmass;
 38:   PetscBool   radiation; /* Either radiation or void growth */
 39:   PetscInt    domain;
 40:   PetscReal   grain; /* some bogus factor that controls the "strength" of the grain boundaries */
 41:   PetscInt    darefine;
 42:   PetscInt    dagrid;
 43: } AppCtx;

 45: PetscErrorCode GetParams(AppCtx*);
 46: PetscErrorCode SetRandomVectors(AppCtx*,PetscReal);
 47: PetscErrorCode SetVariableBounds(DM,Vec,Vec);
 48: PetscErrorCode SetUpMatrices(AppCtx*);
 49: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
 50: PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
 51: PetscErrorCode SetInitialGuess(Vec,AppCtx*);
 52: PetscErrorCode Update_q(AppCtx*);
 53: PetscErrorCode Update_u(Vec,AppCtx*);
 54: PetscErrorCode DPsi(AppCtx*);
 55: PetscErrorCode Llog(Vec,Vec);
 56: PetscErrorCode CheckRedundancy(SNES,IS,IS*,DM);
 57: PetscErrorCode Phi(AppCtx*);
 58: PetscErrorCode Phi_read(AppCtx*);

 62: int main(int argc, char **argv)
 63: {
 65:   Vec            x,r;  /* solution and residual vectors */
 66:   SNES           snes; /* Nonlinear solver context */
 67:   AppCtx         user; /* Application context */
 68:   Vec            xl,xu; /* Upper and lower bounds on variables */
 69:   Mat            J;
 70:   PetscScalar    t=0.0;
 71:   /* PetscViewer    view_out, view_p, view_q, view_psi, view_mat; */
 72:   PetscReal      PETSC_UNUSED bounds[] = {1000.0,-1000.,0.0,1.0,1000.0,-1000.0,0.0,1.0,1000.0,-1000.0};


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

 77:   /* Get physics and time parameters */
 78:   GetParams(&user);

 80:   if (user.periodic) {
 81:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE, 3, 1,NULL,NULL,&user.da1);
 82:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE, 3, 1,NULL,NULL,&user.da1_clone);
 83:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE, 1, 1,NULL,NULL,&user.da2);
 84:   } else {
 85:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE, 3, 1,NULL,NULL,&user.da1);
 86:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE, 3, 1,NULL,NULL,&user.da1_clone);
 87:     DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,-4,-4,PETSC_DECIDE,PETSC_DECIDE, 1, 1,NULL,NULL,&user.da2);
 88:   }
 89:   /* Set Element type (triangular) */
 90:   DMDASetElementType(user.da1,DMDA_ELEMENT_P1);
 91:   DMDASetElementType(user.da2,DMDA_ELEMENT_P1);

 93:   /* Set x and y coordinates */
 94:   DMDASetUniformCoordinates(user.da1,user.xmin,user.xmax,user.ymin,user.ymax,0.0,0.0);
 95:   DMDASetUniformCoordinates(user.da2,user.xmin,user.xmax,user.ymin,user.ymax,0.0,0.0);
 96:   /* Get global vector x from DM (da1) and duplicate vectors r,xl,xu */
 97:   DMCreateGlobalVector(user.da1,&x);
 98:   VecDuplicate(x,&r);
 99:   VecDuplicate(x,&xl);
100:   VecDuplicate(x,&xu);
101:   VecDuplicate(x,&user.q);

103:   /* Get global vector user->wv from da2 and duplicate other vectors */
104:   DMCreateGlobalVector(user.da2,&user.wv);
105:   VecDuplicate(user.wv,&user.cv);
106:   VecDuplicate(user.wv,&user.eta);
107:   VecDuplicate(user.wv,&user.DPsiv);
108:   VecDuplicate(user.wv,&user.DPsieta);
109:   VecDuplicate(user.wv,&user.Pv);
110:   VecDuplicate(user.wv,&user.Pi);
111:   VecDuplicate(user.wv,&user.Piv);
112:   VecDuplicate(user.wv,&user.logcv);
113:   VecDuplicate(user.wv,&user.logcv2);
114:   VecDuplicate(user.wv,&user.work1);
115:   VecDuplicate(user.wv,&user.work2);
116:   /* for twodomain modeling */
117:   VecDuplicate(user.wv,&user.phi1);
118:   VecDuplicate(user.wv,&user.phi2);
119:   VecDuplicate(user.wv,&user.Phi2D_V);
120:   VecDuplicate(user.wv,&user.Sv);
121:   VecDuplicate(user.wv,&user.Si);
122:   VecDuplicate(user.wv,&user.work3);
123:   VecDuplicate(user.wv,&user.work4);

125:   /* Get Jacobian matrix structure from the da for the entire thing, da1 */
126:   DMSetMatType(user.da1,MATAIJ);
127:   DMCreateMatrix(user.da1,&user.M);
128:   /* Get the (usual) mass matrix structure from da2 */
129:   DMSetMatType(user.da2,MATAIJ);
130:   DMCreateMatrix(user.da2,&user.M_0);
131:   /* Form the jacobian matrix and M_0 */
132:   SetUpMatrices(&user);
133:   SetInitialGuess(x,&user);
134:   MatDuplicate(user.M,MAT_DO_NOT_COPY_VALUES,&J);

136:   SNESCreate(PETSC_COMM_WORLD,&snes);
137:   SNESSetDM(snes,user.da1);

139:   SNESSetFunction(snes,r,FormFunction,(void*)&user);
140:   SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);

142:   SetVariableBounds(user.da1,xl,xu);
143:   SNESVISetVariableBounds(snes,xl,xu);
144:   SNESSetFromOptions(snes);
145:   /* SNESVISetRedundancyCheck(snes,(PetscErrorCode (*)(SNES,IS,IS*,void*))CheckRedundancy,user.da1_clone); */
146:   /*
147:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_out",FILE_MODE_WRITE,&view_out);
148:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_p",FILE_MODE_WRITE,&view_p);
149:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_q",FILE_MODE_WRITE,&view_q);
150:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_psi",FILE_MODE_WRITE,&view_psi);
151:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_mat",FILE_MODE_WRITE,&view_mat);
152:    */
153:   /* PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),5,bounds); */


156:   if (user.graphics) {
157:     /*PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),5,bounds);*/

159:     VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
160:   }


163:   /* multidomain modeling */
164:   if (user.domain) {
165:     switch (user.domain) {
166:     case 1:
167:       Phi(&user);
168:       break;
169:     case 2:
170:       Phi_read(&user);
171:       break;
172:     }
173:   }

175:   while (t<user.T) {

177:     SNESSetFunction(snes,r,FormFunction,(void*)&user);
178:     SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);

180:     if (user.radiation) {
181:       SetRandomVectors(&user,t);
182:     }
183:     DPsi(&user);
184:     /*VecView(user.DPsiv,view_psi);*/
185:     /*VecView(user.DPsieta,view_psi);*/

187:     Update_q(&user);
188:     /*VecView(user.q,view_q);*/
189:     /*MatView(user.M,view_mat);*/

191:     SNESSolve(snes,NULL,x);
192:     /*VecView(x,view_out);*/


195:     if (user.graphics) {
196:       VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
197:     }

199:     PetscInt its;
200:     SNESGetIterationNumber(snes,&its);
201:     PetscPrintf(PETSC_COMM_WORLD,"SNESVI solver converged at t = %5.4g in %d iterations\n",t,its);

203:     Update_u(x,&user);
204:     /*
205:     for (i=0; i < (int)(user.T/a) ; i++) {
206:       if (t/a > i - user.dt/a && t/a < i + user.dt/a) {
207:         sprintf(filename,"output_%f",t);
208:         PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&view);
209:         VecView(user.cv,view);
210:         VecView(user.eta,view);
211:         PetscViewerDestroy(&view);
212:       }
213:     }
214:      */

216:     t = t + user.dt;

218:   }

220:   /*
221:   PetscViewerDestroy(&view_out);
222:   PetscViewerDestroy(&view_p);
223:   PetscViewerDestroy(&view_q);
224:   PetscViewerDestroy(&view_psi);
225:   PetscViewerDestroy(&view_mat);
226:    */
227:   VecDestroy(&x);
228:   VecDestroy(&r);
229:   VecDestroy(&xl);
230:   VecDestroy(&xu);
231:   VecDestroy(&user.q);
232:   VecDestroy(&user.wv);
233:   VecDestroy(&user.cv);
234:   VecDestroy(&user.eta);
235:   VecDestroy(&user.DPsiv);
236:   VecDestroy(&user.DPsieta);
237:   VecDestroy(&user.Pv);
238:   VecDestroy(&user.Pi);
239:   VecDestroy(&user.Piv);
240:   VecDestroy(&user.logcv);
241:   VecDestroy(&user.logcv2);
242:   VecDestroy(&user.work1);
243:   VecDestroy(&user.work2);
244:   /* for twodomain modeling */
245:   VecDestroy(&user.phi1);
246:   VecDestroy(&user.phi2);
247:   VecDestroy(&user.Phi2D_V);
248:   VecDestroy(&user.Sv);
249:   VecDestroy(&user.Si);
250:   VecDestroy(&user.work3);
251:   VecDestroy(&user.work4);

253:   MatDestroy(&user.M);
254:   MatDestroy(&user.M_0);
255:   DMDestroy(&user.da1);
256:   DMDestroy(&user.da1_clone);
257:   DMDestroy(&user.da2);
258:   SNESDestroy(&snes);
259:   PetscFinalize();
260:   return 0;
261: }

265: PetscErrorCode Update_u(Vec X,AppCtx *user)
266: {
268:   PetscInt       i,n;
269:   PetscScalar    *xx,*wv_p,*cv_p,*eta_p;

272:   VecGetLocalSize(user->wv,&n);

274:   VecGetArray(X,&xx);
275:   VecGetArray(user->wv,&wv_p);
276:   VecGetArray(user->cv,&cv_p);
277:   VecGetArray(user->eta,&eta_p);


280:   for (i=0; i<n; i++) {
281:     wv_p[i]  = xx[3*i];
282:     cv_p[i]  = xx[3*i+1];
283:     eta_p[i] = xx[3*i+2];
284:   }
285:   VecRestoreArray(X,&xx);
286:   VecRestoreArray(user->wv,&wv_p);
287:   VecRestoreArray(user->cv,&cv_p);
288:   VecRestoreArray(user->eta,&eta_p);
289:   return(0);
290: }

294: PetscErrorCode Update_q(AppCtx *user)
295: {
297:   PetscScalar    *q_p, *w1, *w2;
298:   PetscInt       n,i;

301:   VecGetArray(user->q,&q_p);
302:   VecGetArray(user->work1,&w1);
303:   VecGetArray(user->work2,&w2);
304:   VecGetLocalSize(user->wv,&n);


307:   VecSet(user->work1,0.0);
308:   if (user->radiation) {
309:     VecAXPY(user->work1,-1.0,user->Pv);
310:   }
311:   if (user->domain) {
312:     VecCopy(user->cv,user->work3);
313:     VecShift(user->work3,-1.0*user->cv_eq);
314:     VecCopy(user->Phi2D_V,user->work4);
315:     VecScale(user->work4,-1.0);
316:     VecShift(user->work4,1.0);
317:     VecPointwiseMult(user->work4,user->work4,user->work3);
318:     VecScale(user->work4,user->Svr);
319:     /* Parameter tuning: user->grain */
320:     /* 5000.0 worked for refinement 2 and 5 */
321:     /*10000.0 worked for refinement 2, 3, 4, 5 */
322:     VecAXPY(user->work1,user->grain,user->work4);
323:   }
324:   VecScale(user->work1,user->dt);
325:   VecAXPY(user->work1,-1.0,user->cv);
326:   MatMult(user->M_0,user->work1,user->work2);
327:   for (i=0; i<n; i++) q_p[3*i]=w2[i];

329:   MatMult(user->M_0,user->DPsiv,user->work1);
330:   for (i=0; i<n; i++) q_p[3*i+1]=w1[i];

332:   VecCopy(user->DPsieta,user->work1);
333:   VecScale(user->work1,user->L*user->dt);
334:   VecAXPY(user->work1,-1.0,user->eta);
335:   if (user->radiation) {
336:     VecAXPY(user->work1,-1.0,user->Piv);
337:   }
338:   MatMult(user->M_0,user->work1,user->work2);
339:   for (i=0; i<n; i++) q_p[3*i+2]=w2[i];

341:   VecRestoreArray(user->q,&q_p);
342:   VecRestoreArray(user->work1,&w1);
343:   VecRestoreArray(user->work2,&w2);
344:   return(0);
345: }

349: PetscErrorCode DPsi(AppCtx *user)
350: {
352:   PetscScalar    Evf=user->Evf,A=user->A,B=user->B,cv0=user->cv0;
353:   PetscScalar    *cv_p,*eta_p,*logcv_p,*logcv2_p,*DPsiv_p,*DPsieta_p;
354:   PetscInt       n,i;

357:   VecGetLocalSize(user->cv,&n);
358:   VecGetArray(user->cv,&cv_p);
359:   VecGetArray(user->eta,&eta_p);
360:   VecGetArray(user->logcv,&logcv_p);
361:   VecGetArray(user->logcv2,&logcv2_p);
362:   VecGetArray(user->DPsiv,&DPsiv_p);
363:   VecGetArray(user->DPsieta,&DPsieta_p);

365:   Llog(user->cv,user->logcv);

367:   VecCopy(user->cv,user->work1);
368:   VecScale(user->work1,-1.0);
369:   VecShift(user->work1,1.0);
370:   Llog(user->work1,user->logcv2);

372:   for (i=0; i<n; i++)
373:   {
374:     DPsiv_p[i] = (eta_p[i]-1.0)*(eta_p[i]-1.0)*(eta_p[i]+1.0)*(eta_p[i]+1.0)*(Evf + logcv_p[i] - logcv2_p[i]) - 2.0*A*(cv_p[i] - cv0)*eta_p[i]*(eta_p[i]+2.0)*(eta_p[i]-1.0)*(eta_p[i]-1.0) + 2.0*B*(cv_p[i] - 1.0)*eta_p[i]*eta_p[i];

376:     DPsieta_p[i] = 4.0*eta_p[i]*(eta_p[i]-1.0)*(eta_p[i]+1.0)*(Evf*cv_p[i] + cv_p[i]*logcv_p[i] + (1.0-cv_p[i])*logcv2_p[i]) - A*(cv_p[i] - cv0)*(cv_p[i] - cv0)*(4.0*eta_p[i]*eta_p[i]*eta_p[i] - 6.0*eta_p[i] + 2.0) + 2.0*B*(cv_p[i]-1.0)*(cv_p[i]-1.0)*eta_p[i];

378:   }

380:   VecRestoreArray(user->cv,&cv_p);
381:   VecRestoreArray(user->eta,&eta_p);
382:   VecGetArray(user->logcv,&logcv_p);
383:   VecGetArray(user->logcv2,&logcv2_p);
384:   VecRestoreArray(user->DPsiv,&DPsiv_p);
385:   VecRestoreArray(user->DPsieta,&DPsieta_p);
386:   return(0);
387: }


392: PetscErrorCode Llog(Vec X, Vec Y)
393: {
395:   PetscScalar    *x,*y;
396:   PetscInt       n,i;

399:   VecGetArray(X,&x);
400:   VecGetArray(Y,&y);
401:   VecGetLocalSize(X,&n);
402:   for (i=0; i<n; i++) {
403:     if (x[i] < 1.0e-12) y[i] = log(1.0e-12);
404:     else y[i] = log(x[i]);
405:   }
406:   return(0);
407: }

411: PetscErrorCode SetInitialGuess(Vec X,AppCtx *user)
412: {
413:   PetscErrorCode    ierr;
414:   PetscInt          n,i,Mda,Nda;
415:   PetscScalar       *xx,*cv_p,*wv_p,*eta_p;
416:   /*PetscViewer      view_out;*/
417:   /* needed for the void growth case */
418:   PetscScalar       xmid,ymid,cv_v=1.0,cv_m=user->Sv_scalar*user->cv0,eta_v=1.0,eta_m=0.0,h,lambda;
419:   PetscInt          nele,nen,idx[3];
420:   const PetscInt    *ele;
421:   PetscScalar       x[3],y[3];
422:   Vec               coords;
423:   const PetscScalar *_coords;
424:   PetscScalar       xwidth = user->xmax - user->xmin;

427:   VecGetLocalSize(X,&n);

429:   if (!user->radiation) {
430:     DMDAGetInfo(user->da2,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
431:     DMGetCoordinatesLocal(user->da2,&coords);
432:     VecGetArrayRead(coords,&_coords);

434:     if (user->periodic) h = (user->xmax-user->xmin)/Mda;
435:     else h = (user->xmax-user->xmin)/(Mda-1.0);

437:     xmid   = (user->xmax + user->xmin)/2.0;
438:     ymid   = (user->ymax + user->ymin)/2.0;
439:     lambda = 4.0*h;

441:     DMDAGetElements(user->da2,&nele,&nen,&ele);
442:     for (i=0; i < nele; i++) {
443:       idx[0] = ele[3*i]; idx[1] = ele[3*i+1]; idx[2] = ele[3*i+2];

445:       x[0] = _coords[2*idx[0]]; y[0] = _coords[2*idx[0]+1];
446:       x[1] = _coords[2*idx[1]]; y[1] = _coords[2*idx[1]+1];
447:       x[2] = _coords[2*idx[2]]; y[2] = _coords[2*idx[2]+1];

449:       PetscInt    k;
450:       PetscScalar vals_DDcv[3],vals_cv[3],vals_eta[3],s,hhr,r;
451:       for (k=0; k < 3; k++) {
452:         /*s = PetscAbs(x[k] - xmid);*/
453:         s = PetscSqrtScalar((x[k] - xmid)*(x[k] - xmid) + (y[k] - ymid)*(y[k] - ymid));
454:         if (s <= xwidth*(5.0/64.0)) {
455:           vals_cv[k]   = cv_v;
456:           vals_eta[k]  = eta_v;
457:           vals_DDcv[k] = 0.0;
458:         } else if (s> xwidth*(5.0/64.0) && s<= xwidth*(7.0/64.0)) {
459:           /*r = (s - xwidth*(6.0/64.0))/(0.5*lambda);*/
460:           r            = (s - xwidth*(6.0/64.0))/(xwidth/64.0);
461:           hhr          = 0.25*(-r*r*r + 3*r + 2);
462:           vals_cv[k]   = cv_m + (1.0 - hhr)*(cv_v - cv_m);
463:           vals_eta[k]  = eta_m + (1.0 - hhr)*(eta_v - eta_m);
464:           vals_DDcv[k] = (cv_v - cv_m)*r*6.0/(lambda*lambda);
465:         } else {
466:           vals_cv[k]   = cv_m;
467:           vals_eta[k]  = eta_m;
468:           vals_DDcv[k] = 0.0;
469:         }
470:       }

472:       VecSetValuesLocal(user->cv,3,idx,vals_cv,INSERT_VALUES);
473:       VecSetValuesLocal(user->eta,3,idx,vals_eta,INSERT_VALUES);
474:       VecSetValuesLocal(user->work2,3,idx,vals_DDcv,INSERT_VALUES);

476:     }
477:     VecAssemblyBegin(user->cv);
478:     VecAssemblyEnd(user->cv);
479:     VecAssemblyBegin(user->eta);
480:     VecAssemblyEnd(user->eta);
481:     VecAssemblyBegin(user->work2);
482:     VecAssemblyEnd(user->work2);

484:     DMDARestoreElements(user->da2,&nele,&nen,&ele);
485:     VecRestoreArrayRead(coords,&_coords);
486:   } else {
487:     /*VecSet(user->cv,user->initv);*/
488:     /*VecSet(user->ci,user->initv);*/
489:     VecSet(user->cv,.05);
490:     VecSet(user->eta,user->initeta);

492:   }
493:   DPsi(user);
494:   VecCopy(user->DPsiv,user->wv);
495:   VecAXPY(user->wv,-2.0*user->kav,user->work2);

497:   VecGetArray(X,&xx);
498:   VecGetArray(user->wv,&wv_p);
499:   VecGetArray(user->cv,&cv_p);
500:   VecGetArray(user->eta,&eta_p);

502:   for (i=0; i<n/3; i++)
503:   {
504:     xx[3*i]  =wv_p[i];
505:     xx[3*i+1]=cv_p[i];
506:     xx[3*i+2]=eta_p[i];
507:   }

509:   /*
510:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_initial",FILE_MODE_WRITE,&view_out);
511:   VecView(user->wv,view_out);
512:   VecView(user->cv,view_out);
513:   VecView(user->eta,view_out);
514:   PetscViewerDestroy(&view_out);
515:    */

517:   VecRestoreArray(X,&xx);
518:   VecRestoreArray(user->wv,&wv_p);
519:   VecRestoreArray(user->cv,&cv_p);
520:   VecRestoreArray(user->eta,&eta_p);
521:   return(0);
522: }

524: typedef struct {
525:   PetscReal dt,x,y,strength;
526: } RandomValues;


531: PetscErrorCode SetRandomVectors(AppCtx *user,PetscReal t)
532: {
533:   PetscErrorCode      ierr;
534:   static RandomValues *randomvalues = 0;
535:   static PetscInt     randindex     = 0,n; /* indicates how far into the randomvalues we have currently used */
536:   static PetscReal    randtime      = 0; /* indicates time of last radiation event */
537:   PetscInt            i,j,M,N,cnt = 0;
538:   PetscInt            xs,ys,xm,ym;

541:   if (!randomvalues) {
542:     PetscViewer viewer;
543:     char        filename[PETSC_MAX_PATH_LEN];
544:     PetscBool   flg;
545:     PetscInt    seed;

547:     PetscOptionsGetInt(NULL,"-random_seed",&seed,&flg);
548:     if (flg) {
549:       sprintf(filename,"ex61.random.%d",(int)seed);
550:     } else {
551:       PetscStrcpy(filename,"ex61.random");
552:     }
553:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
554:     PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);
555:     PetscMalloc1(n,&randomvalues);
556:     PetscViewerBinaryRead(viewer,randomvalues,4*n,PETSC_DOUBLE);
557:     for (i=0; i<n; i++) randomvalues[i].dt = randomvalues[i].dt*user->dtevent;
558:     PetscViewerDestroy(&viewer);
559:   }
560:   user->maxevents = PetscMin(user->maxevents,n);

562:   VecSet(user->Pv,0.0);
563:   DMDAGetInfo(user->da1,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
564:   DMDAGetGhostCorners(user->da1,&xs,&ys,0,&xm,&ym,0);
565:   while (user->maxevents > randindex && randtime + randomvalues[randindex].dt < t + user->dt) {  /* radiation event has occured since last time step */
566:     i = ((PetscInt) (randomvalues[randindex].x*M)) - xs;
567:     j = ((PetscInt) (randomvalues[randindex].y*N)) - ys;
568:     if (i >= 0 && i < xm && j >= 0 && j < ym) { /* point is on this process */

570:       /* need to make sure eta at the given point is not great than .8 */
571:       VecSetValueLocal(user->Pv,i  + xm*(j), randomvalues[randindex].strength*user->VG,INSERT_VALUES);
572:     }
573:     randtime += randomvalues[randindex++].dt;
574:     cnt++;
575:   }
576:   PetscPrintf(PETSC_COMM_WORLD,"Number of radiation events %d\n",cnt);
577:   VecAssemblyBegin(user->Pv);
578:   VecAssemblyEnd(user->Pv);

580:   VecCopy(user->Pv,user->Pi);
581:   VecScale(user->Pi,0.9);
582:   VecPointwiseMult(user->Piv,user->Pi,user->Pv);
583:   return(0);
584: }


589: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ctx)
590: {
592:   AppCtx         *user=(AppCtx*)ctx;

595:   MatMultAdd(user->M,X,user->q,F);
596:   return(0);
597: }

601: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
602: {
604:   AppCtx         *user=(AppCtx*)ctx;

607:   MatCopy(user->M,J,*flg);
608:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
609:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
610:   return(0);
611: }
614: PetscErrorCode SetVariableBounds(DM da,Vec xl,Vec xu)
615: {
617:   PetscScalar    ***l,***u;
618:   PetscInt       xs,xm,ys,ym;
619:   PetscInt       i,j;

622:   DMDAVecGetArrayDOF(da,xl,&l);
623:   DMDAVecGetArrayDOF(da,xu,&u);

625:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

627:   for (j=ys; j < ys+ym; j++) {
628:     for (i=xs; i < xs+xm; i++) {
629:       l[j][i][0] = -PETSC_INFINITY;
630:       l[j][i][1] = 0.0;
631:       l[j][i][2] = 0.0;
632:       u[j][i][0] = PETSC_INFINITY;
633:       u[j][i][1] = 1.0;
634:       u[j][i][2] = 1.0;
635:     }
636:   }

638:   DMDAVecRestoreArrayDOF(da,xl,&l);
639:   DMDAVecRestoreArrayDOF(da,xu,&u);
640:   return(0);
641: }

645: PetscErrorCode GetParams(AppCtx *user)
646: {
648:   PetscBool      flg;

651:   /* Set default parameters */
652:   user->xmin       = 0.0; user->xmax = 128.0;
653:   user->ymin       = 0.0; user->ymax = 128.0;
654:   user->Mv         = 1.0;
655:   user->L          = 1.0;
656:   user->kaeta      = 1.0;
657:   user->kav        = 0.5;
658:   user->Evf        = 9.09;
659:   user->A          = 9.09;
660:   user->B          = 9.09;
661:   user->cv0        = 1.13e-4;
662:   user->Sv_scalar  = 500.0;
663:   user->dt         = 1.0e-5;
664:   user->T          = 1.0e-2;
665:   user->initv      = .00069;
666:   user->initeta    = .0;
667:   user->graphics   = PETSC_FALSE;
668:   user->periodic   = PETSC_FALSE;
669:   user->lumpedmass = PETSC_FALSE;
670:   user->radiation  = PETSC_FALSE;
671:   /* multidomain modeling */
672:   user->domain    = 0;
673:   user->grain     = 5000.0;
674:   user->Svr       = 0.5;
675:   user->Sir       = 0.5;
676:   user->cv_eq     = 6.9e-4;
677:   user->ci_eq     = 6.9e-4;
678:   user->VG        = 100.0;
679:   user->maxevents = 10;

681:   user->dtevent = user->dt;
682:   PetscOptionsReal("-dtevent","Average time between random events\n","None",user->dtevent,&user->dtevent,&flg);


685:   PetscOptionsGetReal(NULL,"-xmin",&user->xmin,&flg);
686:   PetscOptionsGetReal(NULL,"-xmax",&user->xmax,&flg);
687:   PetscOptionsGetReal(NULL,"-T",&user->T,&flg);
688:   PetscOptionsGetReal(NULL,"-dt",&user->dt,&flg);
689:   PetscOptionsBool("-graphics","Contour plot solutions at each timestep\n","None",user->graphics,&user->graphics,&flg);
690:   PetscOptionsBool("-periodic","Use periodic boundary conditions\n","None",user->periodic,&user->periodic,&flg);
691:   PetscOptionsBool("-radiation","Allow radiation\n","None",user->radiation,&user->radiation,&flg);
692:   PetscOptionsBool("-lumpedmass","Use lumped mass matrix\n","None",user->lumpedmass,&user->lumpedmass,&flg);
693:   PetscOptionsInt("-domain","Number of domains (0=one domain, 1=two domains, 2=multidomain\n","None",user->domain,&user->domain,&flg);
694:   PetscOptionsReal("-VG","Maximum increase in vacancy (or interstitial) concentration due to a cascade event","None",user->VG,&user->VG,&flg);
695:   PetscOptionsReal("-grain","Some bogus factor that controls the strength of the grain boundaries, makes the picture more plausible","None",user->grain,&user->grain,&flg);
696:   PetscOptionsInt("-maxevents","Maximum random events allowed\n","None",user->maxevents,&user->maxevents,&flg);
697:   PetscOptionsInt("-da_refine","da refine \n","None",user->darefine,&user->darefine,&flg);
698:   PetscOptionsInt("-da_grid_x","da grid x\n","None",user->dagrid,&user->dagrid,&flg);
699:   return(0);
700: }


705: PetscErrorCode SetUpMatrices(AppCtx *user)
706: {
708:   PetscInt       nele,nen,i,n;
709:   const PetscInt *ele;
710:   PetscScalar    dt=user->dt,h;

712:   PetscInt    idx[3];
713:   PetscScalar eM_0[3][3],eM_2[3][3];
714:   Mat         M  =user->M;
715:   Mat         M_0=user->M_0;
716:   PetscInt    Mda,Nda;


720:   MatGetLocalSize(M,&n,NULL);
721:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

723:   if (Mda!=Nda) printf("Currently different Mda and Nda are not supported");
724:   if (user->periodic) h = (user->xmax-user->xmin)/Mda;
725:   else h = (user->xmax-user->xmin)/(Mda-1.0);

727:   if (user->lumpedmass) {
728:     eM_0[0][0] = eM_0[1][1] = eM_0[2][2] = h*h/6.0;
729:     eM_0[0][1] = eM_0[1][0] = eM_0[0][2] = eM_0[2][0] = eM_0[1][2] = eM_0[2][1] = 0.0;
730:   } else {
731:     eM_0[0][0] = eM_0[1][1] = eM_0[2][2] = h*h/12.0;
732:     eM_0[0][1] = eM_0[0][2] = eM_0[1][0] = eM_0[1][2] = eM_0[2][0] = eM_0[2][1] = h*h/24.0;
733:   }
734:   eM_2[0][0] = 1.0;
735:   eM_2[1][1] = eM_2[2][2] = 0.5;
736:   eM_2[0][1] = eM_2[0][2] = eM_2[1][0]= eM_2[2][0] = -0.5;
737:   eM_2[1][2] = eM_2[2][1] = 0.0;


740:   /* Get local element info */
741:   DMDAGetElements(user->da1,&nele,&nen,&ele);
742:   for (i=0; i < nele; i++) {

744:     idx[0] = ele[3*i]; idx[1] = ele[3*i+1]; idx[2] = ele[3*i+2];

746:     PetscInt    row,cols[6],r,row_M_0,cols3[3];
747:     PetscScalar vals[6],vals_M_0[3],vals3[3];

749:     for (r=0; r<3; r++) {
750:       row_M_0    = idx[r];
751:       vals_M_0[0]=eM_0[r][0];
752:       vals_M_0[1]=eM_0[r][1];
753:       vals_M_0[2]=eM_0[r][2];

755:       MatSetValuesLocal(M_0,1,&row_M_0,3,idx,vals_M_0,ADD_VALUES);

757:       row     = 3*idx[r];
758:       cols[0] = 3*idx[0];     vals[0] = dt*eM_2[r][0]*user->Mv;
759:       cols[1] = 3*idx[1];     vals[1] = dt*eM_2[r][1]*user->Mv;
760:       cols[2] = 3*idx[2];     vals[2] = dt*eM_2[r][2]*user->Mv;
761:       cols[3] = 3*idx[0]+1;   vals[3] = eM_0[r][0];
762:       cols[4] = 3*idx[1]+1;   vals[4] = eM_0[r][1];
763:       cols[5] = 3*idx[2]+1;   vals[5] = eM_0[r][2];

765:       /* Insert values in matrix M for 1st dof */
766:       MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);

768:       row     = 3*idx[r]+1;
769:       cols[0] = 3*idx[0];     vals[0] = -eM_0[r][0];
770:       cols[1] = 3*idx[1];     vals[1] = -eM_0[r][1];
771:       cols[2] = 3*idx[2];     vals[2] = -eM_0[r][2];
772:       cols[3] = 3*idx[0]+1;   vals[3] = 2.0*user->kav*eM_2[r][0];
773:       cols[4] = 3*idx[1]+1;   vals[4] = 2.0*user->kav*eM_2[r][1];
774:       cols[5] = 3*idx[2]+1;   vals[5] = 2.0*user->kav*eM_2[r][2];

776:       /* Insert values in matrix M for 2nd dof */
777:       MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);

779:       row      = 3*idx[r]+2;
780:       cols3[0] = 3*idx[0]+2;   vals3[0] = eM_0[r][0] + user->dt*2.0*user->L*user->kaeta*eM_2[r][0];
781:       cols3[1] = 3*idx[1]+2;   vals3[1] = eM_0[r][1] + user->dt*2.0*user->L*user->kaeta*eM_2[r][1];
782:       cols3[2] = 3*idx[2]+2;   vals3[2] = eM_0[r][2] + user->dt*2.0*user->L*user->kaeta*eM_2[r][2];

784:       MatSetValuesLocal(M,1,&row,3,cols3,vals3,ADD_VALUES);
785:     }
786:   }

788:   MatAssemblyBegin(M_0,MAT_FINAL_ASSEMBLY);
789:   MatAssemblyEnd(M_0,MAT_FINAL_ASSEMBLY);

791:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
792:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

794:   DMDARestoreElements(user->da1,&nele,&nen,&ele);
795:   return(0);
796: }

800: PetscErrorCode CheckRedundancy(SNES snes, IS act, IS *outact, DM da)
801: {
803:   PetscScalar    **uin,**uout;
804:   Vec            UIN, UOUT;
805:   PetscInt       xs,xm,*outindex;
806:   const PetscInt *index;
807:   PetscInt       k,i,l,n,M,cnt=0;

810:   DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
811:   DMGetGlobalVector(da,&UIN);
812:   VecSet(UIN,0.0);
813:   DMGetLocalVector(da,&UOUT);
814:   DMDAVecGetArrayDOF(da,UIN,&uin);
815:   ISGetIndices(act,&index);
816:   ISGetLocalSize(act,&n);
817:   for (k=0; k<n; k++) {
818:     l = index[k]%5;
819:     i = index[k]/5;

821:     uin[i][l]=1.0;
822:   }
823:   printf("Number of active constraints before applying redundancy %d\n",n);
824:   ISRestoreIndices(act,&index);
825:   DMDAVecRestoreArrayDOF(da,UIN,&uin);
826:   DMGlobalToLocalBegin(da,UIN,INSERT_VALUES,UOUT);
827:   DMGlobalToLocalEnd(da,UIN,INSERT_VALUES,UOUT);
828:   DMDAVecGetArrayDOF(da,UOUT,&uout);

830:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

832:   for (i=xs; i < xs+xm;i++) {
833:     if (uout[i-1][1] && uout[i][1] && uout[i+1][1]) uout[i][0] = 1.0;
834:     if (uout[i-1][3] && uout[i][3] && uout[i+1][3]) uout[i][2] = 1.0;
835:   }

837:   for (i=xs; i < xs+xm;i++) {
838:     for (l=0;l<5;l++) {
839:       if (uout[i][l]) cnt++;
840:     }
841:   }

843:   printf("Number of active constraints after applying redundancy %d\n",cnt);


846:   PetscMalloc1(cnt,&outindex);
847:   cnt  = 0;

849:   for (i=xs; i < xs+xm;i++) {
850:     for (l=0;l<5;l++) {
851:       if (uout[i][l]) outindex[cnt++] = 5*(i)+l;
852:     }
853:   }


856:   ISCreateGeneral(PETSC_COMM_WORLD,cnt,outindex,PETSC_OWN_POINTER,outact);
857:   DMDAVecRestoreArrayDOF(da,UOUT,&uout);
858:   DMRestoreGlobalVector(da,&UIN);
859:   DMRestoreLocalVector(da,&UOUT);
860:   return(0);
861: }

865: PetscErrorCode Phi(AppCtx *user)
866: {
867:   PetscErrorCode    ierr;
868:   PetscScalar       xmid, xqu, lambda, h,x[3],y[3];
869:   Vec               coords;
870:   const PetscScalar *_coords;
871:   PetscInt          nele,nen,i,idx[3],Mda,Nda;
872:   const PetscInt    *ele;
873:   /*PetscViewer        view;*/

876:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
877:   DMGetCoordinatesLocal(user->da2,&coords);
878:   VecGetArrayRead(coords,&_coords);

880:   h      = (user->xmax - user->xmin)/Mda;
881:   xmid   = (user->xmin + user->xmax)/2.0;
882:   xqu    = (user->xmin + xmid)/2.0;
883:   lambda = 4.0*h;


886:   DMDAGetElements(user->da2,&nele,&nen,&ele);
887:   for (i=0; i < nele; i++) {
888:     idx[0] = ele[3*i]; idx[1] = ele[3*i+1]; idx[2] = ele[3*i+2];
889:     /*printf("idx[0]=%d,idx[1]=%d,idx[2]=%d\n",idx[0],idx[1],idx[2]);*/

891:     x[0] = _coords[2*idx[0]]; y[0] = _coords[2*idx[0]+1];
892:     x[1] = _coords[2*idx[1]]; y[1] = _coords[2*idx[1]+1];
893:     x[2] = _coords[2*idx[2]]; y[2] = _coords[2*idx[2]+1];

895:     /*printf("x[0]=%f,x[1]=%f,x[2]=%f\n",x[0],x[1],x[2]);*/
896:     /*printf("y[0]=%f,y[1]=%f,y[2]=%f\n",y[0],y[1],y[2]);*/

898:     PetscScalar vals1[3],vals2[3],dist1,dist2,s1,r,hhr,xc1,xc2;
899:     PetscInt    k;

901:     xc1 = user->xmin;
902:     xc2 = xmid;

904:     /*VecSet(user->phi1,0.0);*/
905:     for (k=0;k < 3; k++) {
906:       if (x[k]-xqu > 0) s1 = (x[k] - xqu);
907:       else s1 = -(x[k] - xqu);

909:       if (x[k] - xc1 > 0) dist1 = (x[k] - xc1);
910:       else dist1 = -(x[k] - xc1);

912:       if (x[k] - xc2 > 0) dist2 = (x[k] - xc2);
913:       else dist2 = -(x[k] - xc2);

915:       if (dist1 <= 0.5*lambda) {
916:         r = (x[k]-xc1)/(0.5*lambda);
917:         hhr = 0.25*(-r*r*r + 3.0*r + 2.0);
918:         vals1[k] = hhr;
919:       } else if (dist2 <= 0.5*lambda) {
920:         r = (x[k]-xc2)/(0.5*lambda);
921:         hhr = 0.25*(-r*r*r + 3.0*r + 2.0);
922:         vals1[k] = 1.0 - hhr;
923:       } else if (s1 <= xqu - 2.0*h) {
924:         vals1[k] = 1.0;
925:       } else if ((user->xmax-h)-x[k] < 0.1*h) {
926:       /*else if (abs(x[k]-(user->xmax-h)) < 0.1*h) {*/
927:         vals1[k] = .15625;
928:       } else {
929:         vals1[k] = 0.0;
930:       }
931:     }

933:     VecSetValuesLocal(user->phi1,3,idx,vals1,INSERT_VALUES);

935:     xc1 = xmid;
936:     xc2 = user->xmax;

938:     /*VecSet(user->phi2,0.0);*/
939:     for (k=0;k < 3; k++) {
940:       /*
941:       s1 = abs(x[k] - (xqu+xmid));
942:       dist1 = abs(x[k] - xc1);
943:       dist2 = abs(x[k] - xc2);
944:        */
945:       if (x[k]-(xqu + xmid) > 0) s1 = (x[k] - (xqu + xmid));
946:       else s1 = -(x[k] - (xqu + xmid));

948:       if (x[k] - xc1 > 0) dist1 = (x[k] - xc1);
949:       else dist1 = -(x[k] - xc1);

951:       if (x[k] - xc2 > 0) dist2 = (x[k] - xc2);
952:       else dist2 = -(x[k] - xc2);

954:       if (dist1 <= 0.5*lambda) {
955:         r = (x[k]-xc1)/(0.5*lambda);
956:         hhr = 0.25*(-r*r*r + 3.0*r + 2.0);
957:         vals2[k] = hhr;
958:       } else if (dist2 <= 0.5*lambda) {
959:         r = -(x[k]-xc2)/(0.5*lambda);
960:         hhr = 0.25*(-r*r*r + 3.0*r + 2.0);
961:         vals2[k] = hhr;
962:       } else if (s1 <= xqu - 2.0*h) {
963:         vals2[k] = 1.0;
964:       } else if (x[k]-(user->xmin) < 0.1*h) {
965:         vals2[k] = 0.5;
966:       } else if ((x[k]-(user->xmin+h)) < 0.1*h) {
967:         vals2[k] = .15625;
968:       } else {
969:         vals2[k] = 0.0;
970:       }

972:     }

974:     VecSetValuesLocal(user->phi2,3,idx,vals2,INSERT_VALUES);
975:     /*
976:     for (k=0;k < 3; k++) {
977:       vals_sum[k] = vals1[k]*vals1[k] + vals2[k]*vals2[k];
978:     }
979:      */
980:     /*VecSetValuesLocal(user->Phi2D_V,3,idx,vals_sum,INSERT_VALUES);*/

982:   }

984:   VecAssemblyBegin(user->phi1);
985:   VecAssemblyEnd(user->phi1);
986:   VecAssemblyBegin(user->phi2);
987:   VecAssemblyEnd(user->phi2);

989:   /*
990:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_phi",FILE_MODE_WRITE,&view);
991:   VecView(user->phi1,view);
992:   VecView(user->phi2,view);
993:    */

995:   /*VecView(user->phi1,0);*/
996:   /*VecView(user->phi2,0);*/

998:   VecPointwiseMult(user->phi1,user->phi1,user->phi1);
999:   VecPointwiseMult(user->phi2,user->phi2,user->phi2);
1000:   /*
1001:   VecView(user->phi1,view);
1002:   VecView(user->phi2,view);
1003:    */

1005:   VecCopy(user->phi1,user->Phi2D_V);
1006:   VecAXPY(user->Phi2D_V,1.0,user->phi2);
1007:   /*VecView(user->Phi2D_V,0);*/

1009:   /*VecView(user->Phi2D_V,view);*/
1010:   /*PetscViewerDestroy(&view);*/
1011:   return(0);
1012: }

1016: PetscErrorCode Phi_read(AppCtx *user)
1017: {
1019:   PetscViewer    viewer;
1020:   PetscInt       power;

1023:   power = user->darefine + (PetscInt)(PetscLogScalar(user->dagrid)/PetscLogScalar(2.0));
1024:   switch (power) {
1025:   case 6:
1026:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phi4",FILE_MODE_READ,&viewer);
1027:     VecLoad(user->Phi2D_V,viewer);
1028:     PetscViewerDestroy(&viewer);
1029:     break;
1030:   case 7:
1031:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phi3",FILE_MODE_READ,&viewer);
1032:     VecLoad(user->Phi2D_V,viewer);
1033:     PetscViewerDestroy(&viewer);
1034:     break;
1035:   case 8:
1036:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phi2",FILE_MODE_READ,&viewer);
1037:     VecLoad(user->Phi2D_V,viewer);
1038:     PetscViewerDestroy(&viewer);
1039:     break;
1040:   case 9:
1041:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phi1",FILE_MODE_READ,&viewer);
1042:     VecLoad(user->Phi2D_V,viewer);
1043:     PetscViewerDestroy(&viewer);
1044:     break;
1045:   case 10:
1046:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phi",FILE_MODE_READ,&viewer);
1047:     VecLoad(user->Phi2D_V,viewer);
1048:     PetscViewerDestroy(&viewer);
1049:     break;
1050:   case 11:
1051:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phim1",FILE_MODE_READ,&viewer);
1052:     VecLoad(user->Phi2D_V,viewer);
1053:     PetscViewerDestroy(&viewer);
1054:     break;
1055:   case 12:
1056:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phim2",FILE_MODE_READ,&viewer);
1057:     VecLoad(user->Phi2D_V,viewer);
1058:     PetscViewerDestroy(&viewer);
1059:     break;
1060:   case 13:
1061:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phim3",FILE_MODE_READ,&viewer);
1062:     VecLoad(user->Phi2D_V,viewer);
1063:     PetscViewerDestroy(&viewer);
1064:     break;
1065:   case 14:
1066:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phim4",FILE_MODE_READ,&viewer);
1067:     VecLoad(user->Phi2D_V,viewer);
1068:     PetscViewerDestroy(&viewer);
1069:     break;
1070:   }
1071:   return(0);
1072: }