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 DM da; 44 MatFDColoring matfdcoloring; 45 ISColoring iscoloring; 46 PetscReal dt; 47 PetscReal vbounds[] = {-100000,100000,-1.1,1.1}; 48 SNES snes; 49 UserCtx ctx; 50 51 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 52 Initialize program 53 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 54 PetscFunctionBeginUser; 55 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 56 ctx.kappa = 1.0; 57 PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL)); 58 ctx.cahnhillard = PETSC_FALSE; 59 60 PetscCall(PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL)); 61 PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),2,vbounds)); 62 PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),600,600)); 63 ctx.energy = 1; 64 /*PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));*/ 65 PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL)); 66 ctx.tol = 1.0e-8; 67 PetscCall(PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL)); 68 ctx.theta = .001; 69 ctx.theta_c = 1.0; 70 PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL)); 71 PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL)); 72 73 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74 Create distributed array (DMDA) to manage parallel grid and vectors 75 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 76 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,2,2,NULL,&da)); 77 PetscCall(DMSetFromOptions(da)); 78 PetscCall(DMSetUp(da)); 79 PetscCall(DMDASetFieldName(da,0,"Biharmonic heat equation: w = -kappa*u_xx")); 80 PetscCall(DMDASetFieldName(da,1,"Biharmonic heat equation: u")); 81 PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 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 PetscCall(DMCreateGlobalVector(da,&x)); 89 PetscCall(VecDuplicate(x,&r)); 90 91 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92 Create timestepping solver context 93 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 95 PetscCall(TSSetDM(ts,da)); 96 PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 97 PetscCall(TSSetIFunction(ts,NULL,FormFunction,&ctx)); 98 PetscCall(TSSetMaxTime(ts,.02)); 99 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE)); 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 PetscCall(TSGetSNES(ts,&snes)); 114 PetscCall(DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring)); 115 PetscCall(DMSetMatType(da,MATAIJ)); 116 PetscCall(DMCreateMatrix(da,&J)); 117 PetscCall(MatFDColoringCreate(J,iscoloring,&matfdcoloring)); 118 PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts)); 119 PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 120 PetscCall(MatFDColoringSetUp(J,iscoloring,matfdcoloring)); 121 PetscCall(ISColoringDestroy(&iscoloring)); 122 PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring)); 123 124 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125 Customize nonlinear solver 126 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127 PetscCall(TSSetType(ts,TSBEULER)); 128 129 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130 Set initial conditions 131 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132 PetscCall(FormInitialSolution(da,x,ctx.kappa)); 133 PetscCall(TSSetTimeStep(ts,dt)); 134 PetscCall(TSSetSolution(ts,x)); 135 136 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137 Set runtime options 138 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139 PetscCall(TSSetFromOptions(ts)); 140 141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142 Solve nonlinear system 143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144 PetscCall(TSSolve(ts,x)); 145 PetscCall(TSGetStepNumber(ts,&steps)); 146 147 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148 Free work space. All PETSc objects should be destroyed when they 149 are no longer needed. 150 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 151 PetscCall(MatDestroy(&J)); 152 PetscCall(MatFDColoringDestroy(&matfdcoloring)); 153 PetscCall(VecDestroy(&x)); 154 PetscCall(VecDestroy(&r)); 155 PetscCall(TSDestroy(&ts)); 156 PetscCall(DMDestroy(&da)); 157 158 PetscCall(PetscFinalize()); 159 return 0; 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 PetscInt i,Mx,xs,xm; 179 PetscReal hx,sx; 180 Field *x,*xdot,*f; 181 Vec localX,localXdot; 182 UserCtx *ctx = (UserCtx*)ptr; 183 184 PetscFunctionBegin; 185 PetscCall(TSGetDM(ts,&da)); 186 PetscCall(DMGetLocalVector(da,&localX)); 187 PetscCall(DMGetLocalVector(da,&localXdot)); 188 PetscCall(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 PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 199 PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 200 PetscCall(DMGlobalToLocalBegin(da,Xdot,INSERT_VALUES,localXdot)); 201 PetscCall(DMGlobalToLocalEnd(da,Xdot,INSERT_VALUES,localXdot)); 202 203 /* 204 Get pointers to vector data 205 */ 206 PetscCall(DMDAVecGetArrayRead(da,localX,&x)); 207 PetscCall(DMDAVecGetArrayRead(da,localXdot,&xdot)); 208 PetscCall(DMDAVecGetArray(da,F,&f)); 209 210 /* 211 Get local grid boundaries 212 */ 213 PetscCall(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 (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; 230 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; 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 } 234 } 235 f[i].u = xdot[i].u - (x[i-1].w + x[i+1].w - 2.0*x[i].w)*sx; 236 } 237 238 /* 239 Restore vectors 240 */ 241 PetscCall(DMDAVecRestoreArrayRead(da,localXdot,&xdot)); 242 PetscCall(DMDAVecRestoreArrayRead(da,localX,&x)); 243 PetscCall(DMDAVecRestoreArray(da,F,&f)); 244 PetscCall(DMRestoreLocalVector(da,&localX)); 245 PetscCall(DMRestoreLocalVector(da,&localXdot)); 246 PetscFunctionReturn(0); 247 } 248 249 /* ------------------------------------------------------------------- */ 250 PetscErrorCode FormInitialSolution(DM da,Vec X,PetscReal kappa) 251 { 252 PetscInt i,xs,xm,Mx,xgs,xgm; 253 Field *x; 254 PetscReal hx,xx,r,sx; 255 Vec Xg; 256 257 PetscFunctionBegin; 258 PetscCall(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)); 259 260 hx = 1.0/(PetscReal)Mx; 261 sx = 1.0/(hx*hx); 262 263 /* 264 Get pointers to vector data 265 */ 266 PetscCall(DMCreateLocalVector(da,&Xg)); 267 PetscCall(DMDAVecGetArray(da,Xg,&x)); 268 269 /* 270 Get local grid boundaries 271 */ 272 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 273 PetscCall(DMDAGetGhostCorners(da,&xgs,NULL,NULL,&xgm,NULL,NULL)); 274 275 /* 276 Compute u function over the locally owned part of the grid including ghost points 277 */ 278 for (i=xgs; i<xgs+xgm; 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 /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 284 x[i].w = 0; 285 } 286 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 /* 289 Restore vectors 290 */ 291 PetscCall(DMDAVecRestoreArray(da,Xg,&x)); 292 293 /* Grab only the global part of the vector */ 294 PetscCall(VecSet(X,0)); 295 PetscCall(DMLocalToGlobalBegin(da,Xg,ADD_VALUES,X)); 296 PetscCall(DMLocalToGlobalEnd(da,Xg,ADD_VALUES,X)); 297 PetscCall(VecDestroy(&Xg)); 298 PetscFunctionReturn(0); 299 } 300 301 /*TEST 302 303 build: 304 requires: !complex !single 305 306 test: 307 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 308 requires: x 309 310 TEST*/ 311