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 #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 ierr = PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);CHKERRQ(ierr); 58 ctx.cahnhillard = PETSC_FALSE; 59 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 ierr = PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);CHKERRQ(ierr); 66 ctx.tol = 1.0e-8; 67 ierr = PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);CHKERRQ(ierr); 68 ctx.theta = .001; 69 ctx.theta_c = 1.0; 70 ierr = PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);CHKERRQ(ierr); 71 ierr = PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);CHKERRQ(ierr); 72 73 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74 Create distributed array (DMDA) to manage parallel grid and vectors 75 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 76 ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da);CHKERRQ(ierr); 77 ierr = DMSetFromOptions(da);CHKERRQ(ierr); 78 ierr = DMSetUp(da);CHKERRQ(ierr); 79 ierr = DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx");CHKERRQ(ierr); 80 ierr = DMDASetFieldName(da,1,"Biharmonic heat equation: u");CHKERRQ(ierr); 81 ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 82 dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx); 83 84 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 85 Extract global vectors from DMDA; then duplicate for remaining 86 vectors that are the same types 87 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 88 ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 89 ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 90 91 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92 Create timestepping solver context 93 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 95 ierr = TSSetDM(ts,da);CHKERRQ(ierr); 96 ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 97 ierr = TSSetIFunction(ts,NULL,FormFunction,&ctx);CHKERRQ(ierr); 98 ierr = TSSetMaxTime(ts,.02);CHKERRQ(ierr); 99 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);CHKERRQ(ierr); 100 101 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102 Create matrix data structure; set Jacobian evaluation routine 103 104 < Set Jacobian matrix data structure and default Jacobian evaluation 105 routine. User can override with: 106 -snes_mf : matrix-free Newton-Krylov method with no preconditioning 107 (unless user explicitly sets preconditioner) 108 -snes_mf_operator : form preconditioning matrix as set by the user, 109 but use matrix-free approx for Jacobian-vector 110 products within Newton-Krylov method 111 112 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 113 ierr = TSGetSNES(ts,&snes);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 Field *x,*xdot,*f; 182 Vec localX,localXdot; 183 UserCtx *ctx = (UserCtx*)ptr; 184 185 PetscFunctionBegin; 186 ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 187 ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 188 ierr = DMGetLocalVector(da,&localXdot);CHKERRQ(ierr); 189 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); 190 191 hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 192 193 /* 194 Scatter ghost points to local vector,using the 2-step process 195 DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 196 By placing code between these two statements, computations can be 197 done while messages are in transition. 198 */ 199 ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 200 ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 201 ierr = DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 202 ierr = DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 203 204 /* 205 Get pointers to vector data 206 */ 207 ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr); 208 ierr = DMDAVecGetArrayRead(da,localXdot,&xdot);CHKERRQ(ierr); 209 ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 210 211 /* 212 Get local grid boundaries 213 */ 214 ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 215 216 /* 217 Compute function over the locally owned part of the grid 218 */ 219 for (i=xs; i<xs+xm; i++) { 220 f[i].w = x[i].w + ctx->kappa*(x[i-1].u + x[i+1].u - 2.0*x[i].u)*sx; 221 if (ctx->cahnhillard) { 222 switch (ctx->energy) { 223 case 1: /* double well */ 224 f[i].w += -x[i].u*x[i].u*x[i].u + x[i].u; 225 break; 226 case 2: /* double obstacle */ 227 f[i].w += x[i].u; 228 break; 229 case 3: /* logarithmic */ 230 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; 231 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; 232 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; 233 break; 234 } 235 } 236 f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx; 237 } 238 239 /* 240 Restore vectors 241 */ 242 ierr = DMDAVecRestoreArrayRead(da,localXdot,&xdot);CHKERRQ(ierr); 243 ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr); 244 ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 245 ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 246 ierr = DMRestoreLocalVector(da,&localXdot);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 /* ------------------------------------------------------------------- */ 251 PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa) 252 { 253 PetscErrorCode ierr; 254 PetscInt i,xs,xm,Mx,xgs,xgm; 255 Field *x; 256 PetscReal hx,xx,r,sx; 257 Vec Xg; 258 259 PetscFunctionBegin; 260 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); 261 262 hx = 1.0/(PetscReal)Mx; 263 sx = 1.0/(hx*hx); 264 265 /* 266 Get pointers to vector data 267 */ 268 ierr = DMCreateLocalVector(da,&Xg);CHKERRQ(ierr); 269 ierr = DMDAVecGetArray(da,Xg,&x);CHKERRQ(ierr); 270 271 /* 272 Get local grid boundaries 273 */ 274 ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 275 ierr = DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL);CHKERRQ(ierr); 276 277 /* 278 Compute u function over the locally owned part of the grid including ghost points 279 */ 280 for (i=xgs; i<xgs+xgm; i++) { 281 xx = i*hx; 282 r = PetscSqrtReal((xx-.5)*(xx-.5)); 283 if (r < .125) x[i].u = 1.0; 284 else x[i].u = -.50; 285 /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 286 x[i].w = 0; 287 } 288 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; 289 290 /* 291 Restore vectors 292 */ 293 ierr = DMDAVecRestoreArray(da,Xg,&x);CHKERRQ(ierr); 294 295 /* Grab only the global part of the vector */ 296 ierr = VecSet(X,0);CHKERRQ(ierr); 297 ierr = DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X);CHKERRQ(ierr); 298 ierr = DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X);CHKERRQ(ierr); 299 ierr = VecDestroy(&Xg);CHKERRQ(ierr); 300 PetscFunctionReturn(0); 301 } 302 303 /*TEST 304 305 build: 306 requires: !complex !single 307 308 test: 309 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 310 requires: x 311 312 TEST*/ 313