xref: /petsc/src/ts/tutorials/phasefield/biharmonic2.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
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 ./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 
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: (this fails after a few timesteps 12/17/2017)
22 ---------------
23 ./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 */
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 
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));
78 
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);
89 
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));
96 
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));
106 
107   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108      Create matrix data structure; set Jacobian evaluation routine
109 
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
117 
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));
129 
130   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131      Customize nonlinear solver
132    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133   PetscCall(TSSetType(ts, TSBEULER));
134 
135   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136      Set initial conditions
137    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138   PetscCall(FormInitialSolution(da, x, ctx.kappa));
139   PetscCall(TSSetTimeStep(ts, dt));
140   PetscCall(TSSetSolution(ts, x));
141 
142   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143      Set runtime options
144    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145   PetscCall(TSSetFromOptions(ts));
146 
147   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148      Solve nonlinear system
149      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150   PetscCall(TSSolve(ts, x));
151   PetscCall(TSGetStepNumber(ts, &steps));
152 
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));
163 
164   PetscCall(PetscFinalize());
165   return 0;
166 }
167 
168 typedef struct {
169   PetscScalar w, u;
170 } Field;
171 /* ------------------------------------------------------------------- */
172 /*
173    FormFunction - Evaluates nonlinear function, F(x).
174 
175    Input Parameters:
176 .  ts - the TS context
177 .  X - input vector
178 .  ptr - optional user-defined context, as set by SNESSetFunction()
179 
180    Output Parameter:
181 .  F - function vector
182  */
183 PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec Xdot, Vec F, void *ptr) {
184   DM        da;
185   PetscInt  i, Mx, xs, xm;
186   PetscReal hx, sx;
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 (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;
234         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;
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       }
238     }
239     f[i].u = xdot[i].u - (x[i - 1].w + x[i + 1].w - 2.0 * x[i].w) * sx;
240   }
241 
242   /*
243      Restore vectors
244   */
245   PetscCall(DMDAVecRestoreArrayRead(da, localXdot, &xdot));
246   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
247   PetscCall(DMDAVecRestoreArray(da, F, &f));
248   PetscCall(DMRestoreLocalVector(da, &localX));
249   PetscCall(DMRestoreLocalVector(da, &localXdot));
250   PetscFunctionReturn(0);
251 }
252 
253 /* ------------------------------------------------------------------- */
254 PetscErrorCode FormInitialSolution(DM da, Vec X, PetscReal kappa) {
255   PetscInt  i, xs, xm, Mx, xgs, xgm;
256   Field    *x;
257   PetscReal hx, xx, r, sx;
258   Vec       Xg;
259 
260   PetscFunctionBegin;
261   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));
262 
263   hx = 1.0 / (PetscReal)Mx;
264   sx = 1.0 / (hx * hx);
265 
266   /*
267      Get pointers to vector data
268   */
269   PetscCall(DMCreateLocalVector(da, &Xg));
270   PetscCall(DMDAVecGetArray(da, Xg, &x));
271 
272   /*
273      Get local grid boundaries
274   */
275   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
276   PetscCall(DMDAGetGhostCorners(da, &xgs, NULL, NULL, &xgm, NULL, NULL));
277 
278   /*
279      Compute u function over the locally owned part of the grid including ghost points
280   */
281   for (i = xgs; i < xgs + xgm; i++) {
282     xx = i * hx;
283     r  = PetscSqrtReal((xx - .5) * (xx - .5));
284     if (r < .125) x[i].u = 1.0;
285     else x[i].u = -.50;
286     /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
287     x[i].w = 0;
288   }
289   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;
290 
291   /*
292      Restore vectors
293   */
294   PetscCall(DMDAVecRestoreArray(da, Xg, &x));
295 
296   /* Grab only the global part of the vector */
297   PetscCall(VecSet(X, 0));
298   PetscCall(DMLocalToGlobalBegin(da, Xg, ADD_VALUES, X));
299   PetscCall(DMLocalToGlobalEnd(da, Xg, ADD_VALUES, X));
300   PetscCall(VecDestroy(&Xg));
301   PetscFunctionReturn(0);
302 }
303 
304 /*TEST
305 
306    build:
307      requires: !complex !single
308 
309    test:
310      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
311      requires: x
312 
313 TEST*/
314