Actual source code: ex64.c

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

 11: /*
 12: ./ex63 -ksp_type fgmres -snes_vi_monitor   -snes_atol 1.e-11 -snes_converged_reason -ksp_converged_reason   -snes_linesearch_monitor -VG 10000000 -draw_fields 1,3,4  -pc_type mg -pc_mg_galerkin -log_summary -dt .000001 -mg_coarse_pc_type svd  -ksp_monitor_true_residual -ksp_rtol 1.e-9 -snes_linesearch_type basic -T .0020 -P_casc .0005 -snes_monitor_solution -da_refine 10
 13:  */

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

 18: typedef struct{
 19:   PetscReal   dt,T; /* Time step and end time */
 20:   DM          da1,da1_clone,da2;
 21:   Mat         M;    /* Jacobian matrix */
 22:   Mat         M_0;
 23:   Vec         q,wv,cv,eta,DPsiv,DPsieta,logcv,logcv2;
 24:   Vec         work1,work2;
 25:   PetscScalar Mv,L,kaeta,kav,Evf,A,B,cv0,Sv;
 26:   PetscReal   xmin,xmax;
 27:   PetscInt    nx;
 28:   PetscBool   graphics;
 29:   PetscBool   periodic;
 30:   PetscBool   lumpedmass;
 31: }AppCtx;

 33: PetscErrorCode GetParams(AppCtx*);
 34: PetscErrorCode SetVariableBounds(DM,Vec,Vec);
 35: PetscErrorCode SetUpMatrices(AppCtx*);
 36: PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
 37: PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 38: PetscErrorCode SetInitialGuess(Vec,AppCtx*);
 39: PetscErrorCode Update_q(AppCtx*);
 40: PetscErrorCode Update_u(Vec,AppCtx*);
 41: PetscErrorCode DPsi(AppCtx*);
 42: PetscErrorCode Llog(Vec,Vec);
 43: PetscErrorCode CheckRedundancy(SNES,IS,IS*,DM);

 47: int main(int argc, char **argv)
 48: {
 50:   Vec            x,r;  /* solution and residual vectors */
 51:   SNES           snes; /* Nonlinear solver context */
 52:   AppCtx         user; /* Application context */
 53:   Vec            xl,xu; /* Upper and lower bounds on variables */
 54:   Mat            J;
 55:   PetscScalar    t=0.0;
 56:   PetscViewer    view_out, view_p, view_q, view_psi, view_mat;
 57:   PetscReal      bounds[] = {1000.0,-1000.,0.0,1.0,1000.0,-1000.0,0.0,1.0,1000.0,-1000.0};


 60:   PetscInitialize(&argc,&argv, (char*)0, help);
 61: 
 62:   /* Get physics and time parameters */
 63:   GetParams(&user);

 65:   if (user.periodic) {
 66:     DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_PERIODIC, -4, 3, 1,PETSC_NULL,&user.da1);
 67:     DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_PERIODIC, -4, 3, 1,PETSC_NULL,&user.da1_clone);
 68:     DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_PERIODIC, -4, 1, 1,PETSC_NULL,&user.da2);
 69:   } else {
 70:     DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE, -4, 3, 1,PETSC_NULL,&user.da1);
 71:     DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE, -4, 3, 1,PETSC_NULL,&user.da1_clone);
 72:     DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE, -4, 1, 1,PETSC_NULL,&user.da2);

 74:   }
 75:   /* Set Element type (triangular) */
 76:   DMDASetElementType(user.da1,DMDA_ELEMENT_P1);
 77:   DMDASetElementType(user.da2,DMDA_ELEMENT_P1);
 78: 
 79:   /* Set x and y coordinates */
 80:   DMDASetUniformCoordinates(user.da1,user.xmin,user.xmax,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 81:   DMDASetUniformCoordinates(user.da2,user.xmin,user.xmax,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 82:   /* Get global vector x from DM (da1) and duplicate vectors r,xl,xu */
 83:   DMCreateGlobalVector(user.da1,&x);
 84:   VecDuplicate(x,&r);
 85:   VecDuplicate(x,&xl);
 86:   VecDuplicate(x,&xu);
 87:   VecDuplicate(x,&user.q);
 88: 
 89:   /* Get global vector user->wv from da2 and duplicate other vectors */
 90:   DMCreateGlobalVector(user.da2,&user.wv);
 91:   VecDuplicate(user.wv,&user.cv);
 92:   VecDuplicate(user.wv,&user.eta);
 93:   VecDuplicate(user.wv,&user.DPsiv);
 94:   VecDuplicate(user.wv,&user.DPsieta);
 95:   VecDuplicate(user.wv,&user.logcv);
 96:   VecDuplicate(user.wv,&user.logcv2);
 97:   VecDuplicate(user.wv,&user.work1);
 98:   VecDuplicate(user.wv,&user.work2);

100:   /* Get Jacobian matrix structure from the da for the entire thing, da1 */
101:   DMCreateMatrix(user.da1,MATAIJ,&user.M);
102:   /* Get the (usual) mass matrix structure from da2 */
103:   DMCreateMatrix(user.da2,MATAIJ,&user.M_0);
104:   /* Form the jacobian matrix and M_0 */
105:   SetUpMatrices(&user);
106:   SetInitialGuess(x,&user);
107:   MatDuplicate(user.M,MAT_DO_NOT_COPY_VALUES,&J);
108: 
109:   SNESCreate(PETSC_COMM_WORLD,&snes);
110:   SNESSetDM(snes,user.da1);

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

115:   SetVariableBounds(user.da1,xl,xu);
116:   SNESVISetVariableBounds(snes,xl,xu);
117:   SNESSetFromOptions(snes);
118:   //SNESVISetRedundancyCheck(snes,(PetscErrorCode (*)(SNES,IS,IS*,void*))CheckRedundancy,user.da1_clone);
119:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_out",FILE_MODE_WRITE,&view_out);
120:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_p",FILE_MODE_WRITE,&view_p);
121:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_q",FILE_MODE_WRITE,&view_q);
122:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_psi",FILE_MODE_WRITE,&view_psi);
123:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_mat",FILE_MODE_WRITE,&view_mat);
124:   /* PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),5,bounds); */

126:   if (user.graphics) {
127:     PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),5,bounds);

129:     VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
130:   }
131:   while(t<user.T) {

133:     char         filename[PETSC_MAX_PATH_LEN];
134:     PetscScalar  a = 1.0;
135:     PetscInt     i;
136:     PetscViewer  view;


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

142:     DPsi(&user);
143:     VecView(user.DPsiv,view_psi);
144:     VecView(user.DPsieta,view_psi);

146:     Update_q(&user);
147:     VecView(user.q,view_q);
148:     MatView(user.M,view_mat);

150:     SNESSolve(snes,PETSC_NULL,x);
151:     VecView(x,view_out);


154:     if (user.graphics) {
155:       VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
156:     }
157: 
158:     PetscInt its;
159:     SNESGetIterationNumber(snes,&its);
160:     PetscPrintf(PETSC_COMM_WORLD,"SNESVI solver converged at t = %5.4g in %d iterations\n",t,its);

162:     Update_u(x,&user);

164:     for (i=0; i < (int)(user.T/a) ; i++) {
165:       if (t/a > i - user.dt/a && t/a < i + user.dt/a) {
166:         sprintf(filename,"output_%f",t);
167:         PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&view);
168:         VecView(user.cv,view);
169:         VecView(user.eta,view);
170:         PetscViewerDestroy(&view);
171:       }
172: 
173:     }

175:     t = t + user.dt;
176: 
177:   }

179: 
180:   PetscViewerDestroy(&view_out);
181:   PetscViewerDestroy(&view_p);
182:   PetscViewerDestroy(&view_q);
183:   PetscViewerDestroy(&view_psi);
184:   PetscViewerDestroy(&view_mat);
185:   VecDestroy(&x);
186:   VecDestroy(&r);
187:   VecDestroy(&xl);
188:   VecDestroy(&xu);
189:   VecDestroy(&user.q);
190:   VecDestroy(&user.wv);
191:   VecDestroy(&user.cv);
192:   VecDestroy(&user.eta);
193:   VecDestroy(&user.DPsiv);
194:   VecDestroy(&user.DPsieta);
195:   VecDestroy(&user.logcv);
196:   VecDestroy(&user.logcv2);
197:   VecDestroy(&user.work1);
198:   VecDestroy(&user.work2);
199:   MatDestroy(&user.M);
200:   MatDestroy(&user.M_0);
201:   DMDestroy(&user.da1);
202:   DMDestroy(&user.da1_clone);
203:   DMDestroy(&user.da2);
204:   SNESDestroy(&snes);
205:   PetscFinalize();
206:   return 0;
207: }

211: PetscErrorCode Update_u(Vec X,AppCtx *user)
212: {
214:   PetscInt       i,n;
215:   PetscScalar    *xx,*wv_p,*cv_p,*eta_p;
216: 
218:   VecGetLocalSize(user->wv,&n);
219: 
220:   VecGetArray(X,&xx);
221:   VecGetArray(user->wv,&wv_p);
222:   VecGetArray(user->cv,&cv_p);
223:   VecGetArray(user->eta,&eta_p);
224: 
225: 
226:   for(i=0;i<n;i++) {
227:     wv_p[i] = xx[3*i];
228:     cv_p[i] = xx[3*i+1];
229:     eta_p[i] = xx[3*i+2];
230:   }
231:   VecRestoreArray(X,&xx);
232:   VecRestoreArray(user->wv,&wv_p);
233:   VecRestoreArray(user->cv,&cv_p);
234:   VecRestoreArray(user->eta,&eta_p);
235: 
236:   return(0);
237: }

241: PetscErrorCode Update_q(AppCtx *user)
242: {
244:   PetscScalar    *q_p, *w1, *w2;
245:   PetscInt       n,i;

248: 
249:   VecGetArray(user->q,&q_p);
250:   VecGetArray(user->work1,&w1);
251:   VecGetArray(user->work2,&w2);
252:   VecGetLocalSize(user->wv,&n);

254:   MatMult(user->M_0,user->cv,user->work1);
255:   VecScale(user->work1,-1.0);
256:   for (i=0;i<n;i++) {
257:     q_p[3*i]=w1[i];
258:   }

260:   MatMult(user->M_0,user->DPsiv,user->work1);
261:   for (i=0;i<n;i++) {
262:     q_p[3*i+1]=w1[i];
263:   }

265:   VecCopy(user->DPsieta,user->work1);
266:   VecScale(user->work1,user->L*user->dt);
267:   VecAXPY(user->work1,-1.0,user->eta);
268:   MatMult(user->M_0,user->work1,user->work2);
269:   for (i=0;i<n;i++) {
270:     q_p[3*i+2]=w2[i];
271:   }

273:   VecRestoreArray(user->q,&q_p);
274:   VecRestoreArray(user->work1,&w1);
275:   VecRestoreArray(user->work2,&w2);
276: 
277:   return(0);
278: }

282: PetscErrorCode DPsi(AppCtx* user)
283: {
284:   PetscErrorCode  ierr;
285:   PetscScalar     Evf=user->Evf,A=user->A,B=user->B,cv0=user->cv0;
286:   PetscScalar     *cv_p,*eta_p,*logcv_p,*logcv2_p,*DPsiv_p,*DPsieta_p;
287:   PetscInt        n,i;


291:   VecGetLocalSize(user->cv,&n);
292:   VecGetArray(user->cv,&cv_p);
293:   VecGetArray(user->eta,&eta_p);
294:   VecGetArray(user->logcv,&logcv_p);
295:   VecGetArray(user->logcv2,&logcv2_p);
296:   VecGetArray(user->DPsiv,&DPsiv_p);
297:   VecGetArray(user->DPsieta,&DPsieta_p);

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

301:   VecCopy(user->cv,user->work1);
302:   VecScale(user->work1,-1.0);
303:   VecShift(user->work1,1.0);
304:   Llog(user->work1,user->logcv2);

306:   for (i=0;i<n;i++)
307:   {
308:     DPsiv_p[i] = (eta_p[i]-1.0)*(eta_p[i]-1.0)*(eta_p[i]+1.0)*(eta_p[i]+1.0)*( Evf + logcv_p[i] - logcv2_p[i]) - 2.0*A*(cv_p[i] - cv0)*eta_p[i]*(eta_p[i]+2.0)*(eta_p[i]-1.0)*(eta_p[i]-1.0) + 2.0*B*(cv_p[i] - 1.0)*eta_p[i]*eta_p[i];
309: 
310:     DPsieta_p[i] = 4.0*eta_p[i]*(eta_p[i]-1.0)*(eta_p[i]+1.0)*(Evf*cv_p[i] + cv_p[i]*logcv_p[i] + (1.0-cv_p[i])*logcv2_p[i] ) - A*(cv_p[i] - cv0)*(cv_p[i] - cv0)*(4.0*eta_p[i]*eta_p[i]*eta_p[i] - 6.0*eta_p[i] + 2.0) + 2.0*B*(cv_p[i]-1.0)*(cv_p[i]-1.0)*eta_p[i];
311: 
312:   }

314:   VecRestoreArray(user->cv,&cv_p);
315:   VecRestoreArray(user->eta,&eta_p);
316:   VecGetArray(user->logcv,&logcv_p);
317:   VecGetArray(user->logcv2,&logcv2_p);
318:   VecRestoreArray(user->DPsiv,&DPsiv_p);
319:   VecRestoreArray(user->DPsieta,&DPsieta_p);


322:   return(0);

324: }


329: PetscErrorCode Llog(Vec X, Vec Y)
330: {
331:   PetscErrorCode    ierr;
332:   PetscScalar       *x,*y;
333:   PetscInt          n,i;

336: 
337:   VecGetArray(X,&x);
338:   VecGetArray(Y,&y);
339:   VecGetLocalSize(X,&n);
340:   for (i=0;i<n;i++) {
341:     if (x[i] < 1.0e-12) {
342:       y[i] = log(1.0e-12);
343:     }
344:     else {
345:       y[i] = log(x[i]);
346:     }
347:   }
348: 
349:   return(0);
350: }

354: PetscErrorCode SetInitialGuess(Vec X,AppCtx* user)
355: {
356:   PetscErrorCode    ierr;


359:   PetscInt         n,i,Mda;
360:   PetscScalar           *xx,*cv_p,*wv_p,*eta_p;
361:   PetscViewer      view_out;
362:   /* needed for the void growth case */
363:   PetscScalar       xmid,cv_v=1.0,cv_m=user->Sv*user->cv0,eta_v=1.0,eta_m=0.0,h,lambda;
364:   PetscInt          nele,nen,idx[2];
365:   const PetscInt    *ele;
366:   PetscScalar       x[2];
367:   Vec               coords;
368:   const PetscScalar *_coords;
369:   PetscScalar       xwidth = user->xmax - user->xmin;


373:   VecGetLocalSize(X,&n);
374: 
375: 
376:   DMDAGetInfo(user->da2,PETSC_NULL,&Mda,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
377:   DMDAGetGhostedCoordinates(user->da2,&coords);
378:   VecGetArrayRead(coords,&_coords);

380:   if (user->periodic) {
381:     h = (user->xmax-user->xmin)/Mda;
382:   } else {
383:     h = (user->xmax-user->xmin)/(Mda-1.0);
384:   }

386:   xmid = (user->xmax + user->xmin)/2.0;
387:   lambda = 4.0*h;
388: 
389:   DMDAGetElements(user->da2,&nele,&nen,&ele);
390:   for (i=0;i < nele; i++) {
391:     idx[0] = ele[2*i]; idx[1] = ele[2*i+1];
392: 
393:     x[0] = _coords[idx[0]];
394:     x[1] = _coords[idx[1]];
395: 
396: 
397:     PetscInt k;
398:     PetscScalar vals_DDcv[2],vals_cv[2],vals_eta[2],s,hhr,r;
399:     for (k=0; k < 2 ; k++) {
400:       s = PetscAbs(x[k] - xmid);
401:       if (s <= xwidth*(5.0/64.0)) {
402:         vals_cv[k] = cv_v;
403:         vals_eta[k] = eta_v;
404:         vals_DDcv[k] = 0.0;
405:       } else if (s> xwidth*(5.0/64.0) && s<= xwidth*(7.0/64.0) ) {
406:         //r = (s - xwidth*(6.0/64.0) )/(0.5*lambda);
407:         r = (s - xwidth*(6.0/64.0) )/(xwidth/64.0);
408:         hhr = 0.25*(-r*r*r + 3*r + 2);
409:         vals_cv[k] = cv_m + (1.0 - hhr)*(cv_v - cv_m);
410:         vals_eta[k] = eta_m + (1.0 - hhr)*(eta_v - eta_m);
411:         vals_DDcv[k] = (cv_v - cv_m)*r*6.0/(lambda*lambda);
412:       } else {
413:         vals_cv[k] = cv_m;
414:         vals_eta[k] = eta_m;
415:         vals_DDcv[k] = 0.0;
416:       }
417:     }
418: 
419:     VecSetValuesLocal(user->cv,2,idx,vals_cv,INSERT_VALUES);
420:     VecSetValuesLocal(user->eta,2,idx,vals_eta,INSERT_VALUES);
421:     VecSetValuesLocal(user->work2,2,idx,vals_DDcv,INSERT_VALUES);
422: 
423:   }
424:     DMDARestoreElements(user->da2,&nele,&nen,&ele);
425:     VecRestoreArrayRead(coords,&_coords);

427:     VecAssemblyBegin(user->cv);
428:     VecAssemblyEnd(user->cv);
429:     VecAssemblyBegin(user->eta);
430:     VecAssemblyEnd(user->eta);
431:     VecAssemblyBegin(user->work2);
432:     VecAssemblyEnd(user->work2);

434:   DPsi(user);
435:   VecCopy(user->DPsiv,user->wv);
436:   VecAXPY(user->wv,-2.0*user->kav,user->work2);

438:   VecGetArray(X,&xx);
439:   VecGetArray(user->wv,&wv_p);
440:   VecGetArray(user->cv,&cv_p);
441:   VecGetArray(user->eta,&eta_p);

443:   for (i=0;i<n/3;i++)
444:   {
445:     xx[3*i]=wv_p[i];
446:     xx[3*i+1]=cv_p[i];
447:     xx[3*i+2]=eta_p[i];
448:   }

450:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_initial",FILE_MODE_WRITE,&view_out);
451:   VecView(user->wv,view_out);
452:   VecView(user->cv,view_out);
453:   VecView(user->eta,view_out);
454:   PetscViewerDestroy(&view_out);

456:   VecRestoreArray(X,&xx);
457:   VecRestoreArray(user->wv,&wv_p);
458:   VecRestoreArray(user->cv,&cv_p);
459:   VecRestoreArray(user->eta,&eta_p);
460:   return(0);
461: }

465: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void* ctx)
466: {
468:   AppCtx         *user=(AppCtx*)ctx;
469: 
471:   MatMultAdd(user->M,X,user->q,F);
472:   return(0);
473: }

477: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flg,void *ctx)
478: {
480:   AppCtx         *user=(AppCtx*)ctx;
481: 
483:   *flg = SAME_NONZERO_PATTERN;
484:   MatCopy(user->M,*J,*flg);
485:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
486:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
487:   return(0);
488: }
491: PetscErrorCode SetVariableBounds(DM da,Vec xl,Vec xu)
492: {
494:   PetscScalar    **l,**u;
495:   PetscInt       xs,xm;
496:   PetscInt       i;
497: 
499:   DMDAVecGetArrayDOF(da,xl,&l);
500:   DMDAVecGetArrayDOF(da,xu,&u);
501: 
502:   DMDAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);
503: 
504: 
505:   for(i=xs; i < xs+xm;i++) {
506:     l[i][0] = -SNES_VI_INF;
507:     l[i][1] = 0.0;
508:     l[i][2] = 0.0;
509:     u[i][0] = SNES_VI_INF;
510:     u[i][1] = 1.0;
511:     u[i][2] = 1.0;
512:   }
513: 
514:   DMDAVecRestoreArrayDOF(da,xl,&l);
515:   DMDAVecRestoreArrayDOF(da,xu,&u);
516:   return(0);
517: }

521: PetscErrorCode GetParams(AppCtx* user)
522: {
524:   PetscBool      flg;
525: 
527: 
528:   /* Set default parameters */
529:   user->xmin  = 0.0; user->xmax = 128.0;
530:   user->Mv    = 1.0;
531:   user->L     = 1.0;
532:   user->kaeta = 1.0;
533:   user->kav   = 0.5;
534:   user->Evf   = 9.09;
535:   user->A     = 9.09;
536:   user->B     = 9.09;
537:   user->cv0   = 1.13e-4;
538:   user->Sv    = 500.0;
539:   user->dt    = 1.0e-5;
540:   user->T     = 1.0e-2;
541:   user->graphics = PETSC_TRUE;
542:   user->periodic = PETSC_FALSE;
543:   user->lumpedmass = PETSC_FALSE;
544: 
545:   PetscOptionsGetReal(PETSC_NULL,"-xmin",&user->xmin,&flg);
546:   PetscOptionsGetReal(PETSC_NULL,"-xmax",&user->xmax,&flg);
547:   PetscOptionsGetReal(PETSC_NULL,"-T",&user->T,&flg);
548:   PetscOptionsGetReal(PETSC_NULL,"-dt",&user->dt,&flg);
549:   PetscOptionsBool("-graphics","Contour plot solutions at each timestep\n","None",user->graphics,&user->graphics,&flg);
550:   PetscOptionsBool("-periodic","Use periodic boundary conditions\n","None",user->periodic,&user->periodic,&flg);
551:   PetscOptionsBool("-lumpedmass","Use lumped mass matrix\n","None",user->lumpedmass,&user->lumpedmass,&flg);
552: 
553:   return(0);
554:  }


559: PetscErrorCode SetUpMatrices(AppCtx* user)
560: {
561:   PetscErrorCode    ierr;
562:   PetscInt          nele,nen,i,n;
563:   const PetscInt    *ele;
564:   PetscScalar       dt=user->dt,h;
565: 
566:   PetscInt          idx[2];
567:   PetscScalar       eM_0[2][2],eM_2[2][2];
568:   Mat               M=user->M;
569:   Mat               M_0=user->M_0;
570:   PetscInt          Mda;

572: 


576:   MatGetLocalSize(M,&n,PETSC_NULL);
577:   DMDAGetInfo(user->da1,PETSC_NULL,&Mda,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);

579:   if (user->periodic) {
580:     h = (user->xmax-user->xmin)/Mda;
581:   } else {
582:     h = (user->xmax-user->xmin)/(Mda-1.0);
583:   }
584:   if (user->lumpedmass) {
585:     eM_0[0][0] = h/2.0;
586:     eM_0[1][1] = h/2.0;
587:     eM_0[0][1] = eM_0[1][0] = 0.0;
588:   } else {
589:     eM_0[0][0]=eM_0[1][1]=h/3.0;
590:     eM_0[0][1]=eM_0[1][0]=h/6.0;
591:   }
592:   eM_2[0][0]=eM_2[1][1]=1.0/h;
593:   eM_2[0][1]=eM_2[1][0]=-1.0/h;

595:   /* Get local element info */
596:   DMDAGetElements(user->da1,&nele,&nen,&ele);
597:   for(i=0;i < nele;i++) {
598: 
599:     idx[0] = ele[2*i]; idx[1] = ele[2*i+1];


602:     PetscInt    row,cols[4],r,row_M_0,cols2[2];
603:     PetscScalar vals[4],vals_M_0[2],vals2[2];

605:     for(r=0;r<2;r++) {
606:       row_M_0 = idx[r];
607:       vals_M_0[0]=eM_0[r][0];
608:       vals_M_0[1]=eM_0[r][1];
609:       MatSetValuesLocal(M_0,1,&row_M_0,2,idx,vals_M_0,ADD_VALUES);

611:       row = 3*idx[r];
612:       cols[0] = 3*idx[0];     vals[0] = dt*eM_2[r][0]*user->Mv;
613:       cols[1] = 3*idx[1];     vals[1] = dt*eM_2[r][1]*user->Mv;
614:       cols[2] = 3*idx[0]+1;   vals[2] = eM_0[r][0];
615:       cols[3] = 3*idx[1]+1;   vals[3] = eM_0[r][1];
616: 
617:       /* Insert values in matrix M for 1st dof */
618:       MatSetValuesLocal(M,1,&row,4,cols,vals,ADD_VALUES);
619: 
620:       row = 3*idx[r]+1;
621:       cols[0] = 3*idx[0];     vals[0] = -eM_0[r][0];
622:       cols[1] = 3*idx[1];     vals[1] = -eM_0[r][1];
623:       cols[2] = 3*idx[0]+1;   vals[2] = 2.0*user->kav*eM_2[r][0];
624:       cols[3] = 3*idx[1]+1;   vals[3] = 2.0*user->kav*eM_2[r][1];

626:       /* Insert values in matrix M for 2nd dof */
627:       MatSetValuesLocal(M,1,&row,4,cols,vals,ADD_VALUES);
628: 
629: 
630:       row = 3*idx[r]+2;
631:       cols2[0] = 3*idx[0]+2;   vals2[0] = eM_0[r][0] + user->dt*2.0*user->L*user->kaeta*eM_2[r][0];
632:       cols2[1] = 3*idx[1]+2;   vals2[1] = eM_0[r][1] + user->dt*2.0*user->L*user->kaeta*eM_2[r][1];
633:       MatSetValuesLocal(M,1,&row,2,cols2,vals2,ADD_VALUES);
634:     }
635:   }

637:   MatAssemblyBegin(M_0,MAT_FINAL_ASSEMBLY);
638:   MatAssemblyEnd(M_0,MAT_FINAL_ASSEMBLY);
639: 
640:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
641:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
642: 
643:   DMDARestoreElements(user->da1,&nele,&nen,&ele);
644: 
645: 
646:   return(0);
647: }

651: PetscErrorCode CheckRedundancy(SNES snes, IS act, IS *outact, DM da)
652: {
654:   PetscScalar    **uin,**uout;
655:   Vec            UIN, UOUT;
656:   PetscInt       xs,xm,*outindex;
657:   const PetscInt *index;
658:   PetscInt       k,i,l,n,M,cnt=0;
659: 
661:   DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
662:   DMGetGlobalVector(da,&UIN);
663:   VecSet(UIN,0.0);
664:   DMGetLocalVector(da,&UOUT);
665:   DMDAVecGetArrayDOF(da,UIN,&uin);
666:   ISGetIndices(act,&index);
667:   ISGetLocalSize(act,&n);
668:   for (k=0;k<n;k++){
669:     l = index[k]%5;
670:     i = index[k]/5;
671:     uin[i][l]=1.0;
672:   }
673:   printf("Number of active constraints before applying redundancy %d\n",n);
674:   ISRestoreIndices(act,&index);
675:   DMDAVecRestoreArrayDOF(da,UIN,&uin);
676:   DMGlobalToLocalBegin(da,UIN,INSERT_VALUES,UOUT);
677:   DMGlobalToLocalEnd(da,UIN,INSERT_VALUES,UOUT);
678:   DMDAVecGetArrayDOF(da,UOUT,&uout);
679: 
680:   DMDAGetCorners(da,&xs,PETSC_NULL,PETSC_NULL,&xm,PETSC_NULL,PETSC_NULL);

682:   for(i=xs; i < xs+xm;i++) {
683:     if (uout[i-1][1] && uout[i][1] && uout[i+1][1])
684:              uout[i][0] = 1.0;
685:     if (uout[i-1][3] && uout[i][3] && uout[i+1][3])
686:              uout[i][2] = 1.0;
687:   }

689:   for(i=xs; i < xs+xm;i++) {
690:     for(l=0;l<5;l++) {
691:       if (uout[i][l])
692:         cnt++;
693:     }
694:   }

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

699:   PetscMalloc(cnt*sizeof(PetscInt),&outindex);
700:   cnt = 0;
701: 
702:   for(i=xs; i < xs+xm;i++) {
703:     for(l=0;l<5;l++) {
704:       if (uout[i][l])
705:         outindex[cnt++] = 5*(i)+l;
706:     }
707:   }
708: 
709: 
710:   ISCreateGeneral(PETSC_COMM_WORLD,cnt,outindex,PETSC_OWN_POINTER,outact);
711:   DMDAVecRestoreArrayDOF(da,UOUT,&uout);
712:   DMRestoreGlobalVector(da,&UIN);
713:   DMRestoreLocalVector(da,&UOUT);
714:   return(0);
715: }