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