Actual source code: ex60.c

petsc-3.3-p7 2013-05-11
  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

 27:    
 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);
 75: 
 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,PETSC_NULL,PETSC_NULL,&user.da1);
 80: 
 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,PETSC_NULL,PETSC_NULL,&user.da2);


 85:   /* Set Element type (triangular) */
 86:   DMDASetElementType(user.da1,DMDA_ELEMENT_P1);
 87:   DMDASetElementType(user.da2,DMDA_ELEMENT_P1);
 88: 
 89:   /* Set x and y coordinates */
 90:   DMDASetUniformCoordinates(user.da1,user.xmin,user.xmax,user.ymin,user.ymax,PETSC_NULL,PETSC_NULL);
 91:   DMDASetUniformCoordinates(user.da2,user.xmin,user.xmax,user.ymin,user.ymax,PETSC_NULL,PETSC_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);
 99: 
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);
120: 

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);
130: 
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);
136: 
137:   SetVariableBounds(user.da1,xl,xu);
138:   SNESVISetVariableBounds(snes,xl,xu);
139:   SNESSetFromOptions(snes);
140: 
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);
146: 
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,PETSC_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:   }
179: 
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);
214: 
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;
226: 
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: 
251:   return(0);
252: }

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

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

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

274:   VecCopy(user->cv,user->work1);
275:   VecAXPY(user->work1,1.0,user->Pv);
276:   VecScale(user->work1,-1.0);
277:   MatMult(user->M_0,user->work1,user->work2);
278:   VecGetLocalSize(user->work1,&n);
279: 
280:  for (i=0;i<n;i++) {
281:        q_p[5*i]=w2[i];
282:   }
283: 
284:   MatMult(user->M_0,user->DPsiv,user->work1);
285:   for (i=0;i<n;i++) {
286:        q_p[5*i+1]=w1[i];
287:   }

289:   VecCopy(user->ci,user->work1);
290:   VecAXPY(user->work1,1.0,user->Pi);
291:   VecScale(user->work1,-1.0);
292:   MatMult(user->M_0,user->work1,user->work2);
293:   for (i=0;i<n;i++) {
294:        q_p[5*i+2]=w2[i];
295:   }

297:   MatMult(user->M_0,user->DPsii,user->work1);
298:   for (i=0;i<n;i++) {
299:        q_p[5*i+3]=w1[i];
300:   }

302:   VecCopy(user->eta,user->work1);
303:   VecScale(user->work1,-1.0/user->dt);
304:   VecAXPY(user->work1,user->L,user->DPsieta);
305:   VecAXPY(user->work1,-1.0,user->Piv);
306:   MatMult(user->M_0,user->work1,user->work2);
307:   for (i=0;i<n;i++) {
308:        q_p[5*i+4]=w2[i];
309:   }
310: 
311:   VecRestoreArray(user->q,&q_p);
312:   VecRestoreArray(user->work1,&w1);
313:   VecRestoreArray(user->work2,&w2);
314: 
315:   return(0);
316: }

320: PetscErrorCode DPsi(AppCtx* user)
321: {
322:   PetscErrorCode  ierr;
323:   PetscScalar     Evf=user->Evf,Eif=user->Eif,kBT=user->kBT,A=user->A;
324:   PetscScalar     *cv_p,*ci_p,*eta_p,*logcv_p,*logci_p,*logcvi_p,*DPsiv_p,*DPsii_p,*DPsieta_p;
325:   PetscInt        n,i;


329:   VecGetLocalSize(user->cv,&n);
330:   VecGetArray(user->cv,&cv_p);
331:   VecGetArray(user->ci,&ci_p);
332:   VecGetArray(user->eta,&eta_p);
333:   VecGetArray(user->logcv,&logcv_p);
334:   VecGetArray(user->logci,&logci_p);
335:   VecGetArray(user->logcvi,&logcvi_p);
336:   VecGetArray(user->DPsiv,&DPsiv_p);
337:   VecGetArray(user->DPsii,&DPsii_p);
338:   VecGetArray(user->DPsieta,&DPsieta_p);
339: 
340:   Llog(user->cv,user->logcv);
341: 
342:   Llog(user->ci,user->logci);
343: 

345:   VecCopy(user->cv,user->cvi);
346:   VecAXPY(user->cvi,1.0,user->ci);
347:   VecScale(user->cvi,-1.0);
348:   VecShift(user->cvi,1.0);
349:   Llog(user->cvi,user->logcvi);
350: 
351:   for (i=0;i<n;i++)
352:   {
353:     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);

355:     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] ;
356: 
357:     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]);
358: 
359: 
360:   }

362:   VecRestoreArray(user->cv,&cv_p);
363:   VecRestoreArray(user->ci,&ci_p);
364:   VecRestoreArray(user->eta,&eta_p);
365:   VecRestoreArray(user->logcv,&logcv_p);
366:   VecRestoreArray(user->logci,&logci_p);
367:   VecRestoreArray(user->logcvi,&logcvi_p);
368:   VecRestoreArray(user->DPsiv,&DPsiv_p);
369:   VecRestoreArray(user->DPsii,&DPsii_p);
370:   VecRestoreArray(user->DPsieta,&DPsieta_p);


373:   return(0);

375: }


380: PetscErrorCode Llog(Vec X, Vec Y)
381: {
382:   PetscErrorCode    ierr;
383:   PetscScalar       *x,*y;
384:   PetscInt          n,i;

387: 
388:   VecGetArray(X,&x);
389:   VecGetArray(Y,&y);
390:   VecGetLocalSize(X,&n);
391:   for (i=0;i<n;i++) {
392:     if (x[i] < 1.0e-12) {
393:       y[i] = log(1.0e-12);
394:     }
395:     else {
396:       y[i] = log(x[i]);
397:     }
398:   }
399: 
400:   return(0);
401: }


406: PetscErrorCode SetInitialGuess(Vec X,AppCtx* user)
407: {
408:   PetscErrorCode    ierr;
409:   PetscInt          n,i;
410:   PetscScalar           *xx,*cv_p,*ci_p,*wv_p,*wi_p;
411:   PetscViewer       view;
412:   PetscScalar       initv = .00069;


416:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_initial",FILE_MODE_WRITE,&view);
417:   VecGetLocalSize(X,&n);

419:   VecSet(user->cv,initv);
420:   VecSet(user->ci,initv);
421:   VecSet(user->eta,0.0);

423:   DPsi(user);
424:   VecCopy(user->DPsiv,user->wv);
425:   VecCopy(user->DPsii,user->wi);

427:   VecGetArray(X,&xx);
428:   VecGetArray(user->cv,&cv_p);
429:   VecGetArray(user->ci,&ci_p);
430:   VecGetArray(user->wv,&wv_p);
431:   VecGetArray(user->wi,&wi_p);
432:   for (i=0;i<n/5;i++)
433:   {
434:     xx[5*i]=wv_p[i];
435:     xx[5*i+1]=cv_p[i];
436:     xx[5*i+2]=wi_p[i];
437:     xx[5*i+3]=ci_p[i];
438:     xx[5*i+4]=0.0;
439:   }

441:   VecView(user->wv,view);
442:   VecView(user->cv,view);
443:   VecView(user->wi,view);
444:   VecView(user->ci,view);
445:   PetscViewerDestroy(&view);

447:   VecRestoreArray(X,&xx);
448:   VecRestoreArray(user->cv,&cv_p);
449:   VecRestoreArray(user->ci,&ci_p);
450:   VecRestoreArray(user->wv,&wv_p);
451:   VecRestoreArray(user->wi,&wi_p);
452: 
453:   return(0);
454: }

458: PetscErrorCode SetRandomVectors(AppCtx* user)
459: {
461:   PetscInt       i,n,count=0;
462:   PetscScalar    *w1,*w2,*Pv_p,*eta_p;
463: 

466: 
467:   VecSetRandom(user->work1,PETSC_NULL);
468:   VecSetRandom(user->work2,PETSC_NULL);
469:   VecGetArray(user->work1,&w1);
470:   VecGetArray(user->work2,&w2);
471:   VecGetArray(user->Pv,&Pv_p);
472:   VecGetArray(user->eta,&eta_p);
473:   VecGetLocalSize(user->work1,&n);
474:   for (i=0;i<n;i++) {
475: 
476:     if (eta_p[i]>=0.8 || w1[i]>user->P_casc){
477:       Pv_p[i]=0;
478: 
479:     }
480:     else
481:     {
482:       Pv_p[i]=w2[i]*user->VG;
483:       count=count+1;
484:     }

486:   }
487: 
488:   VecCopy(user->Pv,user->Pi);
489:   VecScale(user->Pi,0.9);
490:   VecPointwiseMult(user->Piv,user->Pi,user->Pv);
491:   VecRestoreArray(user->work1,&w1);
492:   VecRestoreArray(user->work2,&w2);
493:   VecRestoreArray(user->Pv,&Pv_p);
494:   VecRestoreArray(user->eta,&eta_p);

496:   return(0);
497: 
498: }

502: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void* ctx)
503: {
505:   AppCtx         *user=(AppCtx*)ctx;
506: 
508:   MatMultAdd(user->M,X,user->q,F);
509:   return(0);
510: }

514: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flg,void *ctx)
515: {
517:   AppCtx         *user=(AppCtx*)ctx;
518: 
520:   *flg = SAME_NONZERO_PATTERN;
521:   MatCopy(user->M,*J,*flg);
522:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
523:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
524:   return(0);
525: }
528: PetscErrorCode SetVariableBounds(DM da,Vec xl,Vec xu)
529: {
531:   PetscScalar    ***l,***u;
532:   PetscInt       xs,xm,ys,ym;
533:   PetscInt       j,i;
534: 
536:   DMDAVecGetArrayDOF(da,xl,&l);
537:   DMDAVecGetArrayDOF(da,xu,&u);
538: 
539:   DMDAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL);
540: 
541:   for (j=ys; j<ys+ym; j++) {
542:     for(i=xs; i < xs+xm;i++) {
543:       l[j][i][0] = -SNES_VI_INF;
544:       l[j][i][1] = 0.0;
545:       l[j][i][2] = -SNES_VI_INF;
546:       l[j][i][3] = 0.0;
547:       l[j][i][4] = 0.0;
548:       u[j][i][0] = SNES_VI_INF;
549:       u[j][i][1] = 1.0;
550:       u[j][i][2] = SNES_VI_INF;
551:       u[j][i][3] = 1.0;
552:       u[j][i][4] = 1.0;
553:     }
554:   }
555: 
556: 
557:   DMDAVecRestoreArrayDOF(da,xl,&l);
558:   DMDAVecRestoreArrayDOF(da,xu,&u);
559:   return(0);
560: }

564: PetscErrorCode GetParams(AppCtx* user)
565: {
567:   PetscBool      flg;
568: 
570: 
571:   /* Set default parameters */
572:   user->xmin = 0.0; user->xmax = 1.0;
573:   user->ymin = 0.0; user->ymax = 1.0;
574:   user->Dv = 1.0; user->Di=4.0;
575:   user->Evf = 0.8; user->Eif = 1.2;
576:   user->A = 1.0;
577:   user->kBT = 0.11;
578:   user->kav = 1.0; user->kai = 1.0; user->kaeta = 1.0;
579:   user->Rsurf = 10.0; user->Rbulk = 1.0;
580:   user->L = 10.0; user->P_casc = 0.05;
581:   user->T = 1.0e-2;    user->dt = 1.0e-4;
582:   user->VG = 100.0;

584:   PetscOptionsGetReal(PETSC_NULL,"-xmin",&user->xmin,&flg);
585:   PetscOptionsGetReal(PETSC_NULL,"-xmax",&user->xmax,&flg);
586:   PetscOptionsGetReal(PETSC_NULL,"-T",&user->T,&flg);
587:   PetscOptionsGetReal(PETSC_NULL,"-dt",&user->dt,&flg);
588:   PetscOptionsGetReal(PETSC_NULL,"-P_casc",&user->P_casc,&flg);
589:   PetscOptionsGetReal(PETSC_NULL,"-VG",&user->VG,&flg);
590: 

592:   return(0);
593:  }


598: PetscErrorCode SetUpMatrices(AppCtx* user)
599: {
600:   PetscErrorCode    ierr;
601:   PetscInt          nele,nen,i,j,n;
602:   const PetscInt    *ele;
603:   PetscScalar       dt=user->dt,hx,hy;
604: 
605:   PetscInt          idx[3],*nodes, *connect, k;
606:   PetscScalar       eM_0[3][3],eM_2_even[3][3],eM_2_odd[3][3];
607:   PetscScalar       cv_sum, ci_sum;
608:   Mat               M=user->M;
609:   Mat               M_0=user->M_0;
610:   PetscInt          Mda=user->Mda, Nda=user->Nda, ld, rd, ru, lu;
611:   PetscScalar       *cv_p,*ci_p;
612: 
614: 
615:   /*  MatSetOption(M,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
616:    MatSetOption(M_0,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);*/

618:   /* Create the mass matrix M_0 */
619:   VecGetArray(user->cv,&cv_p);
620:   VecGetArray(user->ci,&ci_p);
621:   MatGetLocalSize(M,&n,PETSC_NULL);
622:   DMDAGetInfo(user->da1,PETSC_NULL,&Mda,&Nda,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);

624:   PetscMalloc((Mda+1)*(Nda+1)*sizeof(PetscInt),&nodes);
625:   PetscMalloc(Mda*Nda*2*3*sizeof(PetscInt),&connect);
626:   hx = (user->xmax-user->xmin)/Mda;
627:   hy = (user->ymax-user->ymin)/Nda;
628:   for (j=0;j < Nda;j++) {
629:     for (i=0;i < Mda;i++) {
630:       nodes[j*(Mda+1)+i] = j*Mda+i;
631:     }
632:     nodes[j*(Mda+1)+Mda] = j*Mda;
633:   }
634: 
635:   for (i=0;i < Mda;i++){
636:     nodes[Nda*(Mda+1)+i] = i;
637:   }

639:   nodes[Nda*(Mda+1)+Mda] = 0;
640: 
641: 
642:   k = 0;
643:   for (j=0;j<Nda;j++) {
644:     for (i=0;i<Mda;i++) {
645: 
646:       /* ld = nodes[j][i]; */
647:       ld = nodes[j*(Mda+1)+i];
648:       /* rd = nodes[j+1][i]; */
649:       rd = nodes[(j+1)*(Mda+1)+i];
650:       /* ru = nodes[j+1][i+1]; */
651:       ru = nodes[(j+1)*(Mda+1)+i+1];
652:       /* lu = nodes[j][i+1]; */
653:       lu = nodes[j*(Mda+1)+i+1];

655:       /* connect[k][0]=ld; */
656:       connect[k*6]=ld;
657:       /* connect[k][1]=lu; */
658:       connect[k*6+1]=lu;
659:       /* connect[k][2]=ru; */
660:       connect[k*6+2]=rd;
661:       connect[k*6+3]=lu;
662:       connect[k*6+4]=ru;
663:       connect[k*6+5]=rd;
664: 
665:       k = k+1;
666:     }
667:   }
668: 
669: 
670:   eM_0[0][0]=eM_0[1][1]=eM_0[2][2]=hx*hy/12.0;
671:   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;
672: 
673:   eM_2_odd[0][0] = 1.0;
674:   eM_2_odd[1][1] = eM_2_odd[2][2] = 0.5;
675:   eM_2_odd[0][1] = eM_2_odd[0][2] = eM_2_odd[1][0]= eM_2_odd[2][0] = -0.5;
676:   eM_2_odd[1][2] = eM_2_odd[2][1] = 0.0;

678:   eM_2_even[1][1] = 1.0;
679:   eM_2_even[0][0] = eM_2_even[2][2] = 0.5;
680:   eM_2_even[0][1] = eM_2_even[1][0] = eM_2_even[1][2] = eM_2_even[2][1] = -0.5;
681:   eM_2_even[0][2] = eM_2_even[2][0] = 0.0;


684:   for(k=0;k < Mda*Nda*2;k++) {
685:     idx[0] = connect[k*3];
686:     idx[1] = connect[k*3+1];
687:     idx[2] = connect[k*3+2];
688: 
689:     PetscInt    row,cols[6],r,row_M_0,cols3[3];
690:     PetscScalar vals[6],vals_M_0[3],vals3[3];
691: 
692:     for(r=0;r<3;r++) {
693:       row_M_0 = connect[k*3+r];
694: 
695:       vals_M_0[0]=eM_0[r][0];
696:       vals_M_0[1]=eM_0[r][1];
697:       vals_M_0[2]=eM_0[r][2];
698: 
699: 
700:       MatSetValues(M_0,1,&row_M_0,3,idx,vals_M_0,ADD_VALUES);
701: 
702:       //cv_sum = (cv_p[idx[0]] + cv_p[idx[1]] + cv_p[idx[2]])*user->Dv/(3.0*user->kBT);
703:       //ci_sum = (ci_p[idx[0]] + ci_p[idx[1]] + ci_p[idx[2]])*user->Di/(3.0*user->kBT);
704:       cv_sum = .0000069*user->Dv/user->kBT;
705:       ci_sum = .0000069*user->Di/user->kBT;
706: 
707:       if (k%2 == 0)  {
708: 
709:         row = 5*idx[r];
710:         cols[0] = 5*idx[0];     vals[0] = dt*eM_2_odd[r][0]*cv_sum;
711:         cols[1] = 5*idx[1];     vals[1] = dt*eM_2_odd[r][1]*cv_sum;
712:         cols[2] = 5*idx[2];     vals[2] = dt*eM_2_odd[r][2]*cv_sum;
713:         cols[3] = 5*idx[0]+1;   vals[3] = eM_0[r][0];
714:         cols[4] = 5*idx[1]+1;   vals[4] = eM_0[r][1];
715:         cols[5] = 5*idx[2]+1;   vals[5] = eM_0[r][2];

717: 
718:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
719: 
720: 
721:         row = 5*idx[r]+1;
722:         cols[0] = 5*idx[0];     vals[0] = -1.0*eM_0[r][0];
723:         cols[1] = 5*idx[1];     vals[1] = -1.0*eM_0[r][1];
724:         cols[2] = 5*idx[2];     vals[2] = -1.0*eM_0[r][2];
725:         cols[3] = 5*idx[0]+1;   vals[3] =  user->kav*eM_2_odd[r][0];
726:         cols[4] = 5*idx[1]+1;   vals[4] =  user->kav*eM_2_odd[r][1];
727:         cols[5] = 5*idx[2]+1;   vals[5] =  user->kav*eM_2_odd[r][2];
728: 
729:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
730: 
731:         row = 5*idx[r]+2;
732:         cols[0] = 5*idx[0]+2;   vals[0] =  dt*eM_2_odd[r][0]*ci_sum;
733:         cols[1] = 5*idx[1]+2;   vals[1] =  dt*eM_2_odd[r][1]*ci_sum;
734:         cols[2] = 5*idx[2]+2;   vals[2] =  dt*eM_2_odd[r][2]*ci_sum;
735:         cols[3] = 5*idx[0]+3;   vals[3] =  eM_0[r][0];
736:         cols[4] = 5*idx[1]+3;   vals[4] =  eM_0[r][1];
737:         cols[5] = 5*idx[2]+3;   vals[5] =  eM_0[r][2];
738: 
739:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
740: 
741: 
742:         row = 5*idx[r]+3;
743:         cols[0] = 5*idx[0]+2;   vals[0] = -1.0*eM_0[r][0];
744:         cols[1] = 5*idx[1]+2;   vals[1] = -1.0*eM_0[r][1];
745:         cols[2] = 5*idx[2]+2;   vals[2] = -1.0*eM_0[r][2];
746:         cols[3] = 5*idx[0]+3;   vals[3] =  user->kai*eM_2_odd[r][0];
747:         cols[4] = 5*idx[1]+3;   vals[4] =  user->kai*eM_2_odd[r][1];
748:         cols[5] = 5*idx[2]+3;   vals[5] =  user->kai*eM_2_odd[r][2];
749: 
750:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
751: 
752: 
753:         row = 5*idx[r]+4;
754:         /*
755:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*dt*eM_2_odd[r][0];
756:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*dt*eM_2_odd[r][1];
757:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*dt*eM_2_odd[r][2];
758:          */
759:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*eM_2_odd[r][0];
760:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*eM_2_odd[r][1];
761:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*eM_2_odd[r][2];
762: 
763:         MatSetValuesLocal(M,1,&row,3,cols3,vals3,ADD_VALUES);
764: 
765: 
766:       }
767: 
768:       else {
769: 
770: 
771:         row = 5*idx[r];
772:         cols[0] = 5*idx[0];     vals[0] = dt*eM_2_even[r][0]*cv_sum;
773:         cols[1] = 5*idx[1];     vals[1] = dt*eM_2_even[r][1]*cv_sum;
774:         cols[2] = 5*idx[2];     vals[2] = dt*eM_2_even[r][2]*cv_sum;
775:         cols[3] = 5*idx[0]+1;   vals[3] = eM_0[r][0];
776:         cols[4] = 5*idx[1]+1;   vals[4] = eM_0[r][1];
777:         cols[5] = 5*idx[2]+1;   vals[5] = eM_0[r][2];


780: 
781:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
782: 
783: 
784:         row = 5*idx[r]+1;
785:         cols[0] = 5*idx[0];     vals[0] = -1.0*eM_0[r][0];
786:         cols[1] = 5*idx[1];     vals[1] = -1.0*eM_0[r][1];
787:         cols[2] = 5*idx[2];     vals[2] = -1.0*eM_0[r][2];
788:         cols[3] = 5*idx[0]+1;   vals[3] =  user->kav*eM_2_even[r][0];
789:         cols[4] = 5*idx[1]+1;   vals[4] =  user->kav*eM_2_even[r][1];
790:         cols[5] = 5*idx[2]+1;   vals[5] =  user->kav*eM_2_even[r][2];
791: 
792:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
793: 
794:         row = 5*idx[r]+2;
795:         cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_even[r][0]*ci_sum;
796:         cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_even[r][1]*ci_sum;
797:         cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_even[r][2]*ci_sum;
798:         cols[3] = 5*idx[0]+3;   vals[3] = eM_0[r][0];
799:         cols[4] = 5*idx[1]+3;   vals[4] = eM_0[r][1];
800:         cols[5] = 5*idx[2]+3;   vals[5] = eM_0[r][2];
801: 
802:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
803: 
804:         row = 5*idx[r]+3;
805:         cols[0] = 5*idx[0]+2;   vals[0] = -1.0*eM_0[r][0];
806:         cols[1] = 5*idx[1]+2;   vals[1] = -1.0*eM_0[r][1];
807:         cols[2] = 5*idx[2]+2;   vals[2] = -1.0*eM_0[r][2];
808:         cols[3] = 5*idx[0]+3;   vals[3] =  user->kai*eM_2_even[r][0];
809:         cols[4] = 5*idx[1]+3;   vals[4] =  user->kai*eM_2_even[r][1];
810:         cols[5] = 5*idx[2]+3;   vals[5] =  user->kai*eM_2_even[r][2];
811: 
812:         MatSetValuesLocal(M,1,&row,6,cols,vals,ADD_VALUES);
813: 
814:         row = 5*idx[r]+4;
815:         /*
816:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*dt*eM_2_even[r][0];
817:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*dt*eM_2_even[r][1];
818:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*dt*eM_2_even[r][2];
819:          */
820:         cols3[0] = 5*idx[0]+4;   vals3[0] = eM_0[r][0]/dt + user->L*user->kaeta*eM_2_even[r][0];
821:         cols3[1] = 5*idx[1]+4;   vals3[1] = eM_0[r][1]/dt + user->L*user->kaeta*eM_2_even[r][1];
822:         cols3[2] = 5*idx[2]+4;   vals3[2] = eM_0[r][2]/dt + user->L*user->kaeta*eM_2_even[r][2];
823:         MatSetValuesLocal(M,1,&row,3,cols3,vals3,ADD_VALUES);
824: 
825:       }
826: 
827:     }
828:   }

830:   PetscFree(nodes);
831:   PetscFree(connect);

833:   VecRestoreArray(user->cv,&cv_p);
834:   VecRestoreArray(user->ci,&ci_p);
835:   MatAssemblyBegin(M_0,MAT_FINAL_ASSEMBLY);
836:   MatAssemblyEnd(M_0,MAT_FINAL_ASSEMBLY);
837: 
838:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
839:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
840: 
841:   DMDARestoreElements(user->da1,&nele,&nen,&ele);
842: 

844:   return(0);
845: }


850: PetscErrorCode UpdateMatrices(AppCtx* user)
851: {
852:   PetscErrorCode    ierr;
853:   PetscInt          i,j,n,Mda,Nda;
854: 
855:   PetscInt          idx[3],*nodes,*connect,k;
856:   PetscInt          ld,rd,lu,ru;
857:   PetscScalar       eM_2_odd[3][3],eM_2_even[3][3],h,dt=user->dt;
858:   Mat               M=user->M;
859:   PetscScalar       *cv_p,*ci_p,cv_sum,ci_sum;
861: 
862:   /* Create the mass matrix M_0 */
863:   MatGetLocalSize(M,&n,PETSC_NULL);
864:   VecGetArray(user->cv,&cv_p);
865:   VecGetArray(user->ci,&ci_p);
866:   DMDAGetInfo(user->da1,PETSC_NULL,&Mda,&Nda,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);

868:   PetscMalloc((Mda+1)*(Nda+1)*sizeof(PetscInt),&nodes);
869:   PetscMalloc(Mda*Nda*2*3*sizeof(PetscInt),&connect);
870: 
871:   h = 100.0/Mda;

873:   for (j=0;j < Nda;j++) {
874:     for (i=0;i < Mda;i++) {
875:       nodes[j*(Mda+1)+i] = j*Mda+i;
876:     }
877:     nodes[j*(Mda+1)+Mda] = j*Mda;
878:   }
879:   for (i=0;i < Mda;i++){
880:     nodes[Nda*(Mda+1)+i]=i;
881:   }
882:   nodes[Nda*(Mda+1)+Mda]=0;

884: 
885:   k = 0;
886:   for (j=0;j<Nda;j++) {
887:     for (i=0;i<Mda;i++) {
888:       ld = nodes[j*(Mda+1)+i];
889:       rd = nodes[(j+1)*(Mda+1)+i];
890:       ru = nodes[(j+1)*(Mda+1)+i+1];
891:       lu = nodes[j*(Mda+1)+i+1];
892:       connect[k*6]=ld;
893:       connect[k*6+1]=lu;
894:       connect[k*6+2]=rd;
895:       connect[k*6+3]=lu;
896:       connect[k*6+4]=ru;
897:       connect[k*6+5]=rd;
898:       k = k+1;
899:     }
900:   }
901: 
902:   for(k=0;k < Mda*Nda*2;k++) {
903:     idx[0] = connect[k*3];
904:     idx[1] = connect[k*3+1];
905:     idx[2] = connect[k*3+2];
906: 
907:     PetscInt r,row,cols[3];
908:     PetscScalar vals[3];
909:     for (r=0;r<3;r++) {
910:       row = 5*idx[r];
911:       cols[0] = 5*idx[0];     vals[0] = 0.0;
912:       cols[1] = 5*idx[1];     vals[1] = 0.0;
913:       cols[2] = 5*idx[2];     vals[2] = 0.0;

915:       /* Insert values in matrix M for 1st dof */
916:       MatSetValuesLocal(M,1,&row,3,cols,vals,INSERT_VALUES);
917: 
918:       row = 5*idx[r]+2;
919:       cols[0] = 5*idx[0]+2;   vals[0] = 0.0;
920:       cols[1] = 5*idx[1]+2;   vals[1] = 0.0;
921:       cols[2] = 5*idx[2]+2;   vals[2] = 0.0;

923:       /* Insert values in matrix M for 3nd dof */
924:       MatSetValuesLocal(M,1,&row,3,cols,vals,INSERT_VALUES);
925:     }
926:   }
927: 
928:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
929:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

931:   eM_2_odd[0][0] = 1.0;
932:   eM_2_odd[1][1] = eM_2_odd[2][2] = 0.5;
933:   eM_2_odd[0][1] = eM_2_odd[0][2] = eM_2_odd[1][0]= eM_2_odd[2][0] = -0.5;
934:   eM_2_odd[1][2] = eM_2_odd[2][1] = 0.0;

936:   eM_2_even[1][1] = 1.0;
937:   eM_2_even[0][0] = eM_2_even[2][2] = 0.5;
938:   eM_2_even[0][1] = eM_2_even[1][0] = eM_2_even[1][2] = eM_2_even[2][1] = -0.5;
939:   eM_2_even[0][2] = eM_2_even[2][0] = 0.0;

941: 
942:   /* Get local element info */
943:   for(k=0;k < Mda*Nda*2;k++) {
944:     idx[0] = connect[k*3];
945:       idx[1] = connect[k*3+1];
946:       idx[2] = connect[k*3+2];
947: 
948:       PetscInt    row,cols[3],r;
949:       PetscScalar vals[3];
950: 
951:       for(r=0;r<3;r++) {
952: 
953:         // cv_sum = (1.0e-3+cv_p[idx[0]] + cv_p[idx[1]] + cv_p[idx[2]])*user->Dv/(3.0*user->kBT);
954:         //ci_sum = (1.0e-3+ci_p[idx[0]] + ci_p[idx[1]] + ci_p[idx[2]])*user->Di/(3.0*user->kBT);
955:         cv_sum = .0000069*user->Dv/(user->kBT);
956:         ci_sum = .0000069*user->Di/user->kBT;

958:             if (k%2 == 0) /* odd triangle */ {
959: 
960:                 row = 5*idx[r];
961:                 cols[0] = 5*idx[0];     vals[0] = dt*eM_2_odd[r][0]*cv_sum;
962:                 cols[1] = 5*idx[1];     vals[1] = dt*eM_2_odd[r][1]*cv_sum;
963:                 cols[2] = 5*idx[2];     vals[2] = dt*eM_2_odd[r][2]*cv_sum;
964: 
965:                 /* Insert values in matrix M for 1st dof */
966:                 MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);
967: 
968: 
969:                 row = 5*idx[r]+2;
970:                 cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_odd[r][0]*ci_sum;
971:                 cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_odd[r][1]*ci_sum;
972:                 cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_odd[r][2]*ci_sum;
973: 
974:                 MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);
975: 
976:             }
977: 
978:             else {
979:                 row = 5*idx[r];
980:                 cols[0] = 5*idx[0];     vals[0] = dt*eM_2_even[r][0]*cv_sum;
981:                 cols[1] = 5*idx[1];     vals[1] = dt*eM_2_even[r][1]*cv_sum;
982:                 cols[2] = 5*idx[2];     vals[2] = dt*eM_2_even[r][2]*cv_sum;
983: 
984:                 /* Insert values in matrix M for 1st dof */
985:                 MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);
986: 
987: 
988:                 row = 5*idx[r]+2;
989:                 cols[0] = 5*idx[0]+2;   vals[0] = dt*eM_2_even[r][0]*ci_sum;
990:                 cols[1] = 5*idx[1]+2;   vals[1] = dt*eM_2_even[r][1]*ci_sum;
991:                 cols[2] = 5*idx[2]+2;   vals[2] = dt*eM_2_even[r][2]*ci_sum;
992:                 /* Insert values in matrix M for 3nd dof */
993:                 MatSetValuesLocal(M,1,&row,3,cols,vals,ADD_VALUES);
994: 
995:             }
996:         }
997: 
998:     }

1000:   PetscFree(nodes);
1001:   PetscFree(connect);
1002:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
1003:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
1004:   VecRestoreArray(user->cv,&cv_p);
1005:   VecRestoreArray(user->ci,&ci_p);


1008:   return(0);
1009: }