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 */ 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 ierr = PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);CHKERRQ(ierr); 61 ierr = PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds);CHKERRQ(ierr); 62 ierr = PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600);CHKERRQ(ierr); 63 ctx.energy = 1; 64 ierr = PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);CHKERRQ(ierr); 65 ctx.tol = 1.0e-8; 66 ierr = PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);CHKERRQ(ierr); 67 ctx.theta = .001; 68 ctx.theta_c = 1.0; 69 ierr = PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);CHKERRQ(ierr); 70 ierr = PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);CHKERRQ(ierr); 71 72 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73 Create distributed array (DMDA) to manage parallel grid and vectors 74 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 75 ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da);CHKERRQ(ierr); 76 ierr = DMSetFromOptions(da);CHKERRQ(ierr); 77 ierr = DMSetUp(da);CHKERRQ(ierr); 78 ierr = DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");CHKERRQ(ierr); 79 ierr = DMDASetFieldName(da,1,"Biharmonic heat equation: u");CHKERRQ(ierr); 80 ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 81 dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx); 82 83 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84 Extract global vectors from DMDA; then duplicate for remaining 85 vectors that are the same types 86 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 87 ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 88 ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 89 90 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91 Create timestepping solver context 92 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 93 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 94 ierr = TSSetDM(ts,da);CHKERRQ(ierr); 95 ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 96 ierr = TSSetIFunction(ts,NULL,FormFunction,&ctx);CHKERRQ(ierr); 97 ierr = TSSetMaxTime(ts,.02);CHKERRQ(ierr); 98 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 99 100 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101 Create matrix data structure; set Jacobian evaluation routine 102 103 < Set Jacobian matrix data structure and default Jacobian evaluation 104 routine. User can override with: 105 -snes_mf : matrix-free Newton-Krylov method with no preconditioning 106 (unless user explicitly sets preconditioner) 107 -snes_mf_operator : form preconditioning matrix as set by the user, 108 but use matrix-free approx for Jacobian-vector 109 products within Newton-Krylov method 110 111 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 112 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 113 ierr = SNESSetType(snes,SNESVINEWTONRSLS);CHKERRQ(ierr); 114 ierr = DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr); 115 ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr); 116 ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr); 117 ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr); 118 ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);CHKERRQ(ierr); 119 ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr); 120 ierr = MatFDColoringSetUp(J,iscoloring,matfdcoloring);CHKERRQ(ierr); 121 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 122 ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);CHKERRQ(ierr); 123 124 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125 Customize nonlinear solver 126 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127 ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 128 129 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130 Set initial conditions 131 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132 ierr = FormInitialSolution(da,x,ctx.kappa);CHKERRQ(ierr); 133 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 134 ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 135 136 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137 Set runtime options 138 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 140 141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142 Solve nonlinear system 143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144 ierr = TSSolve(ts,x);CHKERRQ(ierr); 145 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 146 147 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148 Free work space. All PETSc objects should be destroyed when they 149 are no longer needed. 150 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 151 ierr = MatDestroy(&J);CHKERRQ(ierr); 152 ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr); 153 ierr = VecDestroy(&x);CHKERRQ(ierr); 154 ierr = VecDestroy(&r);CHKERRQ(ierr); 155 ierr = TSDestroy(&ts);CHKERRQ(ierr); 156 ierr = DMDestroy(&da);CHKERRQ(ierr); 157 158 ierr = PetscFinalize(); 159 return ierr; 160 } 161 162 typedef struct {PetscScalar w,u;} Field; 163 /* ------------------------------------------------------------------- */ 164 /* 165 FormFunction - Evaluates nonlinear function, F(x). 166 167 Input Parameters: 168 . ts - the TS context 169 . X - input vector 170 . ptr - optional user-defined context, as set by SNESSetFunction() 171 172 Output Parameter: 173 . F - function vector 174 */ 175 PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void *ptr) 176 { 177 DM da; 178 PetscErrorCode ierr; 179 PetscInt i,Mx,xs,xm; 180 PetscReal hx,sx; 181 PetscScalar r,l; 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 (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; 232 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; 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 case 4: 236 break; 237 } 238 } 239 f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx; 240 if (ctx->energy==4) { 241 f[i].u = xdot[i].u; 242 /* approximation of \grad (M(u) \grad w), where M(u) = (1-u^2) */ 243 r = (1.0 - x[i+1].u*x[i+1].u)*(x[i+2].w-x[i].w)*.5/hx; 244 l = (1.0 - x[i-1].u*x[i-1].u)*(x[i].w-x[i-2].w)*.5/hx; 245 f[i].u -= (r - l)*.5/hx; 246 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; 247 } 248 } 249 250 /* 251 Restore vectors 252 */ 253 ierr = DMDAVecRestoreArrayRead(da,localXdot,&xdot);CHKERRQ(ierr); 254 ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr); 255 ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 256 ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 257 ierr = DMRestoreLocalVector(da,&localXdot);CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 262 /* ------------------------------------------------------------------- */ 263 PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa) 264 { 265 PetscErrorCode ierr; 266 PetscInt i,xs,xm,Mx,xgs,xgm; 267 Field *x; 268 PetscReal hx,xx,r,sx; 269 Vec Xg; 270 271 PetscFunctionBegin; 272 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); 273 274 hx = 1.0/(PetscReal)Mx; 275 sx = 1.0/(hx*hx); 276 277 /* 278 Get pointers to vector data 279 */ 280 ierr = DMCreateLocalVector(da,&Xg);CHKERRQ(ierr); 281 ierr = DMDAVecGetArray(da,Xg,&x);CHKERRQ(ierr); 282 283 /* 284 Get local grid boundaries 285 */ 286 ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 287 ierr = DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL);CHKERRQ(ierr); 288 289 /* 290 Compute u function over the locally owned part of the grid including ghost points 291 */ 292 for (i=xgs; i<xgs+xgm; i++) { 293 xx = i*hx; 294 r = PetscSqrtReal((xx-.5)*(xx-.5)); 295 if (r < .125) x[i].u = 1.0; 296 else x[i].u = -.50; 297 /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 298 x[i].w = 0; 299 } 300 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; 301 302 /* 303 Restore vectors 304 */ 305 ierr = DMDAVecRestoreArray(da,Xg,&x);CHKERRQ(ierr); 306 307 /* Grab only the global part of the vector */ 308 ierr = VecSet(X,0);CHKERRQ(ierr); 309 ierr = DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X);CHKERRQ(ierr); 310 ierr = DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X);CHKERRQ(ierr); 311 ierr = VecDestroy(&Xg);CHKERRQ(ierr); 312 PetscFunctionReturn(0); 313 } 314 315 /*TEST 316 317 build: 318 requires: !complex !single 319 320 test: 321 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 322 requires: x 323 324 TEST*/ 325