Actual source code: biharmonic3.c

petsc-3.7.3 2016-07-24
Report Typos and Errors
  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: ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason --wait   -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: ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason  --wait   -ts_type beuler    -da_refine 6 -vi  -draw_fields 1  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard


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

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

 40: int main(int argc,char **argv)
 41: {
 42:   TS             ts;                           /* nonlinear solver */
 43:   Vec            x,r;                          /* solution, residual vectors */
 44:   Mat            J;                            /* Jacobian matrix */
 45:   PetscInt       steps,Mx,maxsteps = 10000000;
 47:   DM             da;
 48:   MatFDColoring  matfdcoloring;
 49:   ISColoring     iscoloring;
 50:   PetscReal      dt;
 51:   PetscReal      vbounds[] = {-100000,100000,-1.1,1.1};
 52:   PetscBool      wait;
 53:   Vec            ul,uh;
 54:   SNES           snes;
 55:   UserCtx        ctx;

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

 77:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 78:      Create distributed array (DMDA) to manage parallel grid and vectors
 79:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 80:   DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,2,2,NULL,&da);
 81:   DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
 82:   DMDASetFieldName(da,1,"Biharmonic heat equation: u");
 83:   DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
 84:   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);

 86:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:      Extract global vectors from DMDA; then duplicate for remaining
 88:      vectors that are the same types
 89:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:   DMCreateGlobalVector(da,&x);
 91:   VecDuplicate(x,&r);

 93:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 94:      Create timestepping solver context
 95:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 96:   TSCreate(PETSC_COMM_WORLD,&ts);
 97:   TSSetDM(ts,da);
 98:   TSSetProblemType(ts,TS_NONLINEAR);
 99:   TSSetIFunction(ts,NULL,FormFunction,&ctx);
100:   TSSetDuration(ts,maxsteps,.02);
101:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:      Create matrix data structure; set Jacobian evaluation routine

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

114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115:   TSGetSNES(ts,&snes);
116:   DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
117:   DMSetMatType(da,MATAIJ);
118:   DMCreateMatrix(da,&J);
119:   MatFDColoringCreate(J,iscoloring,&matfdcoloring);
120:   ISColoringDestroy(&iscoloring);
121:   MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
122:   MatFDColoringSetFromOptions(matfdcoloring);
123:   MatFDColoringSetUp(J,iscoloring,matfdcoloring);
124:   SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);

126:   {
127:     VecDuplicate(x,&ul);
128:     VecDuplicate(x,&uh);
129:     VecStrideSet(ul,0,PETSC_NINFINITY);
130:     VecStrideSet(ul,1,-1.0);
131:     VecStrideSet(uh,0,PETSC_INFINITY);
132:     VecStrideSet(uh,1,1.0);
133:     TSVISetVariableBounds(ts,ul,uh);
134:   }

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:      Customize nonlinear solver
138:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:   TSSetType(ts,TSBEULER);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:      Set initial conditions
143:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144:   FormInitialSolution(da,x,ctx.kappa);
145:   TSSetInitialTimeStep(ts,0.0,dt);
146:   TSSetSolution(ts,x);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:      Set runtime options
150:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151:   TSSetFromOptions(ts);

153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154:      Solve nonlinear system
155:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156:   TSSolve(ts,x);
157:   wait = PETSC_FALSE;
158:   PetscOptionsGetBool(NULL,NULL,"-wait",&wait,NULL);
159:   if (wait) {
160:     PetscSleep(-1);
161:   }
162:   TSGetTimeStepNumber(ts,&steps);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Free work space.  All PETSc objects should be destroyed when they
166:      are no longer needed.
167:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168:   {
169:     VecDestroy(&ul);
170:     VecDestroy(&uh);
171:   }
172:   MatDestroy(&J);
173:   MatFDColoringDestroy(&matfdcoloring);
174:   VecDestroy(&x);
175:   VecDestroy(&r);
176:   TSDestroy(&ts);
177:   DMDestroy(&da);

179:   PetscFinalize();
180:   return(0);
181: }

183: typedef struct {PetscScalar w,u;} Field;
184: /* ------------------------------------------------------------------- */
187: /*
188:    FormFunction - Evaluates nonlinear function, F(x).

190:    Input Parameters:
191: .  ts - the TS context
192: .  X - input vector
193: .  ptr - optional user-defined context, as set by SNESSetFunction()

195:    Output Parameter:
196: .  F - function vector
197:  */
198: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
199: {
200:   DM             da;
202:   PetscInt       i,Mx,xs,xm;
203:   PetscReal      hx,sx;
204:   PetscScalar    r,l;
205:   Field          *x,*xdot,*f;
206:   Vec            localX,localXdot;
207:   UserCtx        *ctx = (UserCtx*)ptr;

210:   TSGetDM(ts,&da);
211:   DMGetLocalVector(da,&localX);
212:   DMGetLocalVector(da,&localXdot);
213:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
214:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

218:   /*
219:      Scatter ghost points to local vector,using the 2-step process
220:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
221:      By placing code between these two statements, computations can be
222:      done while messages are in transition.
223:   */
224:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
225:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
226:   DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
227:   DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);

229:   /*
230:      Get pointers to vector data
231:   */
232:   DMDAVecGetArrayRead(da,localX,&x);
233:   DMDAVecGetArrayRead(da,localXdot,&xdot);
234:   DMDAVecGetArray(da,F,&f);

236:   /*
237:      Get local grid boundaries
238:   */
239:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

241:   /*
242:      Compute function over the locally owned part of the grid
243:   */
244:   for (i=xs; i<xs+xm; i++) {
245:     f[i].w =  x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
246:     if (ctx->cahnhillard) {
247:       switch (ctx->energy) {
248:       case 1: /* double well */
249:         f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
250:         break;
251:       case 2: /* double obstacle */
252:         f[i].w += x[i].u;
253:         break;
254:       case 3: /* logarithmic */
255:         if (x[i].u < -1.0 + 2.0*ctx->tol)      f[i].w += .5*ctx->theta*(-log(ctx->tol) + log((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
256:         else if (x[i].u > 1.0 - 2.0*ctx->tol)  f[i].w += .5*ctx->theta*(-log((1.0+x[i].u)/2.0) + log(ctx->tol)) + ctx->theta_c*x[i].u;
257:         else                                   f[i].w += .5*ctx->theta*(-log((1.0+x[i].u)/2.0) + log((1.0-x[i].u)/2.0)) + ctx->theta_c*x[i].u;
258:         break;
259:       case 4:
260:         break;
261:       }
262:     }
263:     f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
264:     if (ctx->energy==4) {
265:       f[i].u = xdot[i].u;
266:       /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
267:       r       = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
268:       l       = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
269:       f[i].u -= (r - l)*.5/hx;
270:       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;
271:     }
272:   }

274:   /*
275:      Restore vectors
276:   */
277:   DMDAVecRestoreArrayRead(da,localXdot,&xdot);
278:   DMDAVecRestoreArrayRead(da,localX,&x);
279:   DMDAVecRestoreArray(da,F,&f);
280:   DMRestoreLocalVector(da,&localX);
281:   DMRestoreLocalVector(da,&localXdot);
282:   return(0);
283: }

285: /* ------------------------------------------------------------------- */
288: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
289: {
291:   PetscInt       i,xs,xm,Mx;
292:   Field          *x;
293:   PetscReal      hx,xx,r,sx;

296:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
297:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

302:   /*
303:      Get pointers to vector data
304:   */
305:   DMDAVecGetArray(da,X,&x);

307:   /*
308:      Get local grid boundaries
309:   */
310:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

312:   /*
313:      Compute function over the locally owned part of the grid
314:   */
315:   for (i=xs; i<xs+xm; i++) {
316:     xx = i*hx;
317:     r  = PetscSqrtScalar((xx-.5)*(xx-.5));
318:     if (r < .125) x[i].u = 1.0;
319:     else          x[i].u = -.50;
320:     /*  u[i] = PetscPowScalar(x - .5,4.0); */
321:   }
322:   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;

324:   /*
325:      Restore vectors
326:   */
327:   DMDAVecRestoreArray(da,X,&x);
328:   return(0);
329: }