1 2 static char help[] = "Solves biharmonic equation in 1d.\n"; 3 4 /* 5 Solves the equation biharmonic equation in split form 6 7 w = -kappa \Delta u 8 u_t = \Delta w 9 -1 <= u <= 1 10 Periodic boundary conditions 11 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 15 16 w = -kappa \Delta u + u^3 - u 17 u_t = \Delta w 18 -1 <= u <= 1 19 Periodic boundary conditions 20 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 24 25 */ 26 #include <petscdm.h> 27 #include <petscdmda.h> 28 #include <petscts.h> 29 #include <petscdraw.h> 30 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; 43 44 int main(int argc, char **argv) { 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; 56 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; 65 PetscCall(PetscOptionsGetBool(NULL, NULL, "-cahn-hillard", &ctx.cahnhillard, NULL)); 66 PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 2, vbounds)); 67 PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 600, 600)); 68 ctx.energy = 1; 69 PetscCall(PetscOptionsGetInt(NULL, NULL, "-energy", &ctx.energy, NULL)); 70 ctx.tol = 1.0e-8; 71 PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &ctx.tol, NULL)); 72 ctx.theta = .001; 73 ctx.theta_c = 1.0; 74 PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta", &ctx.theta, NULL)); 75 PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta_c", &ctx.theta_c, NULL)); 76 77 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78 Create distributed array (DMDA) to manage parallel grid and vectors 79 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 80 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, 2, 2, NULL, &da)); 81 PetscCall(DMSetFromOptions(da)); 82 PetscCall(DMSetUp(da)); 83 PetscCall(DMDASetFieldName(da, 0, "Biharmonic heat equation: w = -kappa*u_xx")); 84 PetscCall(DMDASetFieldName(da, 1, "Biharmonic heat equation: u")); 85 PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 86 dt = 1.0 / (10. * ctx.kappa * Mx * Mx * Mx * Mx); 87 88 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89 Extract global vectors from DMDA; then duplicate for remaining 90 vectors that are the same types 91 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 92 PetscCall(DMCreateGlobalVector(da, &x)); 93 PetscCall(VecDuplicate(x, &r)); 94 95 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96 Create timestepping solver context 97 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 99 PetscCall(TSSetDM(ts, da)); 100 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 101 PetscCall(TSSetIFunction(ts, NULL, FormFunction, &ctx)); 102 PetscCall(TSSetMaxTime(ts, .02)); 103 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 104 105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106 Create matrix data structure; set Jacobian evaluation routine 107 108 < Set Jacobian matrix data structure and default Jacobian evaluation 109 routine. User can override with: 110 -snes_mf : matrix-free Newton-Krylov method with no preconditioning 111 (unless user explicitly sets preconditioner) 112 -snes_mf_operator : form preconditioning matrix as set by the user, 113 but use matrix-free approx for Jacobian-vector 114 products within Newton-Krylov method 115 116 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 117 PetscCall(TSGetSNES(ts, &snes)); 118 PetscCall(SNESSetType(snes, SNESVINEWTONRSLS)); 119 PetscCall(DMCreateColoring(da, IS_COLORING_GLOBAL, &iscoloring)); 120 PetscCall(DMSetMatType(da, MATAIJ)); 121 PetscCall(DMCreateMatrix(da, &J)); 122 PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring)); 123 PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))SNESTSFormFunction, ts)); 124 PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 125 PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring)); 126 PetscCall(ISColoringDestroy(&iscoloring)); 127 PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring)); 128 129 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130 Customize nonlinear solver 131 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132 PetscCall(TSSetType(ts, TSBEULER)); 133 134 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135 Set initial conditions 136 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137 PetscCall(FormInitialSolution(da, x, ctx.kappa)); 138 PetscCall(TSSetTimeStep(ts, dt)); 139 PetscCall(TSSetSolution(ts, x)); 140 141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142 Set runtime options 143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144 PetscCall(TSSetFromOptions(ts)); 145 146 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147 Solve nonlinear system 148 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 149 PetscCall(TSSolve(ts, x)); 150 PetscCall(TSGetStepNumber(ts, &steps)); 151 152 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153 Free work space. All PETSc objects should be destroyed when they 154 are no longer needed. 155 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 156 PetscCall(MatDestroy(&J)); 157 PetscCall(MatFDColoringDestroy(&matfdcoloring)); 158 PetscCall(VecDestroy(&x)); 159 PetscCall(VecDestroy(&r)); 160 PetscCall(TSDestroy(&ts)); 161 PetscCall(DMDestroy(&da)); 162 163 PetscCall(PetscFinalize()); 164 return 0; 165 } 166 167 typedef struct { 168 PetscScalar w, u; 169 } Field; 170 /* ------------------------------------------------------------------- */ 171 /* 172 FormFunction - Evaluates nonlinear function, F(x). 173 174 Input Parameters: 175 . ts - the TS context 176 . X - input vector 177 . ptr - optional user-defined context, as set by SNESSetFunction() 178 179 Output Parameter: 180 . F - function vector 181 */ 182 PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec Xdot, Vec F, void *ptr) { 183 DM da; 184 PetscInt i, Mx, xs, xm; 185 PetscReal hx, sx; 186 PetscScalar r, l; 187 Field *x, *xdot, *f; 188 Vec localX, localXdot; 189 UserCtx *ctx = (UserCtx *)ptr; 190 191 PetscFunctionBegin; 192 PetscCall(TSGetDM(ts, &da)); 193 PetscCall(DMGetLocalVector(da, &localX)); 194 PetscCall(DMGetLocalVector(da, &localXdot)); 195 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)); 196 197 hx = 1.0 / (PetscReal)Mx; 198 sx = 1.0 / (hx * hx); 199 200 /* 201 Scatter ghost points to local vector,using the 2-step process 202 DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 203 By placing code between these two statements, computations can be 204 done while messages are in transition. 205 */ 206 PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 207 PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 208 PetscCall(DMGlobalToLocalBegin(da, Xdot, INSERT_VALUES, localXdot)); 209 PetscCall(DMGlobalToLocalEnd(da, Xdot, INSERT_VALUES, localXdot)); 210 211 /* 212 Get pointers to vector data 213 */ 214 PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 215 PetscCall(DMDAVecGetArrayRead(da, localXdot, &xdot)); 216 PetscCall(DMDAVecGetArray(da, F, &f)); 217 218 /* 219 Get local grid boundaries 220 */ 221 PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 222 223 /* 224 Compute function over the locally owned part of the grid 225 */ 226 for (i = xs; i < xs + xm; i++) { 227 f[i].w = x[i].w + ctx->kappa * (x[i - 1].u + x[i + 1].u - 2.0 * x[i].u) * sx; 228 if (ctx->cahnhillard) { 229 switch (ctx->energy) { 230 case 1: /* double well */ f[i].w += -x[i].u * x[i].u * x[i].u + x[i].u; break; 231 case 2: /* double obstacle */ f[i].w += x[i].u; break; 232 case 3: /* logarithmic */ 233 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; 234 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; 235 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; 236 break; 237 case 4: break; 238 } 239 } 240 f[i].u = xdot[i].u - (x[i - 1].w + x[i + 1].w - 2.0 * x[i].w) * sx; 241 if (ctx->energy == 4) { 242 f[i].u = xdot[i].u; 243 /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */ 244 r = (1.0 - x[i + 1].u * x[i + 1].u) * (x[i + 2].w - x[i].w) * .5 / hx; 245 l = (1.0 - x[i - 1].u * x[i - 1].u) * (x[i].w - x[i - 2].w) * .5 / hx; 246 f[i].u -= (r - l) * .5 / hx; 247 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; 248 } 249 } 250 251 /* 252 Restore vectors 253 */ 254 PetscCall(DMDAVecRestoreArrayRead(da, localXdot, &xdot)); 255 PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 256 PetscCall(DMDAVecRestoreArray(da, F, &f)); 257 PetscCall(DMRestoreLocalVector(da, &localX)); 258 PetscCall(DMRestoreLocalVector(da, &localXdot)); 259 PetscFunctionReturn(0); 260 } 261 262 /* ------------------------------------------------------------------- */ 263 PetscErrorCode FormInitialSolution(DM da, Vec X, PetscReal kappa) { 264 PetscInt i, xs, xm, Mx, xgs, xgm; 265 Field *x; 266 PetscReal hx, xx, r, sx; 267 Vec Xg; 268 269 PetscFunctionBegin; 270 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)); 271 272 hx = 1.0 / (PetscReal)Mx; 273 sx = 1.0 / (hx * hx); 274 275 /* 276 Get pointers to vector data 277 */ 278 PetscCall(DMCreateLocalVector(da, &Xg)); 279 PetscCall(DMDAVecGetArray(da, Xg, &x)); 280 281 /* 282 Get local grid boundaries 283 */ 284 PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 285 PetscCall(DMDAGetGhostCorners(da, &xgs, NULL, NULL, &xgm, NULL, NULL)); 286 287 /* 288 Compute u function over the locally owned part of the grid including ghost points 289 */ 290 for (i = xgs; i < xgs + xgm; i++) { 291 xx = i * hx; 292 r = PetscSqrtReal((xx - .5) * (xx - .5)); 293 if (r < .125) x[i].u = 1.0; 294 else x[i].u = -.50; 295 /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 296 x[i].w = 0; 297 } 298 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; 299 300 /* 301 Restore vectors 302 */ 303 PetscCall(DMDAVecRestoreArray(da, Xg, &x)); 304 305 /* Grab only the global part of the vector */ 306 PetscCall(VecSet(X, 0)); 307 PetscCall(DMLocalToGlobalBegin(da, Xg, ADD_VALUES, X)); 308 PetscCall(DMLocalToGlobalEnd(da, Xg, ADD_VALUES, X)); 309 PetscCall(VecDestroy(&Xg)); 310 PetscFunctionReturn(0); 311 } 312 313 /*TEST 314 315 build: 316 requires: !complex !single 317 318 test: 319 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 320 requires: x 321 322 TEST*/ 323