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 {
 36:   PetscBool cahnhillard;
 37:   PetscReal kappa;
 38:   PetscInt  energy;
 39:   PetscReal tol;
 40:   PetscReal theta;
 41:   PetscReal theta_c;
 42: } UserCtx;

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

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

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

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

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:      Create timestepping solver context
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 99:   TSCreate(PETSC_COMM_WORLD, &ts);
100:   TSSetDM(ts, da);
101:   TSSetProblemType(ts, TS_NONLINEAR);
102:   TSSetIFunction(ts, NULL, FormFunction, &ctx);
103:   TSSetMaxTime(ts, .02);
104:   TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);

106:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107:      Create matrix data structure; set Jacobian evaluation routine

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

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

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:      Customize nonlinear solver
132:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133:   TSSetType(ts, TSBEULER);

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:      Set initial conditions
137:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:   FormInitialSolution(da, x, ctx.kappa);
139:   TSSetTimeStep(ts, dt);
140:   TSSetSolution(ts, x);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Set runtime options
144:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   TSSetFromOptions(ts);

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      Solve nonlinear system
149:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150:   TSSolve(ts, x);
151:   TSGetStepNumber(ts, &steps);

153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154:      Free work space.  All PETSc objects should be destroyed when they
155:      are no longer needed.
156:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   MatDestroy(&J);
158:   MatFDColoringDestroy(&matfdcoloring);
159:   VecDestroy(&x);
160:   VecDestroy(&r);
161:   TSDestroy(&ts);
162:   DMDestroy(&da);

164:   PetscFinalize();
165:   return 0;
166: }

168: typedef struct {
169:   PetscScalar w, u;
170: } Field;
171: /* ------------------------------------------------------------------- */
172: /*
173:    FormFunction - Evaluates nonlinear function, F(x).

175:    Input Parameters:
176: .  ts - the TS context
177: .  X - input vector
178: .  ptr - optional user-defined context, as set by SNESSetFunction()

180:    Output Parameter:
181: .  F - function vector
182:  */
183: PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec Xdot, Vec F, void *ptr)
184: {
185:   DM          da;
186:   PetscInt    i, Mx, xs, xm;
187:   PetscReal   hx, sx;
188:   PetscScalar r, l;
189:   Field      *x, *xdot, *f;
190:   Vec         localX, localXdot;
191:   UserCtx    *ctx = (UserCtx *)ptr;

193:   TSGetDM(ts, &da);
194:   DMGetLocalVector(da, &localX);
195:   DMGetLocalVector(da, &localXdot);
196:   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);

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

201:   /*
202:      Scatter ghost points to local vector,using the 2-step process
203:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
204:      By placing code between these two statements, computations can be
205:      done while messages are in transition.
206:   */
207:   DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX);
208:   DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX);
209:   DMGlobalToLocalBegin(da, Xdot, INSERT_VALUES, localXdot);
210:   DMGlobalToLocalEnd(da, Xdot, INSERT_VALUES, localXdot);

212:   /*
213:      Get pointers to vector data
214:   */
215:   DMDAVecGetArrayRead(da, localX, &x);
216:   DMDAVecGetArrayRead(da, localXdot, &xdot);
217:   DMDAVecGetArray(da, F, &f);

219:   /*
220:      Get local grid boundaries
221:   */
222:   DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL);

224:   /*
225:      Compute function over the locally owned part of the grid
226:   */
227:   for (i = xs; i < xs + xm; i++) {
228:     f[i].w = x[i].w + ctx->kappa * (x[i - 1].u + x[i + 1].u - 2.0 * x[i].u) * sx;
229:     if (ctx->cahnhillard) {
230:       switch (ctx->energy) {
231:       case 1: /* double well */
232:         f[i].w += -x[i].u * x[i].u * x[i].u + x[i].u;
233:         break;
234:       case 2: /* double obstacle */
235:         f[i].w += x[i].u;
236:         break;
237:       case 3: /* logarithmic */
238:         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;
239:         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;
240:         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;
241:         break;
242:       case 4:
243:         break;
244:       }
245:     }
246:     f[i].u = xdot[i].u - (x[i - 1].w + x[i + 1].w - 2.0 * x[i].w) * sx;
247:     if (ctx->energy == 4) {
248:       f[i].u = xdot[i].u;
249:       /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
250:       r = (1.0 - x[i + 1].u * x[i + 1].u) * (x[i + 2].w - x[i].w) * .5 / hx;
251:       l = (1.0 - x[i - 1].u * x[i - 1].u) * (x[i].w - x[i - 2].w) * .5 / hx;
252:       f[i].u -= (r - l) * .5 / hx;
253:       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;
254:     }
255:   }

257:   /*
258:      Restore vectors
259:   */
260:   DMDAVecRestoreArrayRead(da, localXdot, &xdot);
261:   DMDAVecRestoreArrayRead(da, localX, &x);
262:   DMDAVecRestoreArray(da, F, &f);
263:   DMRestoreLocalVector(da, &localX);
264:   DMRestoreLocalVector(da, &localXdot);
265:   return 0;
266: }

268: /* ------------------------------------------------------------------- */
269: PetscErrorCode FormInitialSolution(DM da, Vec X, PetscReal kappa)
270: {
271:   PetscInt  i, xs, xm, Mx, xgs, xgm;
272:   Field    *x;
273:   PetscReal hx, xx, r, sx;
274:   Vec       Xg;

276:   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);

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

281:   /*
282:      Get pointers to vector data
283:   */
284:   DMCreateLocalVector(da, &Xg);
285:   DMDAVecGetArray(da, Xg, &x);

287:   /*
288:      Get local grid boundaries
289:   */
290:   DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL);
291:   DMDAGetGhostCorners(da, &xgs, NULL, NULL, &xgm, NULL, NULL);

293:   /*
294:      Compute u function over the locally owned part of the grid including ghost points
295:   */
296:   for (i = xgs; i < xgs + xgm; i++) {
297:     xx = i * hx;
298:     r  = PetscSqrtReal((xx - .5) * (xx - .5));
299:     if (r < .125) x[i].u = 1.0;
300:     else x[i].u = -.50;
301:     /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
302:     x[i].w = 0;
303:   }
304:   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;

306:   /*
307:      Restore vectors
308:   */
309:   DMDAVecRestoreArray(da, Xg, &x);

311:   /* Grab only the global part of the vector */
312:   VecSet(X, 0);
313:   DMLocalToGlobalBegin(da, Xg, ADD_VALUES, X);
314:   DMLocalToGlobalEnd(da, Xg, ADD_VALUES, X);
315:   VecDestroy(&Xg);
316:   return 0;
317: }

319: /*TEST

321:    build:
322:      requires: !complex !single

324:    test:
325:      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
326:      requires: x

328: TEST*/