Actual source code: biharmonic2.c

petsc-3.7.1 2016-05-15
Report Typos and Errors
  2: static char help[] = "Solves biharmonic equation in 1d.\n";

  4: /*
  5:   Solves the equation biharmonic equation in split form

  7:     w = -kappa \Delta u
  8:     u_t =  \Delta w
  9:     -1  <= u <= 1
 10:     Periodic boundary conditions

 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 --wait   -ts_type beuler  -da_refine 5 -draw_fields 1 -ts_dt 9.53674e-9

 16:     w = -kappa \Delta u  + u^3  - u
 17:     u_t =  \Delta w
 18:     -1  <= u <= 1
 19:     Periodic boundary conditions

 21: Evolve the Cahn-Hillard equations:
 22: ---------------
 23: ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason  --wait   -ts_type beuler    -da_refine 6 -vi  -draw_fields 1  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard


 26: */
 27: #include <petscdm.h>
 28: #include <petscdmda.h>
 29: #include <petscts.h>
 30: #include <petscdraw.h>

 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;

 40: int main(int argc,char **argv)
 41: {
 42:   TS             ts;                           /* nonlinear solver */
 43:   Vec            x,r;                          /* solution, residual vectors */
 44:   Mat            J;                            /* Jacobian matrix */
 45:   PetscInt       steps,Mx,maxsteps = 10000000;
 47:   DM             da;
 48:   MatFDColoring  matfdcoloring;
 49:   ISColoring     iscoloring;
 50:   PetscReal      dt;
 51:   PetscReal      vbounds[] = {-100000,100000,-1.1,1.1};
 52:   PetscBool      wait;
 53:   Vec            ul,uh;
 54:   SNES           snes;
 55:   UserCtx        ctx;

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:      Initialize program
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   PetscInitialize(&argc,&argv,(char*)0,help);
 61:   ctx.kappa = 1.0;
 62:   PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
 63:   ctx.cahnhillard = PETSC_FALSE;

 65:   PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
 66:   PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);
 67:   PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);
 68:   ctx.energy = 1;
 69:   /*PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);*/
 70:   PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
 71:   ctx.tol     = 1.0e-8;
 72:   PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
 73:   ctx.theta   = .001;
 74:   ctx.theta_c = 1.0;
 75:   PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
 76:   PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:      Create distributed array (DMDA) to manage parallel grid and vectors
 80:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 81:   DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,2,2,NULL,&da);
 82:   DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");
 83:   DMDASetFieldName(da,1,"Biharmonic heat equation: u");
 84:   DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
 85:   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);

 87:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:      Extract global vectors from DMDA; then duplicate for remaining
 89:      vectors that are the same types
 90:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 91:   DMCreateGlobalVector(da,&x);
 92:   VecDuplicate(x,&r);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Create timestepping solver context
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   TSCreate(PETSC_COMM_WORLD,&ts);
 98:   TSSetDM(ts,da);
 99:   TSSetProblemType(ts,TS_NONLINEAR);
100:   TSSetIFunction(ts,NULL,FormFunction,&ctx);
101:   TSSetDuration(ts,maxsteps,.02);
102:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Create matrix data structure; set Jacobian evaluation routine

107: <     Set Jacobian matrix data structure and default Jacobian evaluation
108:      routine. User can override with:
109:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
110:                 (unless user explicitly sets preconditioner)
111:      -snes_mf_operator : form preconditioning matrix as set by the user,
112:                          but use matrix-free approx for Jacobian-vector
113:                          products within Newton-Krylov method

115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   TSGetSNES(ts,&snes);
117:   DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
118:   DMSetMatType(da,MATAIJ);
119:   DMCreateMatrix(da,&J);
120:   MatFDColoringCreate(J,iscoloring,&matfdcoloring);
121:   ISColoringDestroy(&iscoloring);
122:   MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
123:   MatFDColoringSetFromOptions(matfdcoloring);
124:   MatFDColoringSetUp(J,iscoloring,matfdcoloring);
125:   SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);

127:   {
128:     VecDuplicate(x,&ul);
129:     VecDuplicate(x,&uh);
130:     VecStrideSet(ul,0,PETSC_NINFINITY);
131:     VecStrideSet(ul,1,-1.0);
132:     VecStrideSet(uh,0,PETSC_INFINITY);
133:     VecStrideSet(uh,1,1.0);
134:     TSVISetVariableBounds(ts,ul,uh);
135:   }

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Customize nonlinear solver
139:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   TSSetType(ts,TSBEULER);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Set initial conditions
144:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   FormInitialSolution(da,x,ctx.kappa);
146:   TSSetInitialTimeStep(ts,0.0,dt);
147:   TSSetSolution(ts,x);

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:      Set runtime options
151:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152:   TSSetFromOptions(ts);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:      Solve nonlinear system
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   TSSolve(ts,x);
158:   wait = PETSC_FALSE;
159:   PetscOptionsGetBool(NULL,NULL,"-wait",&wait,NULL);
160:   if (wait) {
161:     PetscSleep(-1);
162:   }
163:   TSGetTimeStepNumber(ts,&steps);

165:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:      Free work space.  All PETSc objects should be destroyed when they
167:      are no longer needed.
168:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   {
170:     VecDestroy(&ul);
171:     VecDestroy(&uh);
172:   }
173:   MatDestroy(&J);
174:   MatFDColoringDestroy(&matfdcoloring);
175:   VecDestroy(&x);
176:   VecDestroy(&r);
177:   TSDestroy(&ts);
178:   DMDestroy(&da);

180:   PetscFinalize();
181:   return(0);
182: }

184: typedef struct {PetscScalar w,u;} Field;
185: /* ------------------------------------------------------------------- */
188: /*
189:    FormFunction - Evaluates nonlinear function, F(x).

191:    Input Parameters:
192: .  ts - the TS context
193: .  X - input vector
194: .  ptr - optional user-defined context, as set by SNESSetFunction()

196:    Output Parameter:
197: .  F - function vector
198:  */
199: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr)
200: {
201:   DM             da;
203:   PetscInt       i,Mx,xs,xm;
204:   PetscReal      hx,sx;
205:   Field          *x,*xdot,*f;
206:   Vec            localX,localXdot;
207:   UserCtx        *ctx = (UserCtx*)ptr;

210:   TSGetDM(ts,&da);
211:   DMGetLocalVector(da,&localX);
212:   DMGetLocalVector(da,&localXdot);
213:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
214:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

216:   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);

218:   /*
219:      Scatter ghost points to local vector,using the 2-step process
220:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
221:      By placing code between these two statements, computations can be
222:      done while messages are in transition.
223:   */
224:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
225:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
226:   DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
227:   DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);

229:   /*
230:      Get pointers to vector data
231:   */
232:   DMDAVecGetArrayRead(da,localX,&x);
233:   DMDAVecGetArrayRead(da,localXdot,&xdot);
234:   DMDAVecGetArray(da,F,&f);

236:   /*
237:      Get local grid boundaries
238:   */
239:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

241:   /*
242:      Compute function over the locally owned part of the grid
243:   */
244:   for (i=xs; i<xs+xm; i++) {
245:     f[i].w =  x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx;
246:     if (ctx->cahnhillard) {
247:       switch (ctx->energy) {
248:       case 1: /* double well */
249:         f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u;
250:         break;
251:       case 2: /* double obstacle */
252:         f[i].w += x[i].u;
253:         break;
254:       case 3: /* logarithmic */
255:         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;
256:         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;
257:         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;
258:         break;
259:       }
260:     }
261:     f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx;
262:   }

264:   /*
265:      Restore vectors
266:   */
267:   DMDAVecRestoreArrayRead(da,localXdot,&xdot);
268:   DMDAVecRestoreArrayRead(da,localX,&x);
269:   DMDAVecRestoreArray(da,F,&f);
270:   DMRestoreLocalVector(da,&localX);
271:   DMRestoreLocalVector(da,&localXdot);
272:   return(0);
273: }

275: /* ------------------------------------------------------------------- */
278: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
279: {
281:   PetscInt       i,xs,xm,Mx;
282:   Field          *x;
283:   PetscReal      hx,xx,r,sx;

286:   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);

288:   hx = 1.0/(PetscReal)Mx;
289:   sx = 1.0/(hx*hx);

291:   /*
292:      Get pointers to vector data
293:   */
294:   DMDAVecGetArray(da,X,&x);

296:   /*
297:      Get local grid boundaries
298:   */
299:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

301:   /*
302:      Compute function over the locally owned part of the grid
303:   */
304:   for (i=xs; i<xs+xm; i++) {
305:     xx = i*hx;
306:     r = PetscSqrtReal((xx-.5)*(xx-.5));
307:     if (r < .125) x[i].u = 1.0;
308:     else          x[i].u = -.50;
309:     /*  u[i] = PetscPowScalar(x - .5,4.0); */
310:   }
311:   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;

313:   /*
314:      Restore vectors
315:   */
316:   DMDAVecRestoreArray(da,X,&x);
317:   return(0);
318: }