Actual source code: biharmonic2.c

petsc-3.12.0 2019-09-29
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  -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: (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


 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;

 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;
 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;

 53:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 54:      Initialize program
 55:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 56:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 57:   ctx.kappa = 1.0;
 58:   PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
 59:   ctx.cahnhillard = PETSC_FALSE;

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

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

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

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

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:      Create matrix data structure; set Jacobian evaluation routine

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

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

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126:      Customize nonlinear solver
127:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128:   TSSetType(ts,TSBEULER);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:      Set initial conditions
132:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133:   FormInitialSolution(da,x,ctx.kappa);
134:   TSSetTimeStep(ts,dt);
135:   TSSetSolution(ts,x);

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Set runtime options
139:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   TSSetFromOptions(ts);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Solve nonlinear system
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   TSSolve(ts,x);
146:   TSGetStepNumber(ts,&steps);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:      Free work space.  All PETSc objects should be destroyed when they
150:      are no longer needed.
151:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152:   MatDestroy(&J);
153:   MatFDColoringDestroy(&matfdcoloring);
154:   VecDestroy(&x);
155:   VecDestroy(&r);
156:   TSDestroy(&ts);
157:   DMDestroy(&da);

159:   PetscFinalize();
160:   return ierr;
161: }

163: typedef struct {PetscScalar w,u;} Field;
164: /* ------------------------------------------------------------------- */
165: /*
166:    FormFunction - Evaluates nonlinear function, F(x).

168:    Input Parameters:
169: .  ts - the TS context
170: .  X - input vector
171: .  ptr - optional user-defined context, as set by SNESSetFunction()

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;
180:   PetscInt       i,Mx,xs,xm;
181:   PetscReal      hx,sx;
182:   Field          *x,*xdot,*f;
183:   Vec            localX,localXdot;
184:   UserCtx        *ctx = (UserCtx*)ptr;

187:   TSGetDM(ts,&da);
188:   DMGetLocalVector(da,&localX);
189:   DMGetLocalVector(da,&localXdot);
190:   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);

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

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:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
201:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
202:   DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);
203:   DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);

205:   /*
206:      Get pointers to vector data
207:   */
208:   DMDAVecGetArrayRead(da,localX,&x);
209:   DMDAVecGetArrayRead(da,localXdot,&xdot);
210:   DMDAVecGetArray(da,F,&f);

212:   /*
213:      Get local grid boundaries
214:   */
215:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

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:   }

240:   /*
241:      Restore vectors
242:   */
243:   DMDAVecRestoreArrayRead(da,localXdot,&xdot);
244:   DMDAVecRestoreArrayRead(da,localX,&x);
245:   DMDAVecRestoreArray(da,F,&f);
246:   DMRestoreLocalVector(da,&localX);
247:   DMRestoreLocalVector(da,&localXdot);
248:   return(0);
249: }

251: /* ------------------------------------------------------------------- */
252: PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa)
253: {
255:   PetscInt       i,xs,xm,Mx;
256:   Field          *x;
257:   PetscReal      hx,xx,r,sx;

260:   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:   hx = 1.0/(PetscReal)Mx;
263:   sx = 1.0/(hx*hx);

265:   /*
266:      Get pointers to vector data
267:   */
268:   DMDAVecGetArray(da,X,&x);

270:   /*
271:      Get local grid boundaries
272:   */
273:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

275:   /*
276:      Compute function over the locally owned part of the grid
277:   */
278:   for (i=xs; i<xs+xm; i++) {
279:     xx = i*hx;
280:     r = PetscSqrtReal((xx-.5)*(xx-.5));
281:     if (r < .125) x[i].u = 1.0;
282:     else          x[i].u = -.50;
283:     /*  u[i] = PetscPowScalar(x - .5,4.0); */
284:   }
285:   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;

287:   /*
288:      Restore vectors
289:   */
290:   DMDAVecRestoreArray(da,X,&x);
291:   return(0);
292: }

294: /*TEST

296:    build:
297:      requires: !complex !single

299:    test:
300:      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
301:      requires: x

303: TEST*/