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