Actual source code: ex64.c

petsc-master 2015-02-28
Report Typos and Errors
  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,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);

 62:   /* Get physics and time parameters */
 63:   GetParams(&user);

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

 74:   }
 75:   /* Set Element type (triangular) */
 76:   DMDASetElementType(user.da1,DMDA_ELEMENT_P1);
 77:   DMDASetElementType(user.da2,DMDA_ELEMENT_P1);

 79:   /* Set x and y coordinates */
 80:   DMDASetUniformCoordinates(user.da1,user.xmin,user.xmax,NULL,NULL,NULL,NULL);
 81:   DMDASetUniformCoordinates(user.da2,user.xmin,user.xmax,NULL,NULL,NULL,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);

 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:   DMSetMatType(user.da1,MATAIJ);
102:   DMCreateMatrix(user.da1,&user.M);
103:   /* Get the (usual) mass matrix structure from da2 */
104:   DMSetMatType(user.da2,MATAIJ);
105:   DMCreateMatrix(user.da2,&user.M_0);
106:   /* Form the jacobian matrix and M_0 */
107:   SetUpMatrices(&user);
108:   SetInitialGuess(x,&user);
109:   MatDuplicate(user.M,MAT_DO_NOT_COPY_VALUES,&J);

111:   SNESCreate(PETSC_COMM_WORLD,&snes);
112:   SNESSetDM(snes,user.da1);

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

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

128:   if (user.graphics) {
129:     PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),5,bounds);

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

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


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

144:     DPsi(&user);
145:     VecView(user.DPsiv,view_psi);
146:     VecView(user.DPsieta,view_psi);

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

152:     SNESSolve(snes,NULL,x);
153:     VecView(x,view_out);


156:     if (user.graphics) {
157:       VecView(x,PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD));
158:     }

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

164:     Update_u(x,&user);

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

176:     t = t + user.dt;

178:   }


181:   PetscViewerDestroy(&view_out);
182:   PetscViewerDestroy(&view_p);
183:   PetscViewerDestroy(&view_q);
184:   PetscViewerDestroy(&view_psi);
185:   PetscViewerDestroy(&view_mat);
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.eta);
194:   VecDestroy(&user.DPsiv);
195:   VecDestroy(&user.DPsieta);
196:   VecDestroy(&user.logcv);
197:   VecDestroy(&user.logcv2);
198:   VecDestroy(&user.work1);
199:   VecDestroy(&user.work2);
200:   MatDestroy(&user.M);
201:   MatDestroy(&user.M_0);
202:   DMDestroy(&user.da1);
203:   DMDestroy(&user.da1_clone);
204:   DMDestroy(&user.da2);
205:   SNESDestroy(&snes);
206:   PetscFinalize();
207:   return 0;
208: }

212: PetscErrorCode Update_u(Vec X,AppCtx *user)
213: {
214:   PetscErrorCode    ierr;
215:   PetscInt          i,n;
216:   PetscScalar       *wv_p,*cv_p,*eta_p;
217:   const PetscScalar *xx;
218: 
220:   VecGetLocalSize(user->wv,&n);

222:   VecGetArrayRead(X,&xx);
223:   VecGetArray(user->wv,&wv_p);
224:   VecGetArray(user->cv,&cv_p);
225:   VecGetArray(user->eta,&eta_p);


228:   for (i=0; i<n; i++) {
229:     wv_p[i]  = xx[3*i];
230:     cv_p[i]  = xx[3*i+1];
231:     eta_p[i] = xx[3*i+2];
232:   }
233:   VecRestoreArrayRead(X,&xx);
234:   VecRestoreArray(user->wv,&wv_p);
235:   VecRestoreArray(user->cv,&cv_p);
236:   VecRestoreArray(user->eta,&eta_p);
237:   return(0);
238: }

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

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++) q_p[3*i]=w1[i];

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

261:   VecCopy(user->DPsieta,user->work1);
262:   VecScale(user->work1,user->L*user->dt);
263:   VecAXPY(user->work1,-1.0,user->eta);
264:   MatMult(user->M_0,user->work1,user->work2);
265:   for (i=0; i<n; i++) q_p[3*i+2]=w2[i];

267:   VecRestoreArray(user->q,&q_p);
268:   VecRestoreArray(user->work1,&w1);
269:   VecRestoreArray(user->work2,&w2);
270:   return(0);
271: }

275: PetscErrorCode DPsi(AppCtx *user)
276: {
278:   PetscScalar    Evf=user->Evf,A=user->A,B=user->B,cv0=user->cv0;
279:   PetscScalar    *cv_p,*eta_p,*logcv_p,*logcv2_p,*DPsiv_p,*DPsieta_p;
280:   PetscInt       n,i;

283:   VecGetLocalSize(user->cv,&n);
284:   VecGetArray(user->cv,&cv_p);
285:   VecGetArray(user->eta,&eta_p);
286:   VecGetArray(user->logcv,&logcv_p);
287:   VecGetArray(user->logcv2,&logcv2_p);
288:   VecGetArray(user->DPsiv,&DPsiv_p);
289:   VecGetArray(user->DPsieta,&DPsieta_p);

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

293:   VecCopy(user->cv,user->work1);
294:   VecScale(user->work1,-1.0);
295:   VecShift(user->work1,1.0);
296:   Llog(user->work1,user->logcv2);

298:   for (i=0; i<n; i++)
299:   {
300:     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];

302:     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];

304:   }

306:   VecRestoreArray(user->cv,&cv_p);
307:   VecRestoreArray(user->eta,&eta_p);
308:   VecGetArray(user->logcv,&logcv_p);
309:   VecGetArray(user->logcv2,&logcv2_p);
310:   VecRestoreArray(user->DPsiv,&DPsiv_p);
311:   VecRestoreArray(user->DPsieta,&DPsieta_p);
312:   return(0);

314: }


319: PetscErrorCode Llog(Vec X, Vec Y)
320: {
321:   PetscErrorCode    ierr;
322:   PetscScalar       *y;
323:   const PetscScalar *x;
324:   PetscInt          n,i;

327:   VecGetArrayRead(X,&x);
328:   VecGetArray(Y,&y);
329:   VecGetLocalSize(X,&n);
330:   for (i=0; i<n; i++) {
331:     if (x[i] < 1.0e-12) y[i] = log(1.0e-12);
332:     else y[i] = log(x[i]);
333:   }
334:   VecRestoreArrayRead(X,&x);
335:   VecRestoreArray(Y,&y);
336:   return(0);
337: }

341: PetscErrorCode SetInitialGuess(Vec X,AppCtx *user)
342: {

345:   PetscInt          n,i,Mda;
346:   PetscScalar       *xx,*cv_p,*wv_p,*eta_p;
347:   PetscViewer       view_out;

349:   /* needed for the void growth case */
350:   PetscScalar       xmid,cv_v=1.0,cv_m=user->Sv*user->cv0,eta_v=1.0,eta_m=0.0,h,lambda;
351:   PetscInt          nele,nen,idx[2];
352:   const PetscInt    *ele;
353:   PetscScalar       x[2];
354:   Vec               coords;
355:   const PetscScalar *_coords;
356:   PetscScalar       xwidth = user->xmax - user->xmin;

359:   VecGetLocalSize(X,&n);

361:   DMDAGetInfo(user->da2,NULL,&Mda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
362:   DMGetCoordinatesLocal(user->da2,&coords);
363:   VecGetArrayRead(coords,&_coords);

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

368:   xmid   = (user->xmax + user->xmin)/2.0;
369:   lambda = 4.0*h;

371:   DMDAGetElements(user->da2,&nele,&nen,&ele);
372:   for (i=0; i < nele; i++) {
373:     idx[0] = ele[2*i]; idx[1] = ele[2*i+1];

375:     x[0] = _coords[idx[0]];
376:     x[1] = _coords[idx[1]];


379:     PetscInt    k;
380:     PetscScalar vals_DDcv[2],vals_cv[2],vals_eta[2],s,hhr,r;
381:     for (k=0; k < 2; k++) {
382:       s = PetscAbs(x[k] - xmid);
383:       if (s <= xwidth*(5.0/64.0)) {
384:         vals_cv[k]   = cv_v;
385:         vals_eta[k]  = eta_v;
386:         vals_DDcv[k] = 0.0;
387:       } else if (s> xwidth*(5.0/64.0) && s<= xwidth*(7.0/64.0)) {
388:         /*r = (s - xwidth*(6.0/64.0))/(0.5*lambda);*/
389:         r            = (s - xwidth*(6.0/64.0))/(xwidth/64.0);
390:         hhr          = 0.25*(-r*r*r + 3*r + 2);
391:         vals_cv[k]   = cv_m + (1.0 - hhr)*(cv_v - cv_m);
392:         vals_eta[k]  = eta_m + (1.0 - hhr)*(eta_v - eta_m);
393:         vals_DDcv[k] = (cv_v - cv_m)*r*6.0/(lambda*lambda);
394:       } else {
395:         vals_cv[k]   = cv_m;
396:         vals_eta[k]  = eta_m;
397:         vals_DDcv[k] = 0.0;
398:       }
399:     }

401:     VecSetValuesLocal(user->cv,2,idx,vals_cv,INSERT_VALUES);
402:     VecSetValuesLocal(user->eta,2,idx,vals_eta,INSERT_VALUES);
403:     VecSetValuesLocal(user->work2,2,idx,vals_DDcv,INSERT_VALUES);

405:   }
406:   DMDARestoreElements(user->da2,&nele,&nen,&ele);
407:   VecRestoreArrayRead(coords,&_coords);

409:   VecAssemblyBegin(user->cv);
410:   VecAssemblyEnd(user->cv);
411:   VecAssemblyBegin(user->eta);
412:   VecAssemblyEnd(user->eta);
413:   VecAssemblyBegin(user->work2);
414:   VecAssemblyEnd(user->work2);

416:   DPsi(user);
417:   VecCopy(user->DPsiv,user->wv);
418:   VecAXPY(user->wv,-2.0*user->kav,user->work2);

420:   VecGetArray(X,&xx);
421:   VecGetArray(user->wv,&wv_p);
422:   VecGetArray(user->cv,&cv_p);
423:   VecGetArray(user->eta,&eta_p);

425:   for (i=0; i<n/3; i++) {
426:     xx[3*i]  =wv_p[i];
427:     xx[3*i+1]=cv_p[i];
428:     xx[3*i+2]=eta_p[i];
429:   }

431:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"file_initial",FILE_MODE_WRITE,&view_out);
432:   VecView(user->wv,view_out);
433:   VecView(user->cv,view_out);
434:   VecView(user->eta,view_out);
435:   PetscViewerDestroy(&view_out);

437:   VecRestoreArray(X,&xx);
438:   VecRestoreArray(user->wv,&wv_p);
439:   VecRestoreArray(user->cv,&cv_p);
440:   VecRestoreArray(user->eta,&eta_p);
441:   return(0);
442: }

446: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ctx)
447: {
449:   AppCtx         *user=(AppCtx*)ctx;

452:   MatMultAdd(user->M,X,user->q,F);
453:   return(0);
454: }

458: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
459: {
461:   AppCtx         *user=(AppCtx*)ctx;

464:   MatCopy(user->M,J,*flg);
465:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
466:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
467:   return(0);
468: }
471: PetscErrorCode SetVariableBounds(DM da,Vec xl,Vec xu)
472: {
474:   PetscScalar    **l,**u;
475:   PetscInt       xs,xm;
476:   PetscInt       i;

479:   DMDAVecGetArrayDOF(da,xl,&l);
480:   DMDAVecGetArrayDOF(da,xu,&u);

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


485:   for (i=xs; i < xs+xm; i++) {
486:     l[i][0] = -PETSC_INFINITY;
487:     l[i][1] = 0.0;
488:     l[i][2] = 0.0;
489:     u[i][0] = PETSC_INFINITY;
490:     u[i][1] = 1.0;
491:     u[i][2] = 1.0;
492:   }

494:   DMDAVecRestoreArrayDOF(da,xl,&l);
495:   DMDAVecRestoreArrayDOF(da,xu,&u);
496:   return(0);
497: }

501: PetscErrorCode GetParams(AppCtx *user)
502: {
504:   PetscBool      flg;

507:   /* Set default parameters */
508:   user->xmin       = 0.0; user->xmax = 128.0;
509:   user->Mv         = 1.0;
510:   user->L          = 1.0;
511:   user->kaeta      = 1.0;
512:   user->kav        = 0.5;
513:   user->Evf        = 9.09;
514:   user->A          = 9.09;
515:   user->B          = 9.09;
516:   user->cv0        = 1.13e-4;
517:   user->Sv         = 500.0;
518:   user->dt         = 1.0e-5;
519:   user->T          = 1.0e-2;
520:   user->graphics   = PETSC_TRUE;
521:   user->periodic   = PETSC_FALSE;
522:   user->lumpedmass = PETSC_FALSE;

524:   PetscOptionsGetReal(NULL,"-xmin",&user->xmin,&flg);
525:   PetscOptionsGetReal(NULL,"-xmax",&user->xmax,&flg);
526:   PetscOptionsGetReal(NULL,"-T",&user->T,&flg);
527:   PetscOptionsGetReal(NULL,"-dt",&user->dt,&flg);
528:   PetscOptionsBool("-graphics","Contour plot solutions at each timestep\n","None",user->graphics,&user->graphics,&flg);
529:   PetscOptionsBool("-periodic","Use periodic boundary conditions\n","None",user->periodic,&user->periodic,&flg);
530:   PetscOptionsBool("-lumpedmass","Use lumped mass matrix\n","None",user->lumpedmass,&user->lumpedmass,&flg);
531:   return(0);
532: }


537: PetscErrorCode SetUpMatrices(AppCtx *user)
538: {
540:   PetscInt       nele,nen,i,n;
541:   const PetscInt *ele;
542:   PetscScalar    dt=user->dt,h;

544:   PetscInt    idx[2];
545:   PetscScalar eM_0[2][2],eM_2[2][2];
546:   Mat         M  =user->M;
547:   Mat         M_0=user->M_0;
548:   PetscInt    Mda;


552:   MatGetLocalSize(M,&n,NULL);
553:   DMDAGetInfo(user->da1,NULL,&Mda,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

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

558:   if (user->lumpedmass) {
559:     eM_0[0][0] = h/2.0;
560:     eM_0[1][1] = h/2.0;
561:     eM_0[0][1] = eM_0[1][0] = 0.0;
562:   } else {
563:     eM_0[0][0]=eM_0[1][1]=h/3.0;
564:     eM_0[0][1]=eM_0[1][0]=h/6.0;
565:   }
566:   eM_2[0][0]=eM_2[1][1]=1.0/h;
567:   eM_2[0][1]=eM_2[1][0]=-1.0/h;

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

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

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

578:     for (r=0; r<2; r++) {
579:       row_M_0    = idx[r];
580:       vals_M_0[0]=eM_0[r][0];
581:       vals_M_0[1]=eM_0[r][1];

583:       MatSetValuesLocal(M_0,1,&row_M_0,2,idx,vals_M_0,ADD_VALUES);

585:       row     = 3*idx[r];
586:       cols[0] = 3*idx[0];     vals[0] = dt*eM_2[r][0]*user->Mv;
587:       cols[1] = 3*idx[1];     vals[1] = dt*eM_2[r][1]*user->Mv;
588:       cols[2] = 3*idx[0]+1;   vals[2] = eM_0[r][0];
589:       cols[3] = 3*idx[1]+1;   vals[3] = eM_0[r][1];

591:       /* Insert values in matrix M for 1st dof */
592:       MatSetValuesLocal(M,1,&row,4,cols,vals,ADD_VALUES);

594:       row     = 3*idx[r]+1;
595:       cols[0] = 3*idx[0];     vals[0] = -eM_0[r][0];
596:       cols[1] = 3*idx[1];     vals[1] = -eM_0[r][1];
597:       cols[2] = 3*idx[0]+1;   vals[2] = 2.0*user->kav*eM_2[r][0];
598:       cols[3] = 3*idx[1]+1;   vals[3] = 2.0*user->kav*eM_2[r][1];

600:       /* Insert values in matrix M for 2nd dof */
601:       MatSetValuesLocal(M,1,&row,4,cols,vals,ADD_VALUES);


604:       row      = 3*idx[r]+2;
605:       cols2[0] = 3*idx[0]+2;   vals2[0] = eM_0[r][0] + user->dt*2.0*user->L*user->kaeta*eM_2[r][0];
606:       cols2[1] = 3*idx[1]+2;   vals2[1] = eM_0[r][1] + user->dt*2.0*user->L*user->kaeta*eM_2[r][1];

608:       MatSetValuesLocal(M,1,&row,2,cols2,vals2,ADD_VALUES);
609:     }
610:   }

612:   MatAssemblyBegin(M_0,MAT_FINAL_ASSEMBLY);
613:   MatAssemblyEnd(M_0,MAT_FINAL_ASSEMBLY);

615:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
616:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

618:   DMDARestoreElements(user->da1,&nele,&nen,&ele);
619:   return(0);
620: }

624: PetscErrorCode CheckRedundancy(SNES snes, IS act, IS *outact, DM da)
625: {
627:   PetscScalar    **uin,**uout;
628:   Vec            UIN, UOUT;
629:   PetscInt       xs,xm,*outindex;
630:   const PetscInt *index;
631:   PetscInt       k,i,l,n,M,cnt=0;

634:   DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
635:   DMGetGlobalVector(da,&UIN);
636:   VecSet(UIN,0.0);
637:   DMGetLocalVector(da,&UOUT);
638:   DMDAVecGetArrayDOF(da,UIN,&uin);
639:   ISGetIndices(act,&index);
640:   ISGetLocalSize(act,&n);
641:   for (k=0; k<n; k++) {
642:     l = index[k]%5;
643:     i = index[k]/5;

645:     uin[i][l]=1.0;
646:   }
647:   printf("Number of active constraints before applying redundancy %d\n",n);
648:   ISRestoreIndices(act,&index);
649:   DMDAVecRestoreArrayDOF(da,UIN,&uin);
650:   DMGlobalToLocalBegin(da,UIN,INSERT_VALUES,UOUT);
651:   DMGlobalToLocalEnd(da,UIN,INSERT_VALUES,UOUT);
652:   DMDAVecGetArrayDOF(da,UOUT,&uout);

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

656:   for (i=xs; i < xs+xm;i++) {
657:     if (uout[i-1][1] && uout[i][1] && uout[i+1][1]) uout[i][0] = 1.0;
658:     if (uout[i-1][3] && uout[i][3] && uout[i+1][3]) uout[i][2] = 1.0;
659:   }

661:   for (i=xs; i < xs+xm; i++) {
662:     for (l=0; l < 5; l++) {
663:       if (uout[i][l]) cnt++;
664:     }
665:   }

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

669:   PetscMalloc1(cnt,&outindex);
670:   cnt  = 0;

672:   for (i=xs; i < xs+xm;i++) {
673:     for (l=0;l<5;l++) {
674:       if (uout[i][l]) outindex[cnt++] = 5*(i)+l;
675:     }
676:   }


679:   ISCreateGeneral(PETSC_COMM_WORLD,cnt,outindex,PETSC_OWN_POINTER,outact);
680:   DMDAVecRestoreArrayDOF(da,UOUT,&uout);
681:   DMRestoreGlobalVector(da,&UIN);
682:   DMRestoreLocalVector(da,&UOUT);
683:   return(0);
684: }