Actual source code: biharmonic3.c


  2: static char help[] = "Solves biharmonic equation in 1d.\n";

  4: /*
  5:   Solves the equation biharmonic equation in split form

  7:     w = -kappa \Delta u
  8:     u_t =  \Delta w
  9:     -1  <= u <= 1
 10:     Periodic boundary conditions

 12: Evolve the biharmonic heat equation with bounds:  (same as biharmonic)
 13: ---------------
 14: ./biharmonic3 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason -ts_type beuler  -da_refine 5 -draw_fields 1 -ts_dt 9.53674e-9

 16:     w = -kappa \Delta u  + u^3  - u
 17:     u_t =  \Delta w
 18:     -1  <= u <= 1
 19:     Periodic boundary conditions

 21: Evolve the Cahn-Hillard equations:
 22: ---------------
 23: ./biharmonic3 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason  -ts_type beuler    -da_refine 6  -draw_fields 1  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard

 25: */
 26: #include <petscdm.h>
 27: #include <petscdmda.h>
 28: #include <petscts.h>
 29: #include <petscdraw.h>

 31: /*
 32:    User-defined routines
 33: */
 34: extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,Vec,void*),FormInitialSolution(DM,Vec,PetscReal);
 35: typedef struct {PetscBool cahnhillard;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta;PetscReal theta_c;} UserCtx;

 37: int main(int argc,char **argv)
 38: {
 39:   TS             ts;                           /* nonlinear solver */
 40:   Vec            x,r;                          /* solution, residual vectors */
 41:   Mat            J;                            /* Jacobian matrix */
 42:   PetscInt       steps,Mx;
 43:   DM             da;
 44:   MatFDColoring  matfdcoloring;
 45:   ISColoring     iscoloring;
 46:   PetscReal      dt;
 47:   PetscReal      vbounds[] = {-100000,100000,-1.1,1.1};
 48:   SNES           snes;
 49:   UserCtx        ctx;

 51:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52:      Initialize program
 53:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 54:   PetscInitialize(&argc,&argv,(char*)0,help);
 55:   ctx.kappa       = 1.0;
 56:   PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
 57:   ctx.cahnhillard = PETSC_FALSE;
 58:   PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
 59:   PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);
 60:   PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);
 61:   ctx.energy      = 1;
 62:   PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
 63:   ctx.tol     = 1.0e-8;
 64:   PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
 65:   ctx.theta   = .001;
 66:   ctx.theta_c = 1.0;
 67:   PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
 68:   PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:      Create distributed array (DMDA) to manage parallel grid and vectors
 72:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 73:   DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da);
 74:   DMSetFromOptions(da);
 75:   DMSetUp(da);
 76:   DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
 77:   DMDASetFieldName(da,1,"Biharmonic heat equation: u");
 78:   DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
 79:   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);

 81:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:      Extract global vectors from DMDA; then duplicate for remaining
 83:      vectors that are the same types
 84:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 85:   DMCreateGlobalVector(da,&x);
 86:   VecDuplicate(x,&r);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:      Create timestepping solver context
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 91:   TSCreate(PETSC_COMM_WORLD,&ts);
 92:   TSSetDM(ts,da);
 93:   TSSetProblemType(ts,TS_NONLINEAR);
 94:   TSSetIFunction(ts,NULL,FormFunction,&ctx);
 95:   TSSetMaxTime(ts,.02);
 96:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:      Create matrix data structure; set Jacobian evaluation routine

101: <     Set Jacobian matrix data structure and default Jacobian evaluation
102:      routine. User can override with:
103:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
104:                 (unless user explicitly sets preconditioner)
105:      -snes_mf_operator : form preconditioning matrix as set by the user,
106:                          but use matrix-free approx for Jacobian-vector
107:                          products within Newton-Krylov method

109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   TSGetSNES(ts,&snes);
111:   SNESSetType(snes,SNESVINEWTONRSLS);
112:   DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
113:   DMSetMatType(da,MATAIJ);
114:   DMCreateMatrix(da,&J);
115:   MatFDColoringCreate(J,iscoloring,&matfdcoloring);
116:   MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
117:   MatFDColoringSetFromOptions(matfdcoloring);
118:   MatFDColoringSetUp(J,iscoloring,matfdcoloring);
119:   ISColoringDestroy(&iscoloring);
120:   SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123:      Customize nonlinear solver
124:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125:   TSSetType(ts,TSBEULER);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:      Set initial conditions
129:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:   FormInitialSolution(da,x,ctx.kappa);
131:   TSSetTimeStep(ts,dt);
132:   TSSetSolution(ts,x);

134:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135:      Set runtime options
136:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137:   TSSetFromOptions(ts);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:      Solve nonlinear system
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142:   TSSolve(ts,x);
143:   TSGetStepNumber(ts,&steps);

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:      Free work space.  All PETSc objects should be destroyed when they
147:      are no longer needed.
148:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149:   MatDestroy(&J);
150:   MatFDColoringDestroy(&matfdcoloring);
151:   VecDestroy(&x);
152:   VecDestroy(&r);
153:   TSDestroy(&ts);
154:   DMDestroy(&da);

156:   PetscFinalize();
157:   return 0;
158: }

160: typedef struct {PetscScalar w,u;} Field;
161: /* ------------------------------------------------------------------- */
162: /*
163:    FormFunction - Evaluates nonlinear function, F(x).

165:    Input Parameters:
166: .  ts - the TS context
167: .  X - input vector
168: .  ptr - optional user-defined context, as set by SNESSetFunction()

170:    Output Parameter:
171: .  F - function vector
172:  */
173: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
174: {
175:   DM             da;
176:   PetscInt       i,Mx,xs,xm;
177:   PetscReal      hx,sx;
178:   PetscScalar    r,l;
179:   Field          *x,*xdot,*f;
180:   Vec            localX,localXdot;
181:   UserCtx        *ctx = (UserCtx*)ptr;

183:   TSGetDM(ts,&da);
184:   DMGetLocalVector(da,&localX);
185:   DMGetLocalVector(da,&localXdot);
186:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

188:   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);

190:   /*
191:      Scatter ghost points to local vector,using the 2-step process
192:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
193:      By placing code between these two statements, computations can be
194:      done while messages are in transition.
195:   */
196:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
197:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
198:   DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
199:   DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);

201:   /*
202:      Get pointers to vector data
203:   */
204:   DMDAVecGetArrayRead(da,localX,&x);
205:   DMDAVecGetArrayRead(da,localXdot,&xdot);
206:   DMDAVecGetArray(da,F,&f);

208:   /*
209:      Get local grid boundaries
210:   */
211:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

213:   /*
214:      Compute function over the locally owned part of the grid
215:   */
216:   for (i=xs; i<xs+xm; i++) {
217:     f[i].w =  x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
218:     if (ctx->cahnhillard) {
219:       switch (ctx->energy) {
220:       case 1: /* double well */
221:         f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
222:         break;
223:       case 2: /* double obstacle */
224:         f[i].w += x[i].u;
225:         break;
226:       case 3: /* logarithmic */
227:         if (x[i].u < -1.0 + 2.0*ctx->tol)      f[i].w += .5*ctx->theta*(-PetscLogScalar(ctx->tol) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
228:         else if (x[i].u > 1.0 - 2.0*ctx->tol)  f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogScalar(ctx->tol)) + ctx->theta_c*x[i].u;
229:         else                                   f[i].w += .5*ctx->theta*(-PetscLogScalar((1.0+x[i].u)/2.0) + PetscLogScalar((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
230:         break;
231:       case 4:
232:         break;
233:       }
234:     }
235:     f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
236:     if (ctx->energy==4) {
237:       f[i].u = xdot[i].u;
238:       /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
239:       r       = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
240:       l       = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
241:       f[i].u -= (r - l)*.5/hx;
242:       f[i].u += 2.0*ctx->theta_c*x[i].u*(x[i+1].u-x[i-1].u)*(x[i+1].u-x[i-1].u)*.25*sx - (ctx->theta - ctx->theta_c*(1-x[i].u*x[i].u))*(x[i+1].u + x[i-1].u - 2.0*x[i].u)*sx;
243:     }
244:   }

246:   /*
247:      Restore vectors
248:   */
249:   DMDAVecRestoreArrayRead(da,localXdot,&xdot);
250:   DMDAVecRestoreArrayRead(da,localX,&x);
251:   DMDAVecRestoreArray(da,F,&f);
252:   DMRestoreLocalVector(da,&localX);
253:   DMRestoreLocalVector(da,&localXdot);
254:   return 0;
255: }

257: /* ------------------------------------------------------------------- */
258: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
259: {
260:   PetscInt       i,xs,xm,Mx,xgs,xgm;
261:   Field          *x;
262:   PetscReal      hx,xx,r,sx;
263:   Vec            Xg;

265:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

267:   hx = 1.0/(PetscReal)Mx;
268:   sx = 1.0/(hx*hx);

270:   /*
271:      Get pointers to vector data
272:   */
273:   DMCreateLocalVector(da,&Xg);
274:   DMDAVecGetArray(da,Xg,&x);

276:   /*
277:      Get local grid boundaries
278:   */
279:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
280:   DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL);

282:   /*
283:      Compute u function over the locally owned part of the grid including ghost points
284:   */
285:   for (i=xgs; i<xgs+xgm; i++) {
286:     xx = i*hx;
287:     r = PetscSqrtReal((xx-.5)*(xx-.5));
288:     if (r < .125) x[i].u = 1.0;
289:     else          x[i].u = -.50;
290:     /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
291:     x[i].w = 0;
292:   }
293:   for (i=xs; i<xs+xm; i++) x[i].w = -kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;

295:   /*
296:      Restore vectors
297:   */
298:   DMDAVecRestoreArray(da,Xg,&x);

300:   /* Grab only the global part of the vector */
301:   VecSet(X,0);
302:   DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X);
303:   DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X);
304:   VecDestroy(&Xg);
305:   return 0;
306: }

308: /*TEST

310:    build:
311:      requires: !complex !single

313:    test:
314:      args: -ts_monitor -snes_monitor  -pc_type lu   -snes_converged_reason  -ts_type beuler  -da_refine 5 -ts_dt 9.53674e-9 -ts_max_steps 50
315:      requires: x

317: TEST*/