Actual source code: biharmonic2.c

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

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

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

 11: Evolve the biharmonic heat equation with bounds:  (same as biharmonic)
 12: ---------------
 13: ./biharmonic2 -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

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

 20: Evolve the Cahn-Hillard equations: (this fails after a few timesteps 12/17/2017)
 21: ---------------
 22: ./biharmonic2 -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

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

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

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

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:      Initialize program
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   PetscFunctionBeginUser;
 61:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
 62:   ctx.kappa = 1.0;
 63:   PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &ctx.kappa, NULL));
 64:   ctx.cahnhillard = PETSC_FALSE;

 66:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-cahn-hillard", &ctx.cahnhillard, NULL));
 67:   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 2, vbounds));
 68:   PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 600, 600));
 69:   ctx.energy = 1;
 70:   /*PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));*/
 71:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-energy", &ctx.energy, NULL));
 72:   ctx.tol = 1.0e-8;
 73:   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &ctx.tol, NULL));
 74:   ctx.theta   = .001;
 75:   ctx.theta_c = 1.0;
 76:   PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta", &ctx.theta, NULL));
 77:   PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta_c", &ctx.theta_c, NULL));

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

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

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

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

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

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

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

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

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

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

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

164:   PetscCall(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:   Field    *x, *xdot, *f;
189:   Vec       localX, localXdot;
190:   UserCtx  *ctx = (UserCtx *)ptr;

192:   PetscFunctionBegin;
193:   PetscCall(TSGetDM(ts, &da));
194:   PetscCall(DMGetLocalVector(da, &localX));
195:   PetscCall(DMGetLocalVector(da, &localXdot));
196:   PetscCall(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:   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
208:   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
209:   PetscCall(DMGlobalToLocalBegin(da, Xdot, INSERT_VALUES, localXdot));
210:   PetscCall(DMGlobalToLocalEnd(da, Xdot, INSERT_VALUES, localXdot));

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

219:   /*
220:      Get local grid boundaries
221:   */
222:   PetscCall(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 (PetscRealPart(x[i].u) < -1.0 + 2.0 * ctx->tol) f[i].w += .5 * ctx->theta * (-PetscLogReal(ctx->tol) + PetscLogScalar((1.0 - x[i].u) / 2.0)) + ctx->theta_c * x[i].u;
239:         else if (PetscRealPart(x[i].u) > 1.0 - 2.0 * ctx->tol) f[i].w += .5 * ctx->theta * (-PetscLogScalar((1.0 + x[i].u) / 2.0) + PetscLogReal(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:       }
243:     }
244:     f[i].u = xdot[i].u - (x[i - 1].w + x[i + 1].w - 2.0 * x[i].w) * sx;
245:   }

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

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

266:   PetscFunctionBegin;
267:   PetscCall(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));

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

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

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

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

297:   /*
298:      Restore vectors
299:   */
300:   PetscCall(DMDAVecRestoreArray(da, Xg, &x));

302:   /* Grab only the global part of the vector */
303:   PetscCall(VecSet(X, 0));
304:   PetscCall(DMLocalToGlobalBegin(da, Xg, ADD_VALUES, X));
305:   PetscCall(DMLocalToGlobalEnd(da, Xg, ADD_VALUES, X));
306:   PetscCall(VecDestroy(&Xg));
307:   PetscFunctionReturn(PETSC_SUCCESS);
308: }

310: /*TEST

312:    build:
313:      requires: !complex !single

315:    test:
316:      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
317:      requires: x

319: TEST*/