Actual source code: ex61.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: static char help[] = "2D coupled Allen-Cahn and Cahn-Hilliard equation for constant mobility and triangular elements. Use periodic boundary condidtions.\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:  ./ex61 -ksp_type gmres -snes_vi_monitor   -snes_atol 1.e-11  -da_refine 3  -T 0.1   -ksp_monitor_true_residual -pc_type lu -pc_factor_mat_solver_package superlu -snes_converged_reason -ksp_converged_reason  -ksp_rtol 1.e-9  -snes_linesearch_monitor -VG 10 -draw_fields 1,3,4 -snes_linesearch_type basic

 14: ./ex61 -ksp_type gmres -snes_vi_monitor   -snes_atol 1.e-11  -da_refine 4 -T 0.1   -ksp_monitor_true_residual -pc_type sor -snes_converged_reason -ksp_converged_reason  -ksp_rtol 1.e-9  -snes_linesearch_monitor -VG 10 -draw_fields 1,3,4 -snes_linesearch_type basic

 16: ./ex61 -ksp_type fgmres -snes_vi_monitor   -snes_atol 1.e-11  -da_refine 5 -T 0.1   -ksp_monitor_true_residual -snes_converged_reason -ksp_converged_reason  -ksp_rtol 1.e-9  -snes_linesearch_monitor -VG 10 -draw_fields 1,3,4 -snes_linesearch_type basic -pc_type mg -pc_mg_galerkin

 18: ./ex61 -ksp_type fgmres -snes_vi_monitor   -snes_atol 1.e-11  -da_refine 5 -snes_converged_reason -ksp_converged_reason   -snes_linesearch_monitor -VG 1 -draw_fields 1,3,4  -pc_type mg -pc_mg_galerkin -log_summary -dt .0000000000001 -mg_coarse_pc_type svd  -ksp_monitor_true_residual -ksp_rtol 1.e-9

 20: Movie version
 21: ./ex61 -ksp_type fgmres -snes_vi_monitor   -snes_atol 1.e-11  -da_refine 6 -snes_converged_reason -ksp_converged_reason   -snes_linesearch_monitor -VG 10 -draw_fields 1,3,4  -pc_type mg -pc_mg_galerkin -log_summary -dt .000001 -mg_coarse_pc_type redundant -mg_coarse_redundant_pc_type svd  -ksp_monitor_true_residual -ksp_rtol 1.e-9 -snes_linesearch_type basic -T .0020

 23:  */

 25: /*
 26:    Possible additions to the code. At each iteration count the number of solution elements that are at the upper bound and stop the program if large

 28:    Is the solution at time 0 nonsense?  Looks totally different from first time step. Why does cubic line search at beginning screw it up?

 30:    Add command-line option for constant or degenerate mobility
 31:    Add command-line option for graphics at each time step

 33:    Check time-step business; what should it be? How to check that it is good?
 34:    Make random right hand side forcing function proportional to time step so smaller time steps don't mean more radiation
 35:    How does the multigrid linear solver work now?
 36:    What happens when running with degenerate mobility


 39:  */

 41: #include <petscsnes.h>
 42: #include <petscdmda.h>

 44: typedef struct {
 45:   PetscReal   dt,T; /* Time step and end time */
 46:   PetscReal   dtevent;  /* time scale of radiation events, roughly one event per dtevent */
 47:   PetscInt    maxevents; /* once this number of events is reached no more events are generated */
 48:   PetscReal   initv;    /* initial value of phase variables */
 49:   PetscReal   initeta;
 50:   PetscBool   degenerate;  /* use degenerate mobility */
 51:   PetscReal   smallnumber;
 52:   PetscBool   graphics;
 53:   PetscInt    domain;
 54:   PetscBool   radiation;
 55:   PetscBool   voidgrowth; /* use initial conditions for void growth */
 56:   DM          da1,da2;
 57:   Mat         M;    /* Jacobian matrix */
 58:   Mat         M_0;
 59:   Vec         q,wv,cv,wi,ci,eta,cvi,DPsiv,DPsii,DPsieta,Pv,Pi,Piv,logcv,logci,logcvi,Riv;
 60:   Vec         phi1,phi2,Phi2D_V,Sv,Si; /* for twodomain modeling */
 61:   Vec         work1,work2,work3,work4;
 62:   PetscScalar Dv,Di,Evf,Eif,A,kBT,kav,kai,kaeta,Rsurf,Rbulk,L,VG; /* physics parameters */
 63:   PetscScalar Svr,Sir,cv_eq,ci_eq; /* for twodomain modeling */
 64:   PetscReal   asmallnumber; /* gets added to degenerate mobility */
 65:   PetscReal   xmin,xmax,ymin,ymax;
 66:   PetscInt    Mda, Nda;
 67:   PetscViewer graphicsfile;  /* output of solution at each times step */
 68: } AppCtx;

 70: PetscErrorCode GetParams(AppCtx*);
 71: PetscErrorCode SetRandomVectors(AppCtx*,PetscReal);
 72: PetscErrorCode SetVariableBounds(DM,Vec,Vec);
 73: PetscErrorCode SetUpMatrices(AppCtx*);
 74: PetscErrorCode UpdateMatrices(AppCtx*);
 75: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
 76: PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
 77: PetscErrorCode SetInitialGuess(Vec,AppCtx*);
 78: PetscErrorCode Update_q(AppCtx*);
 79: PetscErrorCode Update_u(Vec,AppCtx*);
 80: PetscErrorCode DPsi(AppCtx*);
 81: PetscErrorCode LaplacianFiniteDifference(AppCtx*);
 82: PetscErrorCode Llog(Vec,Vec);
 83: PetscErrorCode Phi(AppCtx*);
 84: PetscErrorCode Phi_read(AppCtx*);

 88: int main(int argc, char **argv)
 89: {
 90:   PetscErrorCode      ierr;
 91:   Vec                 x,r;  /* Solution and residual vectors */
 92:   SNES                snes; /* Nonlinear solver context */
 93:   AppCtx              user; /* Application context */
 94:   Vec                 xl,xu; /* Upper and lower bounds on variables */
 95:   Mat                 J;
 96:   PetscScalar         t=0.0,normq;
 97:   IS                  inactiveconstraints;
 98:   PetscInt            ninactiveconstraints,N;
 99:   SNESConvergedReason reason;
100:   char                cv_filename[80],eta_filename[80];
101:   /*PetscViewer         view_out, view_cv,view_eta,view_vtk_cv,view_vtk_eta;*/
102:   /*PetscReal           bounds[] = {1000.0,-1000.,0.0,1.0,1000.0,-1000.0,0.0,1.0,1000.0,-1000.0}; */
103:   /*  PetscViewer         view_out, view_q, view_psi, view_mat;*/
104:   /*  PetscViewer         view_rand;*/

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

108:   /* Get physics and time parameters */
109:   GetParams(&user);
110:   /* Create a 1D DA with dof = 5; the whole thing */
111:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX, -3,-3,PETSC_DECIDE,PETSC_DECIDE, 5, 1,NULL,NULL,&user.da1);

113:   /* Create a 1D DA with dof = 1; for individual componentes */
114:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX, -3,-3,PETSC_DECIDE,PETSC_DECIDE, 1, 1,NULL,NULL,&user.da2);


117:   /* Set Element type (triangular) */
118:   DMDASetElementType(user.da1,DMDA_ELEMENT_P1);
119:   DMDASetElementType(user.da2,DMDA_ELEMENT_P1);

121:   /* Set x and y coordinates */
122:   DMDASetUniformCoordinates(user.da1,user.xmin,user.xmax,user.ymin,user.ymax,NULL,NULL);
123:   DMDASetUniformCoordinates(user.da2,user.xmin,user.xmax,user.ymin,user.ymax,NULL,NULL);


126:   /* Get global vector x from DM (da1) and duplicate vectors r,xl,xu */
127:   DMCreateGlobalVector(user.da1,&x);
128:   VecGetSize(x,&N);
129:   VecDuplicate(x,&r);
130:   VecDuplicate(x,&xl);
131:   VecDuplicate(x,&xu);
132:   VecDuplicate(x,&user.q);

134:   /* Get global vector user->wv from da2 and duplicate other vectors */
135:   DMCreateGlobalVector(user.da2,&user.wv);
136:   VecDuplicate(user.wv,&user.cv);
137:   VecDuplicate(user.wv,&user.wi);
138:   VecDuplicate(user.wv,&user.ci);
139:   VecDuplicate(user.wv,&user.eta);
140:   VecDuplicate(user.wv,&user.cvi);
141:   VecDuplicate(user.wv,&user.DPsiv);
142:   VecDuplicate(user.wv,&user.DPsii);
143:   VecDuplicate(user.wv,&user.DPsieta);
144:   VecDuplicate(user.wv,&user.Pv);
145:   VecDuplicate(user.wv,&user.Pi);
146:   VecDuplicate(user.wv,&user.Piv);
147:   VecDuplicate(user.wv,&user.logcv);
148:   VecDuplicate(user.wv,&user.logci);
149:   VecDuplicate(user.wv,&user.logcvi);
150:   VecDuplicate(user.wv,&user.work1);
151:   VecDuplicate(user.wv,&user.work2);
152:   VecDuplicate(user.wv,&user.Riv);
153:   VecDuplicate(user.wv,&user.phi1);
154:   VecDuplicate(user.wv,&user.phi2);
155:   VecDuplicate(user.wv,&user.Phi2D_V);
156:   VecDuplicate(user.wv,&user.Sv);
157:   VecDuplicate(user.wv,&user.Si);
158:   VecDuplicate(user.wv,&user.work3);
159:   VecDuplicate(user.wv,&user.work4);
160:   /* for twodomain modeling */
161:   VecDuplicate(user.wv,&user.phi1);
162:   VecDuplicate(user.wv,&user.phi2);
163:   VecDuplicate(user.wv,&user.Phi2D_V);
164:   VecDuplicate(user.wv,&user.Sv);
165:   VecDuplicate(user.wv,&user.Si);
166:   VecDuplicate(user.wv,&user.work3);
167:   VecDuplicate(user.wv,&user.work4);


170:   /* Get Jacobian matrix structure from the da for the entire thing, da1 */
171:   DMSetMatType(user.da1,MATAIJ);
172:   DMCreateMatrix(user.da1,&user.M);
173:   /* Get the (usual) mass matrix structure from da2 */
174:   DMSetMatType(user.da2,MATAIJ);
175:   DMCreateMatrix(user.da2,&user.M_0);
176:   SetInitialGuess(x,&user);
177:   /* twodomain modeling */
178:   if (user.domain) {
179:     switch (user.domain) {
180:     case 1:
181:       Phi(&user);
182:       break;
183:     case 2:
184:       Phi_read(&user);
185:       break;
186:     }
187:   }

189:   /* Form the jacobian matrix and M_0 */
190:   SetUpMatrices(&user);
191:   MatDuplicate(user.M,MAT_DO_NOT_COPY_VALUES,&J);

193:   SNESCreate(PETSC_COMM_WORLD,&snes);
194:   SNESSetDM(snes,user.da1);

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


200:   SetVariableBounds(user.da1,xl,xu);
201:   SNESVISetVariableBounds(snes,xl,xu);
202:   SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,100,PETSC_DEFAULT);
203:   SNESSetFromOptions(snes);

205:   /*  PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_rand",FILE_MODE_WRITE,&view_rand);
206:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_mat2",FILE_MODE_WRITE,&view_mat);
207:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_q",FILE_MODE_WRITE,&view_q);
208:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_out",FILE_MODE_WRITE,&view_out);
209:    PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_psi",FILE_MODE_WRITE,&view_psi);*/

211:   /*PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_out",FILE_MODE_WRITE,&view_out);

213:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_cv",FILE_MODE_WRITE,&view_cv);
214:    PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_eta",FILE_MODE_WRITE,&view_eta);*/

216:   /* PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),5,bounds); */
217:   if (user.graphics) {
218:     VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
219:   }
220:   /*
221:   if (user.graphicsfile) {
222:     DMView(user.da1,user.graphicsfile);
223:     VecView(x,user.graphicsfile);
224:   }
225:    */
226:   while (t<user.T) {
227:     SNESSetFunction(snes,r,FormFunction,(void*)&user);
228:     SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);

230:     SetRandomVectors(&user,t);
231:     /*    VecView(user.Pv,view_rand);
232:     VecView(user.Pi,view_rand);
233:      VecView(user.Piv,view_rand);*/

235:     DPsi(&user);
236:     /*    VecView(user.DPsiv,view_psi);
237:     VecView(user.DPsii,view_psi);
238:      VecView(user.DPsieta,view_psi);*/

240:     Update_q(&user);

242:     /*    VecView(user.q,view_q);*/
243:     /*  MatView(user.M,view_mat);*/



247:     sprintf(cv_filename,"file_cv_%f.vtk",t);
248:     sprintf(eta_filename,"file_eta_%f.vtk",t);
249:     /*    PetscViewerASCIIOpen(PETSC_COMM_WORLD,cv_filename,&view_vtk_cv);
250:     PetscViewerASCIIOpen(PETSC_COMM_WORLD,eta_filename,&view_vtk_eta);
251:     PetscViewerSetFormat(view_vtk_cv, PETSC_VIEWER_ASCII_VTK);
252:     PetscViewerSetFormat(view_vtk_eta, PETSC_VIEWER_ASCII_VTK);
253:     DMView(user.da2,view_vtk_cv);
254:     DMView(user.da2,view_vtk_eta);
255:     VecView(user.cv,view_cv);
256:     VecView(user.eta,view_eta);
257:     VecView(user.cv,view_vtk_cv);
258:     VecView(user.eta,view_vtk_eta);
259:     PetscViewerDestroy(&view_vtk_cv);
260:      PetscViewerDestroy(&view_vtk_eta);*/


263:     VecNorm(user.q,NORM_2,&normq);
264:     printf("2-norm of q = %14.12f\n",normq);
265:     SNESSolve(snes,NULL,x);
266:     SNESGetConvergedReason(snes,&reason);
267:     if (reason < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Nonlinear solver failed");
268:     SNESVIGetInactiveSet(snes,&inactiveconstraints);
269:     ISGetSize(inactiveconstraints,&ninactiveconstraints);
270:     /* if (ninactiveconstraints < .90*N) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP,"To many active constraints, model has become non-physical"); */

272:     /*    VecView(x,view_out);*/
273:     if (user.graphics) {
274:       VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
275:     }
276:     /*    VecView(x,PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD));*/
277:     PetscInt its;
278:     SNESGetIterationNumber(snes,&its);
279:     PetscPrintf(PETSC_COMM_WORLD,"SNESVI solver converged at t = %g in %d iterations\n",t,its);

281:     Update_u(x,&user);
282:     UpdateMatrices(&user);
283:     t    = t + user.dt;
284:     /*
285:     if (user.graphicsfile) {
286:       VecView(x,user.graphicsfile);
287:     }
288:      */
289:   }

291:   /*  PetscViewerDestroy(&view_rand);
292:   PetscViewerDestroy(&view_mat);
293:   PetscViewerDestroy(&view_q);
294:   PetscViewerDestroy(&view_out);
295:    PetscViewerDestroy(&view_psi);
296:   PetscViewerDestroy(&view_out);

298:   PetscViewerDestroy(&view_cv);
299:    PetscViewerDestroy(&view_eta);*/

301:   if (user.graphicsfile) {
302:     PetscViewerDestroy(&user.graphicsfile);
303:   }

305:   VecDestroy(&x);
306:   VecDestroy(&r);
307:   VecDestroy(&xl);
308:   VecDestroy(&xu);
309:   VecDestroy(&user.q);
310:   VecDestroy(&user.wv);
311:   VecDestroy(&user.cv);
312:   VecDestroy(&user.wi);
313:   VecDestroy(&user.ci);
314:   VecDestroy(&user.eta);
315:   VecDestroy(&user.cvi);
316:   VecDestroy(&user.DPsiv);
317:   VecDestroy(&user.DPsii);
318:   VecDestroy(&user.DPsieta);
319:   VecDestroy(&user.Pv);
320:   VecDestroy(&user.Pi);
321:   VecDestroy(&user.Piv);
322:   VecDestroy(&user.logcv);
323:   VecDestroy(&user.logci);
324:   VecDestroy(&user.logcvi);
325:   VecDestroy(&user.work1);
326:   VecDestroy(&user.work2);
327:   VecDestroy(&user.Riv);
328:   MatDestroy(&user.M);
329:   MatDestroy(&user.M_0);
330:   DMDestroy(&user.da1);
331:   DMDestroy(&user.da2);


334:   PetscFinalize();
335:   return 0;
336: }

340: PetscErrorCode Update_u(Vec X,AppCtx *user)
341: {
343:   PetscInt       i,n;
344:   PetscScalar    *xx,*wv_p,*cv_p,*wi_p,*ci_p,*eta_p;

347:   VecGetLocalSize(user->wv,&n);
348:   VecGetArray(X,&xx);
349:   VecGetArray(user->wv,&wv_p);
350:   VecGetArray(user->cv,&cv_p);
351:   VecGetArray(user->wi,&wi_p);
352:   VecGetArray(user->ci,&ci_p);
353:   VecGetArray(user->eta,&eta_p);


356:   for (i=0; i<n; i++) {
357:     wv_p[i]  = xx[5*i];
358:     cv_p[i]  = xx[5*i+1];
359:     wi_p[i]  = xx[5*i+2];
360:     ci_p[i]  = xx[5*i+3];
361:     eta_p[i] = xx[5*i+4];
362:   }
363:   VecRestoreArray(X,&xx);
364:   VecRestoreArray(user->wv,&wv_p);
365:   VecRestoreArray(user->cv,&cv_p);
366:   VecRestoreArray(user->wi,&wi_p);
367:   VecRestoreArray(user->ci,&ci_p);
368:   VecRestoreArray(user->eta,&eta_p);
369:   return(0);
370: }

374: PetscErrorCode Update_q(AppCtx *user)
375: {
377:   PetscScalar    *q_p,*w1,*w2,max1;
378:   PetscInt       i,n;

381:   VecPointwiseMult(user->Riv,user->eta,user->eta);
382:   VecScale(user->Riv,user->Rsurf);
383:   VecShift(user->Riv,user->Rbulk);
384:   VecPointwiseMult(user->Riv,user->ci,user->Riv);
385:   VecPointwiseMult(user->Riv,user->cv,user->Riv);

387:   VecCopy(user->Riv,user->work1);
388:   if (user->radiation) {
389:     VecAXPY(user->work1,-1.0,user->Pv);
390:   }
391:   if (user->domain) {
392:     VecCopy(user->cv,user->work3);
393:     VecShift(user->work3,-1.0*user->cv_eq);
394:     VecCopy(user->Phi2D_V,user->work4);
395:     VecScale(user->work4,-1.0);
396:     VecShift(user->work4,1.0);
397:     VecPointwiseMult(user->work4,user->work4,user->work3);
398:     VecScale(user->work4,user->Svr);
399:     VecAXPY(user->work1,300.0,user->work4);
400:   }
401:   VecScale(user->work1,user->dt);
402:   VecAXPY(user->work1,-1.0,user->cv);
403:   MatMult(user->M_0,user->work1,user->work2);

405:   VecGetArray(user->q,&q_p);
406:   VecGetArray(user->work1,&w1);
407:   VecGetArray(user->work2,&w2);
408:   VecGetLocalSize(user->wv,&n);
409:   for (i=0; i<n; i++) q_p[5*i]=w2[i];

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

414:   VecCopy(user->Riv,user->work1);
415:   if (user->radiation) {
416:     VecAXPY(user->work1,-1.0,user->Pi);
417:   }
418:   if (user->domain) {
419:     VecCopy(user->ci,user->work3);
420:     VecShift(user->work3,-1.0*user->ci_eq);
421:     VecCopy(user->Phi2D_V,user->work4);
422:     VecScale(user->work4,-1.0);
423:     VecShift(user->work4,1.0);
424:     VecPointwiseMult(user->work4,user->work4,user->work3);
425:     VecScale(user->work4,user->Sir);
426:     VecAXPY(user->work1,300.0,user->work4);
427:   }
428:   VecScale(user->work1,user->dt);
429:   VecAXPY(user->work1,-1.0,user->ci);
430:   MatMult(user->M_0,user->work1,user->work2);
431:   for (i=0; i<n; i++) q_p[5*i+2]=w2[i];

433:   MatMult(user->M_0,user->DPsii,user->work1);
434:   for (i=0; i<n; i++) q_p[5*i+3]=w1[i];

436:   VecCopy(user->DPsieta,user->work1);
437:   VecScale(user->work1,user->L);
438:   VecAXPY(user->work1,-1.0/user->dt,user->eta);
439:   if (user->radiation) {
440:     VecAXPY(user->work1,-1.0,user->Piv);
441:   }
442:   MatMult(user->M_0,user->work1,user->work2);
443:   for (i=0; i<n; i++) q_p[5*i+4]=w2[i];

445:   VecRestoreArray(user->q,&q_p);
446:   VecRestoreArray(user->work1,&w1);
447:   VecRestoreArray(user->work2,&w2);
448:   return(0);
449: }

453: PetscErrorCode DPsi(AppCtx *user)
454: {
456:   PetscScalar    Evf=user->Evf,Eif=user->Eif,kBT=user->kBT,A=user->A;
457:   PetscScalar    *cv_p,*ci_p,*eta_p,*logcv_p,*logci_p,*logcvi_p,*DPsiv_p,*DPsii_p,*DPsieta_p;
458:   PetscInt       n,i;

461:   VecGetLocalSize(user->cv,&n);
462:   VecGetArray(user->cv,&cv_p);
463:   VecGetArray(user->ci,&ci_p);
464:   VecGetArray(user->eta,&eta_p);
465:   VecGetArray(user->logcv,&logcv_p);
466:   VecGetArray(user->logci,&logci_p);
467:   VecGetArray(user->logcvi,&logcvi_p);
468:   VecGetArray(user->DPsiv,&DPsiv_p);
469:   VecGetArray(user->DPsii,&DPsii_p);
470:   VecGetArray(user->DPsieta,&DPsieta_p);

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

474:   Llog(user->ci,user->logci);


477:   VecCopy(user->cv,user->cvi);
478:   VecAXPY(user->cvi,1.0,user->ci);
479:   VecScale(user->cvi,-1.0);
480:   VecShift(user->cvi,1.0);
481:   Llog(user->cvi,user->logcvi);

483:   for (i=0; i<n; i++) {
484:     DPsiv_p[i] = (eta_p[i]-1.0)*(eta_p[i]-1.0)*(Evf + kBT*(logcv_p[i] - logcvi_p[i])) + eta_p[i]*eta_p[i]*2*A*(cv_p[i]-1);
485:     DPsii_p[i] = (eta_p[i]-1.0)*(eta_p[i]-1.0)*(Eif + kBT*(logci_p[i] - logcvi_p[i])) + eta_p[i]*eta_p[i]*2*A*ci_p[i];

487:     DPsieta_p[i] = 2.0*(eta_p[i]-1.0)*(Evf*cv_p[i] + Eif*ci_p[i] + kBT*(cv_p[i]* logcv_p[i] + ci_p[i]* logci_p[i] + (1-cv_p[i]-ci_p[i])*logcvi_p[i])) + 2.0*eta_p[i]*A*((cv_p[i]-1.0)*(cv_p[i]-1.0) + ci_p[i]*ci_p[i]);
488:   }

490:   VecRestoreArray(user->cv,&cv_p);
491:   VecRestoreArray(user->ci,&ci_p);
492:   VecRestoreArray(user->eta,&eta_p);
493:   VecRestoreArray(user->logcv,&logcv_p);
494:   VecRestoreArray(user->logci,&logci_p);
495:   VecRestoreArray(user->logcvi,&logcvi_p);
496:   VecRestoreArray(user->DPsiv,&DPsiv_p);
497:   VecRestoreArray(user->DPsii,&DPsii_p);
498:   VecRestoreArray(user->DPsieta,&DPsieta_p);
499:   return(0);
500: }


505: PetscErrorCode Llog(Vec X, Vec Y)
506: {
508:   PetscScalar    *x,*y;
509:   PetscInt       n,i;

512:   VecGetArray(X,&x);
513:   VecGetArray(Y,&y);
514:   VecGetLocalSize(X,&n);
515:   for (i=0; i<n; i++) {
516:     if (x[i] < 1.0e-12) y[i] = log(1.0e-12);
517:     else y[i] = log(x[i]);
518:   }
519:   return(0);
520: }


525: PetscErrorCode SetInitialGuess(Vec X,AppCtx *user)
526: {
528:   PetscInt       n,i,Mda,Nda;
529:   PetscScalar    *xx,*cv_p,*ci_p,*wv_p,*wi_p,*eta;

531:   /* needed for the void growth case */
532:   PetscScalar       xmid,ymid,cv_v=1.0,cv_m=0.122,ci_v=0.0,ci_m=.00069,eta_v=1.0,eta_m=0.0,h,lambda;
533:   PetscInt          nele,nen,idx[3];
534:   const PetscInt    *ele;
535:   PetscScalar       x[3],y[3];
536:   Vec               coords;
537:   const PetscScalar *_coords;
538:   PetscViewer       view;
539:   PetscScalar       xwidth = user->xmax - user->xmin;

542:   VecGetLocalSize(X,&n);

544:   if (user->voidgrowth) {
545:     DMDAGetInfo(user->da2,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
546:     DMGetCoordinatesLocal(user->da2,&coords);
547:     VecGetArrayRead(coords,&_coords);

549:     h      = (user->xmax-user->xmin)/Mda;
550:     xmid   = (user->xmax + user->xmin)/2.0;
551:     ymid   = (user->ymax + user->ymin)/2.0;
552:     lambda = 4.0*h;

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

558:       x[0] = _coords[2*idx[0]]; y[0] = _coords[2*idx[0]+1];
559:       x[1] = _coords[2*idx[1]]; y[1] = _coords[2*idx[1]+1];
560:       x[2] = _coords[2*idx[2]]; y[2] = _coords[2*idx[2]+1];

562:       PetscInt    k;
563:       PetscScalar vals_cv[3],vals_ci[3],vals_eta[3],s,hhr,r;
564:       for (k=0; k < 3; k++) {
565:         s = PetscSqrtScalar((x[k] - xmid)*(x[k] - xmid) + (y[k] - ymid)*(y[k] - ymid));
566:         if (s < xwidth*(5.0/64.0)) {
567:           vals_cv[k]  = cv_v;
568:           vals_ci[k]  = ci_v;
569:           vals_eta[k] = eta_v;
570:         } else if (s>= xwidth*(5.0/64.0) && s<= xwidth*(7.0/64.0)) {
571:           /* r = (s - xwidth*(6.0/64.0))/(0.5*lambda); */
572:           r           = (s - xwidth*(6.0/64.0))/(xwidth/64.0);
573:           hhr         = 0.25*(-r*r*r + 3*r + 2);
574:           vals_cv[k]  = cv_m + (1.0 - hhr)*(cv_v - cv_m);
575:           vals_ci[k]  = ci_m + (1.0 - hhr)*(ci_v - ci_m);
576:           vals_eta[k] = eta_m + (1.0 - hhr)*(eta_v - eta_m);
577:         } else {
578:           vals_cv[k]  = cv_m;
579:           vals_ci[k]  = ci_m;
580:           vals_eta[k] = eta_m;
581:         }
582:       }
583:       VecSetValuesLocal(user->cv,3,idx,vals_cv,INSERT_VALUES);
584:       VecSetValuesLocal(user->ci,3,idx,vals_ci,INSERT_VALUES);
585:       VecSetValuesLocal(user->eta,3,idx,vals_eta,INSERT_VALUES);
586:     }
587:     VecAssemblyBegin(user->cv);
588:     VecAssemblyEnd(user->cv);
589:     VecAssemblyBegin(user->ci);
590:     VecAssemblyEnd(user->ci);
591:     VecAssemblyBegin(user->eta);
592:     VecAssemblyEnd(user->eta);

594:     DMDARestoreElements(user->da2,&nele,&nen,&ele);
595:     VecRestoreArrayRead(coords,&_coords);

597:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_initial",FILE_MODE_WRITE,&view);
598:     VecView(user->cv,view);
599:     VecView(user->ci,view);
600:     VecView(user->eta,view);
601:     PetscViewerDestroy(&view);
602:   } else {
603:     /* VecSet(user->cv,user->initv); */
604:     /* VecSet(user->ci,user->initv); */
605:     VecSet(user->cv,.05);
606:     VecSet(user->ci,.05);
607:     VecSet(user->eta,user->initeta);
608:   }

610:   DPsi(user);
611:   VecCopy(user->DPsiv,user->wv);
612:   VecCopy(user->DPsii,user->wi);

614:   VecGetArray(X,&xx);
615:   VecGetArray(user->cv,&cv_p);
616:   VecGetArray(user->ci,&ci_p);
617:   VecGetArray(user->wv,&wv_p);
618:   VecGetArray(user->wi,&wi_p);
619:   VecGetArray(user->eta,&eta);
620:   for (i=0; i<n/5; i++) {
621:     xx[5*i]  =wv_p[i];
622:     xx[5*i+1]=cv_p[i];
623:     xx[5*i+2]=wi_p[i];
624:     xx[5*i+3]=ci_p[i];
625:     xx[5*i+4]=eta[i];
626:   }

628:   /* VecView(user->wv,view);
629:   VecView(user->cv,view);
630:   VecView(user->wi,view);
631:   VecView(user->ci,view);
632:    PetscViewerDestroy(&view);*/

634:   VecRestoreArray(X,&xx);
635:   VecRestoreArray(user->cv,&cv_p);
636:   VecRestoreArray(user->ci,&ci_p);
637:   VecRestoreArray(user->wv,&wv_p);
638:   VecRestoreArray(user->wi,&wi_p);
639:   VecRestoreArray(user->eta,&eta);
640:   return(0);
641: }

643: typedef struct {
644:   PetscReal dt,x,y,strength;
645: } RandomValues;


650: PetscErrorCode SetRandomVectors(AppCtx *user,PetscReal t)
651: {
652:   PetscErrorCode      ierr;
653:   static RandomValues *randomvalues = 0;
654:   static PetscInt     randindex     = 0,n; /* indicates how far into the randomvalues we have currently used */
655:   static PetscReal    randtime      = 0; /* indicates time of last radiation event */
656:   PetscInt            i,j,M,N,cnt = 0;
657:   PetscInt            xs,ys,xm,ym;

660:   if (!randomvalues) {
661:     PetscViewer viewer;
662:     char        filename[PETSC_MAX_PATH_LEN];
663:     PetscBool   flg;
664:     PetscInt    seed;

666:     PetscOptionsGetInt(NULL,"-random_seed",&seed,&flg);
667:     if (flg) {
668:       sprintf(filename,"ex61.random.%d",(int)seed);
669:     } else {
670:       PetscStrcpy(filename,"ex61.random");
671:     }
672:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
673:     PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);
674:     PetscMalloc1(n,&randomvalues);
675:     PetscViewerBinaryRead(viewer,randomvalues,4*n,PETSC_DOUBLE);
676:     for (i=0; i<n; i++) randomvalues[i].dt = randomvalues[i].dt*user->dtevent;
677:     PetscViewerDestroy(&viewer);
678:   }
679:   user->maxevents = PetscMin(user->maxevents,n);

681:   VecSet(user->Pv,0.0);
682:   DMDAGetInfo(user->da1,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
683:   DMDAGetGhostCorners(user->da1,&xs,&ys,0,&xm,&ym,0);
684:   while (user->maxevents > randindex && randtime + randomvalues[randindex].dt < t + user->dt) {  /* radiation event has occured since last time step */
685:     i = ((PetscInt) (randomvalues[randindex].x*M)) - xs;
686:     j = ((PetscInt) (randomvalues[randindex].y*N)) - ys;
687:     if (i >= 0 && i < xm && j >= 0 && j < ym) { /* point is on this process */
688:       /* need to make sure eta at the given point is not great than .8 */
689:       VecSetValueLocal(user->Pv,i + 1 + xm*(j + 1), randomvalues[randindex].strength*user->VG,INSERT_VALUES);
690:     }
691:     randtime += randomvalues[randindex++].dt;
692:     cnt++;
693:   }
694:   PetscPrintf(PETSC_COMM_WORLD,"Number of radiation events %d\n",cnt);
695:   VecAssemblyBegin(user->Pv);
696:   VecAssemblyEnd(user->Pv);

698:   VecCopy(user->Pv,user->Pi);
699:   VecScale(user->Pi,0.9);
700:   VecPointwiseMult(user->Piv,user->Pi,user->Pv);
701:   return(0);
702: }

706: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ctx)
707: {
709:   AppCtx         *user=(AppCtx*)ctx;

712:   MatMultAdd(user->M,X,user->q,F);
713:   return(0);
714: }

718: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
719: {
721:   AppCtx         *user=(AppCtx*)ctx;

724:   MatCopy(user->M,J,*flg);
725:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
726:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
727:   return(0);
728: }

732: PetscErrorCode SetVariableBounds(DM da,Vec xl,Vec xu)
733: {
735:   PetscScalar    ***l,***u;
736:   PetscInt       xs,xm,ys,ym;
737:   PetscInt       j,i;

740:   DMDAVecGetArrayDOF(da,xl,&l);
741:   DMDAVecGetArrayDOF(da,xu,&u);

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

745:   for (j=ys; j<ys+ym; j++) {
746:     for (i=xs; i < xs+xm; i++) {
747:       l[j][i][0] = -PETSC_INFINITY;
748:       l[j][i][1] = 0.0;
749:       l[j][i][2] = -PETSC_INFINITY;
750:       l[j][i][3] = 0.0;
751:       l[j][i][4] = 0.0;
752:       u[j][i][0] = PETSC_INFINITY;
753:       u[j][i][1] = 1.0;
754:       u[j][i][2] = PETSC_INFINITY;
755:       u[j][i][3] = 1.0;
756:       u[j][i][4] = 1.0;
757:     }
758:   }


761:   DMDAVecRestoreArrayDOF(da,xl,&l);
762:   DMDAVecRestoreArrayDOF(da,xu,&u);
763:   return(0);
764: }

768: PetscErrorCode GetParams(AppCtx *user)
769: {
771:   PetscBool      flg,graphicsfile = PETSC_FALSE;

774:   /* Set default parameters */
775:   user->xmin  = 0.0; user->xmax = 128.0;
776:   user->ymin  = 0.0; user->ymax = 128.0;
777:   user->Dv    = 1.0;
778:   user->Di    = 4.0;
779:   user->Evf   = 0.8;
780:   user->Eif   = 1.2;
781:   user->A     = 1.0;
782:   user->kBT   = 0.11;
783:   user->kav   = 1.0;
784:   user->kai   = 1.0;
785:   user->kaeta = 1.0;
786:   user->Rsurf = 10.0;
787:   user->Rbulk = 1.0;
788:   user->VG    = 100.0;
789:   user->L     = 10.0;

791:   user->T          = 1.0e-2;
792:   user->dt         = 1.0e-4;
793:   user->initv      = .00069;
794:   user->initeta    = 0.0;
795:   user->degenerate = PETSC_FALSE;
796:   user->maxevents  = 10;
797:   user->graphics   = PETSC_TRUE;

799:   /* multidomain modeling */
800:   user->domain =   1;
801:   user->Svr    = 0.5;
802:   user->Sir    = 0.5;
803:   user->cv_eq  = 6.9e-4;
804:   user->ci_eq  = 6.9e-4;
805:   /* void growth */
806:   user->voidgrowth = PETSC_FALSE;

808:   user->radiation = PETSC_FALSE;

810:   /* degenerate mobility */
811:   user->smallnumber = 1.0e-3;

813:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Coupled Cahn-Hillard/Allen-Cahn Equations","Phasefield");
814:   PetscOptionsInt("-domain","Number of domains (0=one domain, 1=two domains, 2=multidomain\n","None",user->domain,&user->domain,&flg);

816:   PetscOptionsReal("-Dv","Vacancy Diffusivity\n","None",user->Dv,&user->Dv,&flg);
817:   PetscOptionsReal("-Di","Interstitial Diffusivity\n","None",user->Di,&user->Di,&flg);
818:   PetscOptionsReal("-Evf","Vacancy Formation Energy\n","None",user->Evf,&user->Evf,&flg);
819:   PetscOptionsReal("-Eif","Interstitial Formation energy\n","None",user->Eif,&user->Eif,&flg);
820:   PetscOptionsReal("-A","???","None",user->A,&user->A,&flg);
821:   PetscOptionsReal("-kBT","Boltzmann's Constant times the Absolute Temperature","None",user->kBT,&user->kBT,&flg);
822:   PetscOptionsReal("-kav","???","None",user->kav,&user->kav,&flg);
823:   PetscOptionsReal("-kai","???","None",user->kai,&user->kai,&flg);
824:   PetscOptionsReal("-kaeta","???","None",user->kaeta,&user->kaeta,&flg);
825:   PetscOptionsReal("-Rsurf","???","None",user->Rsurf,&user->Rsurf,&flg);
826:   PetscOptionsReal("-Rbulk","???","None",user->Rbulk,&user->Rbulk,&flg);
827:   PetscOptionsReal("-VG","Maximum increase in vacancy (or interstitial) concentration due to a cascade event","None",user->VG,&user->VG,&flg);
828:   PetscOptionsReal("-L","???","None",user->L,&user->L,&flg);

830:   PetscOptionsReal("-initv","Initial solution of Cv and Ci","None",user->initv,&user->initv,&flg);
831:   PetscOptionsReal("-initeta","Initial solution of Eta","None",user->initeta,&user->initeta,&flg);
832:   PetscOptionsBool("-degenerate","Run with degenerate mobility\n","None",user->degenerate,&user->degenerate,&flg);
833:   PetscOptionsReal("-smallnumber","Small number added to degenerate mobility\n","None",user->smallnumber,&user->smallnumber,&flg);

835:   PetscOptionsBool("-voidgrowth","Use initial conditions for void growth\n","None",user->voidgrowth,&user->voidgrowth,&flg);
836:   PetscOptionsBool("-radiation","Use initial conditions for void growth\n","None",user->radiation,&user->radiation,&flg);
837:   PetscOptionsReal("-xmin","Lower X coordinate of domain\n","None",user->xmin,&user->xmin,&flg);
838:   PetscOptionsReal("-xmax","Upper X coordinate of domain\n","None",user->xmax,&user->xmax,&flg);
839:   PetscOptionsReal("-T","Total runtime\n","None",user->T,&user->T,&flg);
840:   PetscOptionsReal("-dt","Time step\n","None",user->dt,&user->dt,&flg);

842:   user->dtevent = user->dt;

844:   PetscOptionsReal("-dtevent","Average time between random events\n","None",user->dtevent,&user->dtevent,&flg);
845:   PetscOptionsInt("-maxevents","Maximum random events allowed\n","None",user->maxevents,&user->maxevents,&flg);

847:   PetscOptionsBool("-graphics","Contour plot solutions at each timestep\n","None",user->graphics,&user->graphics,&flg);
848:   PetscOptionsBool("-graphicsfile","Save solution at each timestep\n","None",graphicsfile,&graphicsfile,&flg);
849:   if (graphicsfile) {
850:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex61.data",FILE_MODE_WRITE,&user->graphicsfile);
851:   }
852:   PetscOptionsEnd();
853:   return(0);
854: }


859: PetscErrorCode SetUpMatrices(AppCtx *user)
860: {
862:   PetscInt       nele,nen,i,n;
863:   const PetscInt *ele;
864:   PetscScalar    dt=user->dt,hx,hy;

866:   PetscInt    idx[3];
867:   PetscScalar eM_0[3][3],eM_2_even[3][3],eM_2_odd[3][3];
868:   PetscScalar cv_sum, ci_sum;
869:   Mat         M  =user->M;
870:   Mat         M_0=user->M_0;
871:   PetscInt    Mda=user->Mda, Nda=user->Nda;
872:   PetscScalar *cv_p,*ci_p;
873:   /* newly added */
874:   Vec cvlocal,cilocal;

877:   /*  MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
878:    MatSetOption(M_0,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);*/

880:   /* new stuff */
881:   DMGetLocalVector(user->da2,&cvlocal);
882:   DMGetLocalVector(user->da2,&cilocal);
883:   DMGlobalToLocalBegin(user->da2,user->cv,INSERT_VALUES,cvlocal);
884:   DMGlobalToLocalEnd(user->da2,user->cv,INSERT_VALUES,cvlocal);
885:   DMGlobalToLocalBegin(user->da2,user->ci,INSERT_VALUES,cilocal);
886:   DMGlobalToLocalEnd(user->da2,user->ci,INSERT_VALUES,cilocal);
887:   /* old stuff */
888:   /*
889:   VecGetArray(user->cv,&cv_p);
890:   VecGetArray(user->ci,&ci_p);
891:    */
892:   /* new stuff */
893:   VecGetArray(cvlocal,&cv_p);
894:   VecGetArray(cilocal,&ci_p);

896:   MatGetLocalSize(M,&n,NULL);
897:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
898:   hx   = (user->xmax-user->xmin)/Mda;
899:   hy   = (user->ymax-user->ymin)/Nda;

901:   eM_0[0][0]=eM_0[1][1]=eM_0[2][2]=hx*hy/12.0;
902:   eM_0[0][1]=eM_0[0][2]=eM_0[1][0]=eM_0[1][2]=eM_0[2][0]=eM_0[2][1]=hx*hy/24.0;

904:   eM_2_odd[0][0] = 1.0;
905:   eM_2_odd[1][1] = eM_2_odd[2][2] = 0.5;
906:   eM_2_odd[0][1] = eM_2_odd[0][2] = eM_2_odd[1][0]= eM_2_odd[2][0] = -0.5;
907:   eM_2_odd[1][2] = eM_2_odd[2][1] = 0.0;

909:   eM_2_even[0][0] = 1.0;
910:   eM_2_even[1][1] = eM_2_even[2][2] = 0.5;
911:   eM_2_even[0][1] = eM_2_even[0][2] = eM_2_even[1][0]= eM_2_even[2][0] = -0.5;
912:   eM_2_even[1][2] = eM_2_even[2][1] = 0.0;

914:   /*  eM_2_even[1][1] = 1.0;
915:   eM_2_even[0][0] = eM_2_even[2][2] = 0.5;
916:   eM_2_even[0][1] = eM_2_even[1][0] = eM_2_even[1][2] = eM_2_even[2][1] = -0.5;
917:   eM_2_even[0][2] = eM_2_even[2][0] = 0.0;
918:    */

920:   /*  for (k=0;k < Mda*Nda*2;k++) { */
921:   DMDAGetElements(user->da1,&nele,&nen,&ele);
922:   for (i=0; i < nele; i++) {
923:     /*
924:     idx[0] = connect[k*3];
925:     idx[1] = connect[k*3+1];
926:     idx[2] = connect[k*3+2];
927:      */
928:     idx[0] = ele[3*i];
929:     idx[1] = ele[3*i+1];
930:     idx[2] = ele[3*i+2];

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

935:     for (r=0; r<3; r++) {
936:       /* row_M_0 = connect[k*3+r]; */
937:       row_M_0 = idx[r];

939:       vals_M_0[0]=eM_0[r][0];
940:       vals_M_0[1]=eM_0[r][1];
941:       vals_M_0[2]=eM_0[r][2];


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

946:       if (user->degenerate) {
947:         cv_sum = (cv_p[idx[0]] + cv_p[idx[1]] + cv_p[idx[2]])*user->Dv/(3.0*user->kBT);
948:         ci_sum = (ci_p[idx[0]] + ci_p[idx[1]] + ci_p[idx[2]])*user->Di/(3.0*user->kBT);
949:       } else {
950:         cv_sum = user->initv*user->Dv/user->kBT;
951:         ci_sum = user->initv*user->Di/user->kBT;
952:       }



956:       row     = 5*idx[r];
957:       cols[0] = 5*idx[0];     vals[0] = dt*eM_2_odd[r][0]*cv_sum;
958:       cols[1] = 5*idx[1];     vals[1] = dt*eM_2_odd[r][1]*cv_sum;
959:       cols[2] = 5*idx[2];     vals[2] = dt*eM_2_odd[r][2]*cv_sum;
960:       cols[3] = 5*idx[0]+1;   vals[3] = eM_0[r][0];
961:       cols[4] = 5*idx[1]+1;   vals[4] = eM_0[r][1];
962:       cols[5] = 5*idx[2]+1;   vals[5] = eM_0[r][2];


965:       MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);


968:       row     = 5*idx[r]+1;
969:       cols[0] = 5*idx[0];     vals[0] = -1.0*eM_0[r][0];
970:       cols[1] = 5*idx[1];     vals[1] = -1.0*eM_0[r][1];
971:       cols[2] = 5*idx[2];     vals[2] = -1.0*eM_0[r][2];
972:       cols[3] = 5*idx[0]+1;   vals[3] =  user->kav*eM_2_odd[r][0];
973:       cols[4] = 5*idx[1]+1;   vals[4] =  user->kav*eM_2_odd[r][1];
974:       cols[5] = 5*idx[2]+1;   vals[5] =  user->kav*eM_2_odd[r][2];

976:       MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);

978:       row     = 5*idx[r]+2;
979:       cols[0] = 5*idx[0]+2;   vals[0] =  dt*eM_2_odd[r][0]*ci_sum;
980:       cols[1] = 5*idx[1]+2;   vals[1] =  dt*eM_2_odd[r][1]*ci_sum;
981:       cols[2] = 5*idx[2]+2;   vals[2] =  dt*eM_2_odd[r][2]*ci_sum;
982:       cols[3] = 5*idx[0]+3;   vals[3] =  eM_0[r][0];
983:       cols[4] = 5*idx[1]+3;   vals[4] =  eM_0[r][1];
984:       cols[5] = 5*idx[2]+3;   vals[5] =  eM_0[r][2];

986:       MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);


989:       row     = 5*idx[r]+3;
990:       cols[0] = 5*idx[0]+2;   vals[0] = -1.0*eM_0[r][0];
991:       cols[1] = 5*idx[1]+2;   vals[1] = -1.0*eM_0[r][1];
992:       cols[2] = 5*idx[2]+2;   vals[2] = -1.0*eM_0[r][2];
993:       cols[3] = 5*idx[0]+3;   vals[3] =  user->kai*eM_2_odd[r][0];
994:       cols[4] = 5*idx[1]+3;   vals[4] =  user->kai*eM_2_odd[r][1];
995:       cols[5] = 5*idx[2]+3;   vals[5] =  user->kai*eM_2_odd[r][2];

997:       MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);


1000:       row = 5*idx[r]+4;
1001:       /*
1002:       cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*dt*eM_2_odd[r][0];
1003:       cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*dt*eM_2_odd[r][1];
1004:       cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*dt*eM_2_odd[r][2];
1005:        */
1006:       cols3[0] = 5*idx[0]+4;   vals3[0] = (eM_0[r][0]/dt + user->L*user->kaeta*eM_2_odd[r][0]);
1007:       cols3[1] = 5*idx[1]+4;   vals3[1] = (eM_0[r][1]/dt + user->L*user->kaeta*eM_2_odd[r][1]);
1008:       cols3[2] = 5*idx[2]+4;   vals3[2] = (eM_0[r][2]/dt + user->L*user->kaeta*eM_2_odd[r][2]);

1010:       MatSetValuesLocal(M,1,&row,3,cols3,vals3,ADD_VALUES);

1012:     }
1013:   }

1015:   /* new */
1016:   VecRestoreArray(cvlocal,&cv_p);
1017:   VecRestoreArray(cilocal,&ci_p);
1018:   DMRestoreLocalVector(user->da2,&cvlocal);
1019:   DMRestoreLocalVector(user->da2,&cilocal);
1020:   /* old */
1021:   /*
1022:   VecRestoreArray(user->cv,&cv_p);
1023:   VecRestoreArray(user->ci,&ci_p);
1024:    */
1025:   MatAssemblyBegin(M_0,MAT_FINAL_ASSEMBLY);
1026:   MatAssemblyEnd(M_0,MAT_FINAL_ASSEMBLY);

1028:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
1029:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

1031:   DMDARestoreElements(user->da1,&nele,&nen,&ele);
1032:   return(0);
1033: }


1038: PetscErrorCode UpdateMatrices(AppCtx *user)
1039: {
1041:   PetscInt       i,n,Mda,Nda,nele,nen;
1042:   const PetscInt *ele;

1044:   PetscInt    idx[3];
1045:   PetscScalar eM_2_odd[3][3],eM_2_even[3][3],h,dt=user->dt;
1046:   Mat         M=user->M;
1047:   PetscScalar *cv_p,*ci_p,cv_sum,ci_sum;
1048:   /* newly added */
1049:   Vec cvlocal,cilocal;

1052:   MatGetLocalSize(M,&n,NULL);

1054:   /* new stuff */
1055:   DMGetLocalVector(user->da2,&cvlocal);
1056:   DMGetLocalVector(user->da2,&cilocal);
1057:   DMGlobalToLocalBegin(user->da2,user->cv,INSERT_VALUES,cvlocal);
1058:   DMGlobalToLocalEnd(user->da2,user->cv,INSERT_VALUES,cvlocal);
1059:   DMGlobalToLocalBegin(user->da2,user->ci,INSERT_VALUES,cilocal);
1060:   DMGlobalToLocalEnd(user->da2,user->ci,INSERT_VALUES,cilocal);
1061:   /* new stuff */
1062:   VecGetArray(cvlocal,&cv_p);
1063:   VecGetArray(cilocal,&ci_p);

1065:   /* old stuff */
1066:   /*
1067:   VecGetArray(user->cv,&cv_p);
1068:   VecGetArray(user->ci,&ci_p);
1069:    */
1070:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);



1074:   h = (user->xmax-user->xmin)/Mda;

1076:   DMDAGetElements(user->da1,&nele,&nen,&ele);

1078:   for (i=0; i < nele; i++) {
1079:     /*
1080:     idx[0] = connect[k*3];
1081:     idx[1] = connect[k*3+1];
1082:     idx[2] = connect[k*3+2];
1083:      */
1084:     idx[0] = ele[3*i];
1085:     idx[1] = ele[3*i+1];
1086:     idx[2] = ele[3*i+2];

1088:     PetscInt    r,row,cols[3];
1089:     PetscScalar vals[3];
1090:     for (r=0; r<3; r++) {
1091:       row     = 5*idx[r];
1092:       cols[0] = 5*idx[0];     vals[0] = 0.0;
1093:       cols[1] = 5*idx[1];     vals[1] = 0.0;
1094:       cols[2] = 5*idx[2];     vals[2] = 0.0;

1096:       /* Insert values in matrix M for 1st dof */
1097:       MatSetValuesLocal(M,1,&row,3,cols,vals,INSERT_VALUES);

1099:       row     = 5*idx[r]+2;
1100:       cols[0] = 5*idx[0]+2;   vals[0] = 0.0;
1101:       cols[1] = 5*idx[1]+2;   vals[1] = 0.0;
1102:       cols[2] = 5*idx[2]+2;   vals[2] = 0.0;

1104:       /* Insert values in matrix M for 3nd dof */
1105:       MatSetValuesLocal(M,1,&row,3,cols,vals,INSERT_VALUES);
1106:     }
1107:   }

1109:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
1110:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

1112:   eM_2_odd[0][0] = 1.0;
1113:   eM_2_odd[1][1] = eM_2_odd[2][2] = 0.5;
1114:   eM_2_odd[0][1] = eM_2_odd[0][2] = eM_2_odd[1][0]= eM_2_odd[2][0] = -0.5;
1115:   eM_2_odd[1][2] = eM_2_odd[2][1] = 0.0;

1117:   eM_2_even[0][0] = 1.0;
1118:   eM_2_even[1][1] = eM_2_even[2][2] = 0.5;
1119:   eM_2_even[0][1] = eM_2_even[0][2] = eM_2_even[1][0]= eM_2_even[2][0] = -0.5;
1120:   eM_2_even[1][2] = eM_2_even[2][1] = 0.0;

1122:   /*
1123:   eM_2_even[1][1] = 1.0;
1124:   eM_2_even[0][0] = eM_2_even[2][2] = 0.5;
1125:   eM_2_even[0][1] = eM_2_even[1][0] = eM_2_even[1][2] = eM_2_even[2][1] = -0.5;
1126:   eM_2_even[0][2] = eM_2_even[2][0] = 0.0;
1127:    */

1129:   /* Get local element info */
1130:   /* for (k=0;k < Mda*Nda*2;k++) { */
1131:   for (i=0; i < nele; i++) {
1132:     /*
1133:       idx[0] = connect[k*3];
1134:       idx[1] = connect[k*3+1];
1135:       idx[2] = connect[k*3+2];
1136:      */
1137:     idx[0] = ele[3*i];
1138:     idx[1] = ele[3*i+1];
1139:     idx[2] = ele[3*i+2];

1141:     PetscInt    row,cols[3],r;
1142:     PetscScalar vals[3];

1144:     for (r=0; r<3; r++) {

1146:       if (user->degenerate) {
1147:         printf("smallnumber = %14.12f\n",user->smallnumber);
1148:         cv_sum = (user->smallnumber + cv_p[idx[0]] + cv_p[idx[1]] + cv_p[idx[2]])*user->Dv/(3.0*user->kBT);
1149:         ci_sum = (user->smallnumber + ci_p[idx[0]] + ci_p[idx[1]] + ci_p[idx[2]])*user->Di/(3.0*user->kBT);
1150:       } else {
1151:         cv_sum = user->initv*user->Dv/(user->kBT);
1152:         ci_sum = user->initv*user->Di/user->kBT;
1153:       }

1155:       row     = 5*idx[r];
1156:       cols[0] = 5*idx[0];     vals[0] = dt*eM_2_odd[r][0]*cv_sum;
1157:       cols[1] = 5*idx[1];     vals[1] = dt*eM_2_odd[r][1]*cv_sum;
1158:       cols[2] = 5*idx[2];     vals[2] = dt*eM_2_odd[r][2]*cv_sum;

1160:       /* Insert values in matrix M for 1st dof */
1161:       MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);

1163:       row     = 5*idx[r]+2;
1164:       cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_odd[r][0]*ci_sum;
1165:       cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_odd[r][1]*ci_sum;
1166:       cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_odd[r][2]*ci_sum;

1168:       MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);
1169:     }
1170:   }

1172:   /* new stuff */
1173:   VecRestoreArray(cvlocal,&cv_p);
1174:   VecRestoreArray(cilocal,&ci_p);
1175:   DMRestoreLocalVector(user->da2,&cvlocal);
1176:   DMRestoreLocalVector(user->da2,&cilocal);
1177:   /* old stuff */
1178:   /*
1179:   VecRestoreArray(user->cv,&cv_p);
1180:   VecRestoreArray(user->ci,&ci_p);
1181:    */

1183:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
1184:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
1185:   return(0);
1186: }



1192: PetscErrorCode Phi(AppCtx *user)
1193: {
1194:   PetscErrorCode    ierr;
1195:   PetscScalar       xmid, xqu, lambda, h,x[3],y[3];
1196:   Vec               coords;
1197:   const PetscScalar *_coords;
1198:   PetscInt          nele,nen,i,idx[3],Mda,Nda;
1199:   const PetscInt    *ele;
1200:   PetscViewer       view;

1203:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
1204:   DMGetCoordinatesLocal(user->da2,&coords);
1205:   VecGetArrayRead(coords,&_coords);

1207:   h      = (user->xmax - user->xmin)/Mda;
1208:   xmid   = (user->xmin + user->xmax)/2.0;
1209:   xqu    = (user->xmin + xmid)/2.0;
1210:   lambda = 4.0*h;


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

1218:     x[0] = _coords[2*idx[0]]; y[0] = _coords[2*idx[0]+1];
1219:     x[1] = _coords[2*idx[1]]; y[1] = _coords[2*idx[1]+1];
1220:     x[2] = _coords[2*idx[2]]; y[2] = _coords[2*idx[2]+1];

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

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

1228:     xc1 = user->xmin;
1229:     xc2 = xmid;

1231:     /*VecSet(user->phi1,0.0);*/
1232:     for (k=0; k < 3; k++) {
1233:       if (x[k]-xqu > 0) s1 = (x[k] - xqu);
1234:       else s1 = -(x[k] - xqu);

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

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

1242:       if (dist1 <= 0.5*lambda) {
1243:         r        = (x[k]-xc1)/(0.5*lambda);
1244:         hhr      = 0.25*(-r*r*r + 3.0*r + 2.0);
1245:         vals1[k] = hhr;
1246:       } else if (dist2 <= 0.5*lambda) {
1247:         r        = (x[k]-xc2)/(0.5*lambda);
1248:         hhr      = 0.25*(-r*r*r + 3.0*r + 2.0);
1249:         vals1[k] = 1.0 - hhr;
1250:       } else if (s1 <= xqu - 2.0*h) vals1[k] = 1.0;
1251:       /*else if (abs(x[k]-(user->xmax-h)) < 0.1*h) {*/
1252:       else if ((user->xmax-h)-x[k] < 0.1*h) vals1[k] = .15625;
1253:       else vals1[k] = 0.0;
1254:     }

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

1258:     xc1 = xmid;
1259:     xc2 = user->xmax;

1261:     /*VecSet(user->phi2,0.0);*/
1262:     for (k=0;k < 3; k++) {
1263:       /*
1264:       s1 = abs(x[k] - (xqu+xmid));
1265:       dist1 = abs(x[k] - xc1);
1266:       dist2 = abs(x[k] - xc2);
1267:        */
1268:       if (x[k]-(xqu + xmid) > 0) s1 = (x[k] - (xqu + xmid));
1269:       else s1 = -(x[k] - (xqu + xmid));

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

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

1277:       if (dist1 <= 0.5*lambda) {
1278:         r        = (x[k]-xc1)/(0.5*lambda);
1279:         hhr      = 0.25*(-r*r*r + 3.0*r + 2.0);
1280:         vals2[k] = hhr;
1281:       } else if (dist2 <= 0.5*lambda) {
1282:         r        = -(x[k]-xc2)/(0.5*lambda);
1283:         hhr      = 0.25*(-r*r*r + 3.0*r + 2.0);
1284:         vals2[k] = hhr;
1285:       } else if (s1 <= xqu - 2.0*h)           vals2[k] = 1.0;
1286:       else if (x[k]-(user->xmin) < 0.1*h)     vals2[k] = 0.5;
1287:       else if ((x[k]-(user->xmin+h)) < 0.1*h) vals2[k] = .15625;
1288:       else vals2[k] = 0.0;

1290:     }

1292:     VecSetValuesLocal(user->phi2,3,idx,vals2,INSERT_VALUES);
1293:     /*
1294:     for (k=0;k < 3; k++) {
1295:       vals_sum[k] = vals1[k]*vals1[k] + vals2[k]*vals2[k];
1296:     }
1297:      */
1298:     /*VecSetValuesLocal(user->Phi2D_V,3,idx,vals_sum,INSERT_VALUES);*/

1300:   }

1302:   VecAssemblyBegin(user->phi1);
1303:   VecAssemblyEnd(user->phi1);
1304:   VecAssemblyBegin(user->phi2);
1305:   VecAssemblyEnd(user->phi2);
1306:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_phi",FILE_MODE_WRITE,&view);

1308:   VecView(user->phi1,view);
1309:   VecView(user->phi2,view);


1312:   /*VecView(user->phi1,0);*/
1313:   /*VecView(user->phi2,0);*/

1315:   VecPointwiseMult(user->phi1,user->phi1,user->phi1);
1316:   VecPointwiseMult(user->phi2,user->phi2,user->phi2);
1317:   VecView(user->phi1,view);
1318:   VecView(user->phi2,view);

1320:   VecCopy(user->phi1,user->Phi2D_V);
1321:   VecAXPY(user->Phi2D_V,1.0,user->phi2);
1322:   /*VecView(user->Phi2D_V,0);*/

1324:   VecView(user->Phi2D_V,view);
1325:   PetscViewerDestroy(&view);
1326:   /*  VecNorm(user->Phi2D_V,NORM_INFINITY,&max1);*/
1327:   /*VecMin(user->Phi2D_V,&loc1,&min1);*/
1328:   /*printf("norm phi = %f, min phi = %f\n",max1,min1);*/
1329:   return(0);
1330: }

1334: PetscErrorCode Phi_read(AppCtx *user)
1335: {
1337:   PetscReal      *values;
1338:   PetscViewer    viewer;

1341:   VecGetArray(user->Phi2D_V,&values);
1342:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"phi3",FILE_MODE_READ,&viewer);
1343:   PetscViewerBinaryRead(viewer,values,16384,PETSC_DOUBLE);
1344:   PetscViewerDestroy(&viewer);
1345:   return(0);
1346: }