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 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL)); 58 ctx.cahnhillard = PETSC_FALSE; 59 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL)); 60 CHKERRQ(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds)); 61 CHKERRQ(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600)); 62 ctx.energy = 1; 63 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL)); 64 ctx.tol = 1.0e-8; 65 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL)); 66 ctx.theta = .001; 67 ctx.theta_c = 1.0; 68 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL)); 69 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL)); 70 71 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72 Create distributed array (DMDA) to manage parallel grid and vectors 73 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74 CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da)); 75 CHKERRQ(DMSetFromOptions(da)); 76 CHKERRQ(DMSetUp(da)); 77 CHKERRQ(DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx")); 78 CHKERRQ(DMDASetFieldName(da,1,"Biharmonic heat equation: u")); 79 CHKERRQ(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 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 CHKERRQ(DMCreateGlobalVector(da,&x)); 87 CHKERRQ(VecDuplicate(x,&r)); 88 89 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 90 Create timestepping solver context 91 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 92 CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 93 CHKERRQ(TSSetDM(ts,da)); 94 CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 95 CHKERRQ(TSSetIFunction(ts,NULL,FormFunction,&ctx)); 96 CHKERRQ(TSSetMaxTime(ts,.02)); 97 CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 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 CHKERRQ(TSGetSNES(ts,&snes)); 112 CHKERRQ(SNESSetType(snes,SNESVINEWTONRSLS)); 113 CHKERRQ(DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring)); 114 CHKERRQ(DMSetMatType(da,MATAIJ)); 115 CHKERRQ(DMCreateMatrix(da,&J)); 116 CHKERRQ(MatFDColoringCreate(J,iscoloring,&matfdcoloring)); 117 CHKERRQ(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts)); 118 CHKERRQ(MatFDColoringSetFromOptions(matfdcoloring)); 119 CHKERRQ(MatFDColoringSetUp(J,iscoloring,matfdcoloring)); 120 CHKERRQ(ISColoringDestroy(&iscoloring)); 121 CHKERRQ(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring)); 122 123 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124 Customize nonlinear solver 125 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126 CHKERRQ(TSSetType(ts,TSBEULER)); 127 128 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129 Set initial conditions 130 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131 CHKERRQ(FormInitialSolution(da,x,ctx.kappa)); 132 CHKERRQ(TSSetTimeStep(ts,dt)); 133 CHKERRQ(TSSetSolution(ts,x)); 134 135 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136 Set runtime options 137 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138 CHKERRQ(TSSetFromOptions(ts)); 139 140 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 141 Solve nonlinear system 142 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 143 CHKERRQ(TSSolve(ts,x)); 144 CHKERRQ(TSGetStepNumber(ts,&steps)); 145 146 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147 Free work space. All PETSc objects should be destroyed when they 148 are no longer needed. 149 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 150 CHKERRQ(MatDestroy(&J)); 151 CHKERRQ(MatFDColoringDestroy(&matfdcoloring)); 152 CHKERRQ(VecDestroy(&x)); 153 CHKERRQ(VecDestroy(&r)); 154 CHKERRQ(TSDestroy(&ts)); 155 CHKERRQ(DMDestroy(&da)); 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 PetscInt i,Mx,xs,xm; 178 PetscReal hx,sx; 179 PetscScalar r,l; 180 Field *x,*xdot,*f; 181 Vec localX,localXdot; 182 UserCtx *ctx = (UserCtx*)ptr; 183 184 PetscFunctionBegin; 185 CHKERRQ(TSGetDM(ts,&da)); 186 CHKERRQ(DMGetLocalVector(da,&localX)); 187 CHKERRQ(DMGetLocalVector(da,&localXdot)); 188 CHKERRQ(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)); 189 190 hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 191 192 /* 193 Scatter ghost points to local vector,using the 2-step process 194 DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 195 By placing code between these two statements, computations can be 196 done while messages are in transition. 197 */ 198 CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 199 CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 200 CHKERRQ(DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot)); 201 CHKERRQ(DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot)); 202 203 /* 204 Get pointers to vector data 205 */ 206 CHKERRQ(DMDAVecGetArrayRead(da,localX,&x)); 207 CHKERRQ(DMDAVecGetArrayRead(da,localXdot,&xdot)); 208 CHKERRQ(DMDAVecGetArray(da,F,&f)); 209 210 /* 211 Get local grid boundaries 212 */ 213 CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 214 215 /* 216 Compute function over the locally owned part of the grid 217 */ 218 for (i=xs; i<xs+xm; i++) { 219 f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx; 220 if (ctx->cahnhillard) { 221 switch (ctx->energy) { 222 case 1: /* double well */ 223 f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u; 224 break; 225 case 2: /* double obstacle */ 226 f[i].w += x[i].u; 227 break; 228 case 3: /* logarithmic */ 229 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; 230 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; 231 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; 232 break; 233 case 4: 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 if (ctx->energy==4) { 239 f[i].u = xdot[i].u; 240 /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */ 241 r = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx; 242 l = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx; 243 f[i].u -= (r - l)*.5/hx; 244 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; 245 } 246 } 247 248 /* 249 Restore vectors 250 */ 251 CHKERRQ(DMDAVecRestoreArrayRead(da,localXdot,&xdot)); 252 CHKERRQ(DMDAVecRestoreArrayRead(da,localX,&x)); 253 CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 254 CHKERRQ(DMRestoreLocalVector(da,&localX)); 255 CHKERRQ(DMRestoreLocalVector(da,&localXdot)); 256 PetscFunctionReturn(0); 257 } 258 259 /* ------------------------------------------------------------------- */ 260 PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa) 261 { 262 PetscInt i,xs,xm,Mx,xgs,xgm; 263 Field *x; 264 PetscReal hx,xx,r,sx; 265 Vec Xg; 266 267 PetscFunctionBegin; 268 CHKERRQ(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)); 269 270 hx = 1.0/(PetscReal)Mx; 271 sx = 1.0/(hx*hx); 272 273 /* 274 Get pointers to vector data 275 */ 276 CHKERRQ(DMCreateLocalVector(da,&Xg)); 277 CHKERRQ(DMDAVecGetArray(da,Xg,&x)); 278 279 /* 280 Get local grid boundaries 281 */ 282 CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 283 CHKERRQ(DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL)); 284 285 /* 286 Compute u function over the locally owned part of the grid including ghost points 287 */ 288 for (i=xgs; i<xgs+xgm; i++) { 289 xx = i*hx; 290 r = PetscSqrtReal((xx-.5)*(xx-.5)); 291 if (r < .125) x[i].u = 1.0; 292 else x[i].u = -.50; 293 /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 294 x[i].w = 0; 295 } 296 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; 297 298 /* 299 Restore vectors 300 */ 301 CHKERRQ(DMDAVecRestoreArray(da,Xg,&x)); 302 303 /* Grab only the global part of the vector */ 304 CHKERRQ(VecSet(X,0)); 305 CHKERRQ(DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X)); 306 CHKERRQ(DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X)); 307 CHKERRQ(VecDestroy(&Xg)); 308 PetscFunctionReturn(0); 309 } 310 311 /*TEST 312 313 build: 314 requires: !complex !single 315 316 test: 317 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 318 requires: x 319 320 TEST*/ 321