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