Actual source code: ex60.c

petsc-3.4.5 2014-06-29
  1: static char help[] = "2D coupled Allen-Cahn and Cahn-Hilliard equation for constant mobility and triangular elements.\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:  ./ex60 -ksp_type fgmres -pc_type mg  -snes_vi_monitor   -snes_atol 1.e-11  -da_grid_x 72 -da_grid_y 72 -ksp_rtol 1.e-8  -T 0.1  -VG 100 -pc_type lu -ksp_monitor_true_residual -pc_factor_mat_solver_package superlu -snes_converged_reason -ksp_converged_reason -pc_type sor  -ksp_rtol 1.e-9  -snes_linesearch_monitor -VG 10 -draw_fields 1,3,4 -snes_monitor_solution

 14:  */

 16: /*
 17:    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

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

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


 28:  */

 30:  #include petscsnes.h
 31:  #include petscdmda.h

 33: typedef struct {
 34:   PetscReal   dt,T; /* Time step and end time */
 35:   DM          da1,da2;
 36:   Mat         M;    /* Jacobian matrix */
 37:   Mat         M_0;
 38:   Vec         q,wv,cv,wi,ci,eta,cvi,DPsiv,DPsii,DPsieta,Pv,Pi,Piv,logcv,logci,logcvi,Rr,Riv;
 39:   Vec         work1,work2,work3,work4;
 40:   PetscScalar Dv,Di,Evf,Eif,A,kBT,kav,kai,kaeta,Rsurf,Rbulk,L,P_casc,VG; /* physics parameters */
 41:   PetscReal   xmin,xmax,ymin,ymax;
 42:   PetscInt    Mda, Nda;
 43: } AppCtx;

 45: PetscErrorCode GetParams(AppCtx*);
 46: PetscErrorCode SetRandomVectors(AppCtx*);
 47: PetscErrorCode SetVariableBounds(DM,Vec,Vec);
 48: PetscErrorCode SetUpMatrices(AppCtx*);
 49: PetscErrorCode UpdateMatrices(AppCtx*);
 50: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
 51: PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 52: PetscErrorCode SetInitialGuess(Vec,AppCtx*);
 53: PetscErrorCode Update_q(AppCtx*);
 54: PetscErrorCode Update_u(Vec,AppCtx*);
 55: PetscErrorCode DPsi(AppCtx*);
 56: PetscErrorCode LaplacianFiniteDifference(AppCtx*);
 57: PetscErrorCode Llog(Vec,Vec);
 60: int main(int argc, char **argv)
 61: {
 63:   Vec            x,r;  /* Solution and residual vectors */
 64:   SNES           snes; /* Nonlinear solver context */
 65:   AppCtx         user; /* Application context */
 66:   Vec            xl,xu; /* Upper and lower bounds on variables */
 67:   Mat            J;
 68:   PetscScalar    t=0.0;
 69:   PetscViewer    view_out, view_q, view_psi, view_mat;
 70:   PetscViewer    view_rand;
 71:   IS             inactiveconstraints;
 72:   PetscInt       ninactiveconstraints,N;

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

 76:   /* Get physics and time parameters */
 77:   GetParams(&user);
 78:   /* Create a 1D DA with dof = 5; the whole thing */
 79:   DMDACreate2d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX, -3,-3,PETSC_DECIDE,PETSC_DECIDE, 5, 1,NULL,NULL,&user.da1);

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


 85:   /* Set Element type (triangular) */
 86:   DMDASetElementType(user.da1,DMDA_ELEMENT_P1);
 87:   DMDASetElementType(user.da2,DMDA_ELEMENT_P1);

 89:   /* Set x and y coordinates */
 90:   DMDASetUniformCoordinates(user.da1,user.xmin,user.xmax,user.ymin,user.ymax,NULL,NULL);
 91:   DMDASetUniformCoordinates(user.da2,user.xmin,user.xmax,user.ymin,user.ymax,NULL,NULL);
 92:   /* Get global vector x from DM (da1) and duplicate vectors r,xl,xu */
 93:   DMCreateGlobalVector(user.da1,&x);
 94:   VecGetSize(x,&N);
 95:   VecDuplicate(x,&r);
 96:   VecDuplicate(x,&xl);
 97:   VecDuplicate(x,&xu);
 98:   VecDuplicate(x,&user.q);

100:   /* Get global vector user->wv from da2 and duplicate other vectors */
101:   DMCreateGlobalVector(user.da2,&user.wv);
102:   VecDuplicate(user.wv,&user.cv);
103:   VecDuplicate(user.wv,&user.wi);
104:   VecDuplicate(user.wv,&user.ci);
105:   VecDuplicate(user.wv,&user.eta);
106:   VecDuplicate(user.wv,&user.cvi);
107:   VecDuplicate(user.wv,&user.DPsiv);
108:   VecDuplicate(user.wv,&user.DPsii);
109:   VecDuplicate(user.wv,&user.DPsieta);
110:   VecDuplicate(user.wv,&user.Pv);
111:   VecDuplicate(user.wv,&user.Pi);
112:   VecDuplicate(user.wv,&user.Piv);
113:   VecDuplicate(user.wv,&user.logcv);
114:   VecDuplicate(user.wv,&user.logci);
115:   VecDuplicate(user.wv,&user.logcvi);
116:   VecDuplicate(user.wv,&user.work1);
117:   VecDuplicate(user.wv,&user.work2);
118:   VecDuplicate(user.wv,&user.Rr);
119:   VecDuplicate(user.wv,&user.Riv);


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

131:   SNESCreate(PETSC_COMM_WORLD,&snes);
132:   SNESSetDM(snes,user.da1);

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

137:   SetVariableBounds(user.da1,xl,xu);
138:   SNESVISetVariableBounds(snes,xl,xu);
139:   SNESSetFromOptions(snes);

141:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_rand",FILE_MODE_WRITE,&view_rand);
142:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_mat",FILE_MODE_WRITE,&view_mat);
143:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_q",FILE_MODE_WRITE,&view_q);
144:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_out",FILE_MODE_WRITE,&view_out);
145:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_psi",FILE_MODE_WRITE,&view_psi);

147:   while (t<user.T) {
148:     SNESSetFunction(snes,r,FormFunction,(void*)&user);
149:     SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);

151:     SetRandomVectors(&user);
152:     /*    VecView(user.Pv,view_rand);
153:     VecView(user.Pi,view_rand);
154:      VecView(user.Piv,view_rand);*/

156:     DPsi(&user);
157:     /*    VecView(user.DPsiv,view_psi);
158:     VecView(user.DPsii,view_psi);
159:      VecView(user.DPsieta,view_psi);*/

161:     Update_q(&user);
162:     /*    VecView(user.q,view_q);
163:      MatView(user.M,view_mat);*/
164:     SNESSolve(snes,NULL,x);
165:     SNESVIGetInactiveSet(snes,&inactiveconstraints);
166:     ISGetSize(inactiveconstraints,&ninactiveconstraints);
167:     /* if (ninactiveconstraints < .90*N) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP,"To many active constraints, model has become non-physical"); */

169:     /*    VecView(x,view_out);*/
170:     VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
171:     PetscInt its;
172:     SNESGetIterationNumber(snes,&its);
173:     PetscPrintf(PETSC_COMM_WORLD,"SNESVI solver converged at t = %5.4f in %d iterations\n",t,its);

175:     Update_u(x,&user);
176:     UpdateMatrices(&user);
177:     t    = t + user.dt;
178:   }

180:   PetscViewerDestroy(&view_rand);
181:   PetscViewerDestroy(&view_mat);
182:   PetscViewerDestroy(&view_q);
183:   PetscViewerDestroy(&view_out);
184:   PetscViewerDestroy(&view_psi);

186:   VecDestroy(&x);
187:   VecDestroy(&r);
188:   VecDestroy(&xl);
189:   VecDestroy(&xu);
190:   VecDestroy(&user.q);
191:   VecDestroy(&user.wv);
192:   VecDestroy(&user.cv);
193:   VecDestroy(&user.wi);
194:   VecDestroy(&user.ci);
195:   VecDestroy(&user.eta);
196:   VecDestroy(&user.cvi);
197:   VecDestroy(&user.DPsiv);
198:   VecDestroy(&user.DPsii);
199:   VecDestroy(&user.DPsieta);
200:   VecDestroy(&user.Pv);
201:   VecDestroy(&user.Pi);
202:   VecDestroy(&user.Piv);
203:   VecDestroy(&user.logcv);
204:   VecDestroy(&user.logci);
205:   VecDestroy(&user.logcvi);
206:   VecDestroy(&user.work1);
207:   VecDestroy(&user.work2);
208:   VecDestroy(&user.Rr);
209:   VecDestroy(&user.Riv);
210:   MatDestroy(&user.M);
211:   MatDestroy(&user.M_0);
212:   DMDestroy(&user.da1);
213:   DMDestroy(&user.da2);

215:   PetscFinalize();
216:   return 0;
217: }

221: PetscErrorCode Update_u(Vec X,AppCtx *user)
222: {
224:   PetscInt       i,n;
225:   PetscScalar    *xx,*wv_p,*cv_p,*wi_p,*ci_p,*eta_p;

228:   VecGetLocalSize(user->wv,&n);
229:   VecGetArray(X,&xx);
230:   VecGetArray(user->wv,&wv_p);
231:   VecGetArray(user->cv,&cv_p);
232:   VecGetArray(user->wi,&wi_p);
233:   VecGetArray(user->ci,&ci_p);
234:   VecGetArray(user->eta,&eta_p);


237:   for (i=0; i<n; i++) {
238:     wv_p[i]  = xx[5*i];
239:     cv_p[i]  = xx[5*i+1];
240:     wi_p[i]  = xx[5*i+2];
241:     ci_p[i]  = xx[5*i+3];
242:     eta_p[i] = xx[5*i+4];
243:   }
244:   VecRestoreArray(X,&xx);
245:   VecRestoreArray(user->wv,&wv_p);
246:   VecRestoreArray(user->cv,&cv_p);
247:   VecRestoreArray(user->wi,&wi_p);
248:   VecRestoreArray(user->ci,&ci_p);
249:   VecRestoreArray(user->eta,&eta_p);
250:   return(0);
251: }

255: PetscErrorCode Update_q(AppCtx *user)
256: {
258:   PetscScalar    *q_p,*w1,*w2;
259:   PetscInt       i,n;

262:   VecPointwiseMult(user->Rr,user->eta,user->eta);
263:   VecScale(user->Rr,user->Rsurf);
264:   VecShift(user->Rr,user->Rbulk);
265:   VecPointwiseMult(user->Riv,user->cv,user->ci);
266:   VecPointwiseMult(user->Riv,user->Rr,user->Riv);

268:   VecGetArray(user->q,&q_p);
269:   VecGetArray(user->work1,&w1);
270:   VecGetArray(user->work2,&w2);

272:   VecCopy(user->cv,user->work1);
273:   VecAXPY(user->work1,1.0,user->Pv);
274:   VecScale(user->work1,-1.0);
275:   MatMult(user->M_0,user->work1,user->work2);
276:   VecGetLocalSize(user->work1,&n);

278:   for (i=0; i<n; i++) q_p[5*i]=w2[i];

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

283:   VecCopy(user->ci,user->work1);
284:   VecAXPY(user->work1,1.0,user->Pi);
285:   VecScale(user->work1,-1.0);
286:   MatMult(user->M_0,user->work1,user->work2);
287:   for (i=0; i<n; i++) q_p[5*i+2]=w2[i];

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

292:   VecCopy(user->eta,user->work1);
293:   VecScale(user->work1,-1.0/user->dt);
294:   VecAXPY(user->work1,user->L,user->DPsieta);
295:   VecAXPY(user->work1,-1.0,user->Piv);
296:   MatMult(user->M_0,user->work1,user->work2);
297:   for (i=0; i<n; i++) q_p[5*i+4]=w2[i];

299:   VecRestoreArray(user->q,&q_p);
300:   VecRestoreArray(user->work1,&w1);
301:   VecRestoreArray(user->work2,&w2);
302:   return(0);
303: }

307: PetscErrorCode DPsi(AppCtx *user)
308: {
310:   PetscScalar    Evf=user->Evf,Eif=user->Eif,kBT=user->kBT,A=user->A;
311:   PetscScalar    *cv_p,*ci_p,*eta_p,*logcv_p,*logci_p,*logcvi_p,*DPsiv_p,*DPsii_p,*DPsieta_p;
312:   PetscInt       n,i;

315:   VecGetLocalSize(user->cv,&n);
316:   VecGetArray(user->cv,&cv_p);
317:   VecGetArray(user->ci,&ci_p);
318:   VecGetArray(user->eta,&eta_p);
319:   VecGetArray(user->logcv,&logcv_p);
320:   VecGetArray(user->logci,&logci_p);
321:   VecGetArray(user->logcvi,&logcvi_p);
322:   VecGetArray(user->DPsiv,&DPsiv_p);
323:   VecGetArray(user->DPsii,&DPsii_p);
324:   VecGetArray(user->DPsieta,&DPsieta_p);

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

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


331:   VecCopy(user->cv,user->cvi);
332:   VecAXPY(user->cvi,1.0,user->ci);
333:   VecScale(user->cvi,-1.0);
334:   VecShift(user->cvi,1.0);
335:   Llog(user->cvi,user->logcvi);

337:   for (i=0; i<n; i++)
338:   {
339:     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);
340:     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];

342:     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]);
343:   }

345:   VecRestoreArray(user->cv,&cv_p);
346:   VecRestoreArray(user->ci,&ci_p);
347:   VecRestoreArray(user->eta,&eta_p);
348:   VecRestoreArray(user->logcv,&logcv_p);
349:   VecRestoreArray(user->logci,&logci_p);
350:   VecRestoreArray(user->logcvi,&logcvi_p);
351:   VecRestoreArray(user->DPsiv,&DPsiv_p);
352:   VecRestoreArray(user->DPsii,&DPsii_p);
353:   VecRestoreArray(user->DPsieta,&DPsieta_p);
354:   return(0);
355: }


360: PetscErrorCode Llog(Vec X, Vec Y)
361: {
363:   PetscScalar    *x,*y;
364:   PetscInt       n,i;

367:   VecGetArray(X,&x);
368:   VecGetArray(Y,&y);
369:   VecGetLocalSize(X,&n);
370:   for (i=0; i<n; i++) {
371:     if (x[i] < 1.0e-12) y[i] = log(1.0e-12);
372:     else y[i] = log(x[i]);
373:   }
374:   return(0);
375: }


380: PetscErrorCode SetInitialGuess(Vec X,AppCtx *user)
381: {
383:   PetscInt       n,i;
384:   PetscScalar    *xx,*cv_p,*ci_p,*wv_p,*wi_p;
385:   PetscViewer    view;
386:   PetscScalar    initv = .00069;

389:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_initial",FILE_MODE_WRITE,&view);
390:   VecGetLocalSize(X,&n);

392:   VecSet(user->cv,initv);
393:   VecSet(user->ci,initv);
394:   VecSet(user->eta,0.0);

396:   DPsi(user);
397:   VecCopy(user->DPsiv,user->wv);
398:   VecCopy(user->DPsii,user->wi);

400:   VecGetArray(X,&xx);
401:   VecGetArray(user->cv,&cv_p);
402:   VecGetArray(user->ci,&ci_p);
403:   VecGetArray(user->wv,&wv_p);
404:   VecGetArray(user->wi,&wi_p);
405:   for (i=0; i<n/5; i++) {
406:     xx[5*i]  =wv_p[i];
407:     xx[5*i+1]=cv_p[i];
408:     xx[5*i+2]=wi_p[i];
409:     xx[5*i+3]=ci_p[i];
410:     xx[5*i+4]=0.0;
411:   }

413:   VecView(user->wv,view);
414:   VecView(user->cv,view);
415:   VecView(user->wi,view);
416:   VecView(user->ci,view);
417:   PetscViewerDestroy(&view);

419:   VecRestoreArray(X,&xx);
420:   VecRestoreArray(user->cv,&cv_p);
421:   VecRestoreArray(user->ci,&ci_p);
422:   VecRestoreArray(user->wv,&wv_p);
423:   VecRestoreArray(user->wi,&wi_p);
424:   return(0);
425: }

429: PetscErrorCode SetRandomVectors(AppCtx *user)
430: {
432:   PetscInt       i,n,count=0;
433:   PetscScalar    *w1,*w2,*Pv_p,*eta_p;

436:   VecSetRandom(user->work1,NULL);
437:   VecSetRandom(user->work2,NULL);
438:   VecGetArray(user->work1,&w1);
439:   VecGetArray(user->work2,&w2);
440:   VecGetArray(user->Pv,&Pv_p);
441:   VecGetArray(user->eta,&eta_p);
442:   VecGetLocalSize(user->work1,&n);
443:   for (i=0; i<n; i++) {

445:     if (eta_p[i]>=0.8 || w1[i]>user->P_casc) Pv_p[i]=0;
446:     else {
447:       Pv_p[i]=w2[i]*user->VG;
448:       count  =count+1;
449:     }

451:   }

453:   VecCopy(user->Pv,user->Pi);
454:   VecScale(user->Pi,0.9);
455:   VecPointwiseMult(user->Piv,user->Pi,user->Pv);
456:   VecRestoreArray(user->work1,&w1);
457:   VecRestoreArray(user->work2,&w2);
458:   VecRestoreArray(user->Pv,&Pv_p);
459:   VecRestoreArray(user->eta,&eta_p);
460:   return(0);

462: }

466: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ctx)
467: {
469:   AppCtx         *user=(AppCtx*)ctx;

472:   MatMultAdd(user->M,X,user->q,F);
473:   return(0);
474: }

478: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flg,void *ctx)
479: {
481:   AppCtx         *user=(AppCtx*)ctx;

484:   *flg = SAME_NONZERO_PATTERN;
485:   MatCopy(user->M,*J,*flg);
486:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
487:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
488:   return(0);
489: }
492: PetscErrorCode SetVariableBounds(DM da,Vec xl,Vec xu)
493: {
495:   PetscScalar    ***l,***u;
496:   PetscInt       xs,xm,ys,ym;
497:   PetscInt       j,i;

500:   DMDAVecGetArrayDOF(da,xl,&l);
501:   DMDAVecGetArrayDOF(da,xu,&u);

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

505:   for (j=ys; j<ys+ym; j++) {
506:     for (i=xs; i < xs+xm; i++) {
507:       l[j][i][0] = -SNES_VI_INF;
508:       l[j][i][1] = 0.0;
509:       l[j][i][2] = -SNES_VI_INF;
510:       l[j][i][3] = 0.0;
511:       l[j][i][4] = 0.0;
512:       u[j][i][0] = SNES_VI_INF;
513:       u[j][i][1] = 1.0;
514:       u[j][i][2] = SNES_VI_INF;
515:       u[j][i][3] = 1.0;
516:       u[j][i][4] = 1.0;
517:     }
518:   }


521:   DMDAVecRestoreArrayDOF(da,xl,&l);
522:   DMDAVecRestoreArrayDOF(da,xu,&u);
523:   return(0);
524: }

528: PetscErrorCode GetParams(AppCtx *user)
529: {
531:   PetscBool      flg;

534:   /* Set default parameters */
535:   user->xmin  = 0.0; user->xmax = 1.0;
536:   user->ymin  = 0.0; user->ymax = 1.0;
537:   user->Dv    = 1.0; user->Di=4.0;
538:   user->Evf   = 0.8; user->Eif = 1.2;
539:   user->A     = 1.0;
540:   user->kBT   = 0.11;
541:   user->kav   = 1.0; user->kai = 1.0; user->kaeta = 1.0;
542:   user->Rsurf = 10.0; user->Rbulk = 1.0;
543:   user->L     = 10.0; user->P_casc = 0.05;
544:   user->T     = 1.0e-2;    user->dt = 1.0e-4;
545:   user->VG    = 100.0;

547:   PetscOptionsGetReal(NULL,"-xmin",&user->xmin,&flg);
548:   PetscOptionsGetReal(NULL,"-xmax",&user->xmax,&flg);
549:   PetscOptionsGetReal(NULL,"-T",&user->T,&flg);
550:   PetscOptionsGetReal(NULL,"-dt",&user->dt,&flg);
551:   PetscOptionsGetReal(NULL,"-P_casc",&user->P_casc,&flg);
552:   PetscOptionsGetReal(NULL,"-VG",&user->VG,&flg);
553:   return(0);
554: }


559: PetscErrorCode SetUpMatrices(AppCtx *user)
560: {
562:   PetscInt       nele,nen,i,j,n;
563:   const PetscInt *ele;
564:   PetscScalar    dt=user->dt,hx,hy;

566:   PetscInt    idx[3],*nodes, *connect, k;
567:   PetscScalar eM_0[3][3],eM_2_even[3][3],eM_2_odd[3][3];
568:   PetscScalar cv_sum, ci_sum;
569:   Mat         M  =user->M;
570:   Mat         M_0=user->M_0;
571:   PetscInt    Mda=user->Mda, Nda=user->Nda, ld, rd, ru, lu;
572:   PetscScalar *cv_p,*ci_p;

575:   /*  MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
576:    MatSetOption(M_0,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);*/

578:   /* Create the mass matrix M_0 */
579:   VecGetArray(user->cv,&cv_p);
580:   VecGetArray(user->ci,&ci_p);
581:   MatGetLocalSize(M,&n,NULL);
582:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

584:   PetscMalloc((Mda+1)*(Nda+1)*sizeof(PetscInt),&nodes);
585:   PetscMalloc(Mda*Nda*2*3*sizeof(PetscInt),&connect);
586:   hx   = (user->xmax-user->xmin)/Mda;
587:   hy   = (user->ymax-user->ymin)/Nda;
588:   for (j=0; j < Nda; j++) {
589:     for (i=0; i < Mda; i++) nodes[j*(Mda+1)+i] = j*Mda+i;
590:     nodes[j*(Mda+1)+Mda] = j*Mda;
591:   }

593:   for (i=0; i < Mda; i++) nodes[Nda*(Mda+1)+i] = i;

595:   nodes[Nda*(Mda+1)+Mda] = 0;

597:   k = 0;
598:   for (j=0; j<Nda; j++) {
599:     for (i=0; i<Mda; i++) {

601:       /* ld = nodes[j][i]; */
602:       ld = nodes[j*(Mda+1)+i];
603:       /* rd = nodes[j+1][i]; */
604:       rd = nodes[(j+1)*(Mda+1)+i];
605:       /* ru = nodes[j+1][i+1]; */
606:       ru = nodes[(j+1)*(Mda+1)+i+1];
607:       /* lu = nodes[j][i+1]; */
608:       lu = nodes[j*(Mda+1)+i+1];

610:       /* connect[k][0]=ld; */
611:       connect[k*6]=ld;
612:       /* connect[k][1]=lu; */
613:       connect[k*6+1]=lu;
614:       /* connect[k][2]=ru; */
615:       connect[k*6+2]=rd;
616:       connect[k*6+3]=lu;
617:       connect[k*6+4]=ru;
618:       connect[k*6+5]=rd;

620:       k = k+1;
621:     }
622:   }


625:   eM_0[0][0]=eM_0[1][1]=eM_0[2][2]=hx*hy/12.0;
626:   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;

628:   eM_2_odd[0][0] = 1.0;
629:   eM_2_odd[1][1] = eM_2_odd[2][2] = 0.5;
630:   eM_2_odd[0][1] = eM_2_odd[0][2] = eM_2_odd[1][0]= eM_2_odd[2][0] = -0.5;
631:   eM_2_odd[1][2] = eM_2_odd[2][1] = 0.0;

633:   eM_2_even[1][1] = 1.0;
634:   eM_2_even[0][0] = eM_2_even[2][2] = 0.5;
635:   eM_2_even[0][1] = eM_2_even[1][0] = eM_2_even[1][2] = eM_2_even[2][1] = -0.5;
636:   eM_2_even[0][2] = eM_2_even[2][0] = 0.0;


639:   for (k=0; k < Mda*Nda*2; k++) {
640:     idx[0] = connect[k*3];
641:     idx[1] = connect[k*3+1];
642:     idx[2] = connect[k*3+2];

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

647:     for (r=0; r<3; r++) {
648:       row_M_0 = connect[k*3+r];

650:       vals_M_0[0]=eM_0[r][0];
651:       vals_M_0[1]=eM_0[r][1];
652:       vals_M_0[2]=eM_0[r][2];


655:       MatSetValues(M_0,1,&row_M_0,3,idx,vals_M_0,ADD_VALUES);

657:       /* cv_sum = (cv_p[idx[0]] + cv_p[idx[1]] + cv_p[idx[2]])*user->Dv/(3.0*user->kBT); */
658:       /* ci_sum = (ci_p[idx[0]] + ci_p[idx[1]] + ci_p[idx[2]])*user->Di/(3.0*user->kBT); */
659:       cv_sum = .0000069*user->Dv/user->kBT;
660:       ci_sum = .0000069*user->Di/user->kBT;

662:       if (k%2 == 0) {
663:         row     = 5*idx[r];
664:         cols[0] = 5*idx[0];     vals[0] = dt*eM_2_odd[r][0]*cv_sum;
665:         cols[1] = 5*idx[1];     vals[1] = dt*eM_2_odd[r][1]*cv_sum;
666:         cols[2] = 5*idx[2];     vals[2] = dt*eM_2_odd[r][2]*cv_sum;
667:         cols[3] = 5*idx[0]+1;   vals[3] = eM_0[r][0];
668:         cols[4] = 5*idx[1]+1;   vals[4] = eM_0[r][1];
669:         cols[5] = 5*idx[2]+1;   vals[5] = eM_0[r][2];

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

673:         row     = 5*idx[r]+1;
674:         cols[0] = 5*idx[0];     vals[0] = -1.0*eM_0[r][0];
675:         cols[1] = 5*idx[1];     vals[1] = -1.0*eM_0[r][1];
676:         cols[2] = 5*idx[2];     vals[2] = -1.0*eM_0[r][2];
677:         cols[3] = 5*idx[0]+1;   vals[3] =  user->kav*eM_2_odd[r][0];
678:         cols[4] = 5*idx[1]+1;   vals[4] =  user->kav*eM_2_odd[r][1];
679:         cols[5] = 5*idx[2]+1;   vals[5] =  user->kav*eM_2_odd[r][2];

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

683:         row     = 5*idx[r]+2;
684:         cols[0] = 5*idx[0]+2;   vals[0] =  dt*eM_2_odd[r][0]*ci_sum;
685:         cols[1] = 5*idx[1]+2;   vals[1] =  dt*eM_2_odd[r][1]*ci_sum;
686:         cols[2] = 5*idx[2]+2;   vals[2] =  dt*eM_2_odd[r][2]*ci_sum;
687:         cols[3] = 5*idx[0]+3;   vals[3] =  eM_0[r][0];
688:         cols[4] = 5*idx[1]+3;   vals[4] =  eM_0[r][1];
689:         cols[5] = 5*idx[2]+3;   vals[5] =  eM_0[r][2];

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

693:         row     = 5*idx[r]+3;
694:         cols[0] = 5*idx[0]+2;   vals[0] = -1.0*eM_0[r][0];
695:         cols[1] = 5*idx[1]+2;   vals[1] = -1.0*eM_0[r][1];
696:         cols[2] = 5*idx[2]+2;   vals[2] = -1.0*eM_0[r][2];
697:         cols[3] = 5*idx[0]+3;   vals[3] =  user->kai*eM_2_odd[r][0];
698:         cols[4] = 5*idx[1]+3;   vals[4] =  user->kai*eM_2_odd[r][1];
699:         cols[5] = 5*idx[2]+3;   vals[5] =  user->kai*eM_2_odd[r][2];

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

703:         row = 5*idx[r]+4;
704:         /*
705:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*dt*eM_2_odd[r][0];
706:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*dt*eM_2_odd[r][1];
707:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*dt*eM_2_odd[r][2];
708:          */
709:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*eM_2_odd[r][0];
710:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*eM_2_odd[r][1];
711:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*eM_2_odd[r][2];

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


716:       } else {


719:         row     = 5*idx[r];
720:         cols[0] = 5*idx[0];     vals[0] = dt*eM_2_even[r][0]*cv_sum;
721:         cols[1] = 5*idx[1];     vals[1] = dt*eM_2_even[r][1]*cv_sum;
722:         cols[2] = 5*idx[2];     vals[2] = dt*eM_2_even[r][2]*cv_sum;
723:         cols[3] = 5*idx[0]+1;   vals[3] = eM_0[r][0];
724:         cols[4] = 5*idx[1]+1;   vals[4] = eM_0[r][1];
725:         cols[5] = 5*idx[2]+1;   vals[5] = eM_0[r][2];

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

729:         row     = 5*idx[r]+1;
730:         cols[0] = 5*idx[0];     vals[0] = -1.0*eM_0[r][0];
731:         cols[1] = 5*idx[1];     vals[1] = -1.0*eM_0[r][1];
732:         cols[2] = 5*idx[2];     vals[2] = -1.0*eM_0[r][2];
733:         cols[3] = 5*idx[0]+1;   vals[3] =  user->kav*eM_2_even[r][0];
734:         cols[4] = 5*idx[1]+1;   vals[4] =  user->kav*eM_2_even[r][1];
735:         cols[5] = 5*idx[2]+1;   vals[5] =  user->kav*eM_2_even[r][2];

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

739:         row     = 5*idx[r]+2;
740:         cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_even[r][0]*ci_sum;
741:         cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_even[r][1]*ci_sum;
742:         cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_even[r][2]*ci_sum;
743:         cols[3] = 5*idx[0]+3;   vals[3] = eM_0[r][0];
744:         cols[4] = 5*idx[1]+3;   vals[4] = eM_0[r][1];
745:         cols[5] = 5*idx[2]+3;   vals[5] = eM_0[r][2];

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

749:         row     = 5*idx[r]+3;
750:         cols[0] = 5*idx[0]+2;   vals[0] = -1.0*eM_0[r][0];
751:         cols[1] = 5*idx[1]+2;   vals[1] = -1.0*eM_0[r][1];
752:         cols[2] = 5*idx[2]+2;   vals[2] = -1.0*eM_0[r][2];
753:         cols[3] = 5*idx[0]+3;   vals[3] =  user->kai*eM_2_even[r][0];
754:         cols[4] = 5*idx[1]+3;   vals[4] =  user->kai*eM_2_even[r][1];
755:         cols[5] = 5*idx[2]+3;   vals[5] =  user->kai*eM_2_even[r][2];

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

759:         row = 5*idx[r]+4;
760:         /*
761:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*dt*eM_2_even[r][0];
762:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*dt*eM_2_even[r][1];
763:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*dt*eM_2_even[r][2];
764:          */
765:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*eM_2_even[r][0];
766:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*eM_2_even[r][1];
767:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*eM_2_even[r][2];

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

771:       }

773:     }
774:   }

776:   PetscFree(nodes);
777:   PetscFree(connect);

779:   VecRestoreArray(user->cv,&cv_p);
780:   VecRestoreArray(user->ci,&ci_p);
781:   MatAssemblyBegin(M_0,MAT_FINAL_ASSEMBLY);
782:   MatAssemblyEnd(M_0,MAT_FINAL_ASSEMBLY);

784:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
785:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

787:   DMDARestoreElements(user->da1,&nele,&nen,&ele);
788:   return(0);
789: }


794: PetscErrorCode UpdateMatrices(AppCtx *user)
795: {
797:   PetscInt       i,j,n,Mda,Nda;

799:   PetscInt    idx[3],*nodes,*connect,k;
800:   PetscInt    ld,rd,lu,ru;
801:   PetscScalar eM_2_odd[3][3],eM_2_even[3][3],h,dt=user->dt;
802:   Mat         M=user->M;
803:   PetscScalar *cv_p,*ci_p,cv_sum,ci_sum;

806:   /* Create the mass matrix M_0 */
807:   MatGetLocalSize(M,&n,NULL);
808:   VecGetArray(user->cv,&cv_p);
809:   VecGetArray(user->ci,&ci_p);
810:   DMDAGetInfo(user->da1,NULL,&Mda,&Nda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

812:   PetscMalloc((Mda+1)*(Nda+1)*sizeof(PetscInt),&nodes);
813:   PetscMalloc(Mda*Nda*2*3*sizeof(PetscInt),&connect);

815:   h = 100.0/Mda;

817:   for (j=0; j < Nda; j++) {
818:     for (i=0; i < Mda; i++) nodes[j*(Mda+1)+i] = j*Mda+i;
819:     nodes[j*(Mda+1)+Mda] = j*Mda;
820:   }
821:   for (i=0; i < Mda; i++) nodes[Nda*(Mda+1)+i]=i;
822:   nodes[Nda*(Mda+1)+Mda]=0;


825:   k = 0;
826:   for (j=0; j<Nda; j++) {
827:     for (i=0; i<Mda; i++) {
828:       ld = nodes[j*(Mda+1)+i];
829:       rd = nodes[(j+1)*(Mda+1)+i];
830:       ru = nodes[(j+1)*(Mda+1)+i+1];
831:       lu = nodes[j*(Mda+1)+i+1];

833:       connect[k*6]   = ld;
834:       connect[k*6+1] = lu;
835:       connect[k*6+2] = rd;
836:       connect[k*6+3] = lu;
837:       connect[k*6+4] = ru;
838:       connect[k*6+5] = rd;

840:       k = k+1;
841:     }
842:   }

844:   for (k=0; k < Mda*Nda*2; k++) {
845:     idx[0] = connect[k*3];
846:     idx[1] = connect[k*3+1];
847:     idx[2] = connect[k*3+2];

849:     PetscInt    r,row,cols[3];
850:     PetscScalar vals[3];
851:     for (r=0; r<3; r++) {
852:       row     = 5*idx[r];
853:       cols[0] = 5*idx[0];     vals[0] = 0.0;
854:       cols[1] = 5*idx[1];     vals[1] = 0.0;
855:       cols[2] = 5*idx[2];     vals[2] = 0.0;

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

860:       row     = 5*idx[r]+2;
861:       cols[0] = 5*idx[0]+2;   vals[0] = 0.0;
862:       cols[1] = 5*idx[1]+2;   vals[1] = 0.0;
863:       cols[2] = 5*idx[2]+2;   vals[2] = 0.0;

865:       /* Insert values in matrix M for 3nd dof */
866:       MatSetValuesLocal(M,1,&row,3,cols,vals,INSERT_VALUES);
867:     }
868:   }

870:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
871:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

873:   eM_2_odd[0][0] = 1.0;
874:   eM_2_odd[1][1] = eM_2_odd[2][2] = 0.5;
875:   eM_2_odd[0][1] = eM_2_odd[0][2] = eM_2_odd[1][0]= eM_2_odd[2][0] = -0.5;
876:   eM_2_odd[1][2] = eM_2_odd[2][1] = 0.0;

878:   eM_2_even[1][1] = 1.0;
879:   eM_2_even[0][0] = eM_2_even[2][2] = 0.5;
880:   eM_2_even[0][1] = eM_2_even[1][0] = eM_2_even[1][2] = eM_2_even[2][1] = -0.5;
881:   eM_2_even[0][2] = eM_2_even[2][0] = 0.0;


884:   /* Get local element info */
885:   for (k=0; k < Mda*Nda*2; k++) {
886:     idx[0] = connect[k*3];
887:     idx[1] = connect[k*3+1];
888:     idx[2] = connect[k*3+2];

890:     PetscInt    row,cols[3],r;
891:     PetscScalar vals[3];

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

895:       /* cv_sum = (1.0e-3+cv_p[idx[0]] + cv_p[idx[1]] + cv_p[idx[2]])*user->Dv/(3.0*user->kBT); */
896:       /* ci_sum = (1.0e-3+ci_p[idx[0]] + ci_p[idx[1]] + ci_p[idx[2]])*user->Di/(3.0*user->kBT); */
897:       cv_sum = .0000069*user->Dv/(user->kBT);
898:       ci_sum = .0000069*user->Di/user->kBT;

900:       if (k%2 == 0) { /* odd triangle */

902:         row     = 5*idx[r];
903:         cols[0] = 5*idx[0];     vals[0] = dt*eM_2_odd[r][0]*cv_sum;
904:         cols[1] = 5*idx[1];     vals[1] = dt*eM_2_odd[r][1]*cv_sum;
905:         cols[2] = 5*idx[2];     vals[2] = dt*eM_2_odd[r][2]*cv_sum;

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

910:         row     = 5*idx[r]+2;
911:         cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_odd[r][0]*ci_sum;
912:         cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_odd[r][1]*ci_sum;
913:         cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_odd[r][2]*ci_sum;

915:         MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);

917:       } else {
918:         row     = 5*idx[r];
919:         cols[0] = 5*idx[0];     vals[0] = dt*eM_2_even[r][0]*cv_sum;
920:         cols[1] = 5*idx[1];     vals[1] = dt*eM_2_even[r][1]*cv_sum;
921:         cols[2] = 5*idx[2];     vals[2] = dt*eM_2_even[r][2]*cv_sum;

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

926:         row     = 5*idx[r]+2;
927:         cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_even[r][0]*ci_sum;
928:         cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_even[r][1]*ci_sum;
929:         cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_even[r][2]*ci_sum;
930:         /* Insert values in matrix M for 3nd dof */
931:         MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);
932:       }
933:     }
934:   }

936:   PetscFree(nodes);
937:   PetscFree(connect);
938:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
939:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
940:   VecRestoreArray(user->cv,&cv_p);
941:   VecRestoreArray(user->ci,&ci_p);
942:   return(0);
943: }