xref: /petsc/src/ts/tutorials/phasefield/biharmonic3.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
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 {PetscBool cahnhillard;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta;PetscReal theta_c;} UserCtx;
36 
37 int main(int argc,char **argv)
38 {
39   TS             ts;                           /* nonlinear solver */
40   Vec            x,r;                          /* solution, residual vectors */
41   Mat            J;                            /* Jacobian matrix */
42   PetscInt       steps,Mx;
43   DM             da;
44   MatFDColoring  matfdcoloring;
45   ISColoring     iscoloring;
46   PetscReal      dt;
47   PetscReal      vbounds[] = {-100000,100000,-1.1,1.1};
48   SNES           snes;
49   UserCtx        ctx;
50 
51   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52      Initialize program
53      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
55   ctx.kappa       = 1.0;
56   PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL));
57   ctx.cahnhillard = PETSC_FALSE;
58   PetscCall(PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL));
59   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds));
60   PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600));
61   ctx.energy      = 1;
62   PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));
63   ctx.tol     = 1.0e-8;
64   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL));
65   ctx.theta   = .001;
66   ctx.theta_c = 1.0;
67   PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL));
68   PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL));
69 
70   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71      Create distributed array (DMDA) to manage parallel grid and vectors
72   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da));
74   PetscCall(DMSetFromOptions(da));
75   PetscCall(DMSetUp(da));
76   PetscCall(DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx"));
77   PetscCall(DMDASetFieldName(da,1,"Biharmonic heat equation: u"));
78   PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0));
79   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
80 
81   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82      Extract global vectors from DMDA; then duplicate for remaining
83      vectors that are the same types
84    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85   PetscCall(DMCreateGlobalVector(da,&x));
86   PetscCall(VecDuplicate(x,&r));
87 
88   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89      Create timestepping solver context
90      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
92   PetscCall(TSSetDM(ts,da));
93   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
94   PetscCall(TSSetIFunction(ts,NULL,FormFunction,&ctx));
95   PetscCall(TSSetMaxTime(ts,.02));
96   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
97 
98   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99      Create matrix data structure; set Jacobian evaluation routine
100 
101 <     Set Jacobian matrix data structure and default Jacobian evaluation
102      routine. User can override with:
103      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
104                 (unless user explicitly sets preconditioner)
105      -snes_mf_operator : form preconditioning matrix as set by the user,
106                          but use matrix-free approx for Jacobian-vector
107                          products within Newton-Krylov method
108 
109      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110   PetscCall(TSGetSNES(ts,&snes));
111   PetscCall(SNESSetType(snes,SNESVINEWTONRSLS));
112   PetscCall(DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring));
113   PetscCall(DMSetMatType(da,MATAIJ));
114   PetscCall(DMCreateMatrix(da,&J));
115   PetscCall(MatFDColoringCreate(J,iscoloring,&matfdcoloring));
116   PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts));
117   PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
118   PetscCall(MatFDColoringSetUp(J,iscoloring,matfdcoloring));
119   PetscCall(ISColoringDestroy(&iscoloring));
120   PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring));
121 
122   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123      Customize nonlinear solver
124    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125   PetscCall(TSSetType(ts,TSBEULER));
126 
127   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128      Set initial conditions
129    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130   PetscCall(FormInitialSolution(da,x,ctx.kappa));
131   PetscCall(TSSetTimeStep(ts,dt));
132   PetscCall(TSSetSolution(ts,x));
133 
134   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135      Set runtime options
136    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137   PetscCall(TSSetFromOptions(ts));
138 
139   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140      Solve nonlinear system
141      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142   PetscCall(TSSolve(ts,x));
143   PetscCall(TSGetStepNumber(ts,&steps));
144 
145   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146      Free work space.  All PETSc objects should be destroyed when they
147      are no longer needed.
148    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149   PetscCall(MatDestroy(&J));
150   PetscCall(MatFDColoringDestroy(&matfdcoloring));
151   PetscCall(VecDestroy(&x));
152   PetscCall(VecDestroy(&r));
153   PetscCall(TSDestroy(&ts));
154   PetscCall(DMDestroy(&da));
155 
156   PetscCall(PetscFinalize());
157   return 0;
158 }
159 
160 typedef struct {PetscScalar w,u;} Field;
161 /* ------------------------------------------------------------------- */
162 /*
163    FormFunction - Evaluates nonlinear function, F(x).
164 
165    Input Parameters:
166 .  ts - the TS context
167 .  X - input vector
168 .  ptr - optional user-defined context, as set by SNESSetFunction()
169 
170    Output Parameter:
171 .  F - function vector
172  */
173 PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
174 {
175   DM             da;
176   PetscInt       i,Mx,xs,xm;
177   PetscReal      hx,sx;
178   PetscScalar    r,l;
179   Field          *x,*xdot,*f;
180   Vec            localX,localXdot;
181   UserCtx        *ctx = (UserCtx*)ptr;
182 
183   PetscFunctionBegin;
184   PetscCall(TSGetDM(ts,&da));
185   PetscCall(DMGetLocalVector(da,&localX));
186   PetscCall(DMGetLocalVector(da,&localXdot));
187   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));
188 
189   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
190 
191   /*
192      Scatter ghost points to local vector,using the 2-step process
193         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
194      By placing code between these two statements, computations can be
195      done while messages are in transition.
196   */
197   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
198   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
199   PetscCall(DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot));
200   PetscCall(DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot));
201 
202   /*
203      Get pointers to vector data
204   */
205   PetscCall(DMDAVecGetArrayRead(da,localX,&x));
206   PetscCall(DMDAVecGetArrayRead(da,localXdot,&xdot));
207   PetscCall(DMDAVecGetArray(da,F,&f));
208 
209   /*
210      Get local grid boundaries
211   */
212   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
213 
214   /*
215      Compute function over the locally owned part of the grid
216   */
217   for (i=xs; i<xs+xm; i++) {
218     f[i].w =  x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
219     if (ctx->cahnhillard) {
220       switch (ctx->energy) {
221       case 1: /* double well */
222         f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
223         break;
224       case 2: /* double obstacle */
225         f[i].w += x[i].u;
226         break;
227       case 3: /* logarithmic */
228         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;
229         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;
230         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;
231         break;
232       case 4:
233         break;
234       }
235     }
236     f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
237     if (ctx->energy==4) {
238       f[i].u = xdot[i].u;
239       /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */
240       r       = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx;
241       l       = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx;
242       f[i].u -= (r - l)*.5/hx;
243       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;
244     }
245   }
246 
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(0);
256 }
257 
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;
265 
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));
268 
269   hx = 1.0/(PetscReal)Mx;
270   sx = 1.0/(hx*hx);
271 
272   /*
273      Get pointers to vector data
274   */
275   PetscCall(DMCreateLocalVector(da,&Xg));
276   PetscCall(DMDAVecGetArray(da,Xg,&x));
277 
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));
283 
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;
296 
297   /*
298      Restore vectors
299   */
300   PetscCall(DMDAVecRestoreArray(da,Xg,&x));
301 
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(0);
308 }
309 
310 /*TEST
311 
312    build:
313      requires: !complex !single
314 
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
318 
319 TEST*/
320