1 static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n"; 2 /* 3 u_t + a1*u_x = -k1*u + k2*v + s1 4 v_t + a2*v_x = k1*u - k2*v + s2 5 0 < x < 1; 6 a1 = 1, k1 = 10^6, s1 = 0, 7 a2 = 0, k2 = 2*k1, s2 = 1 8 9 Initial conditions: 10 u(x,0) = 1 + s2*x 11 v(x,0) = k0/k1*u(x,0) + s1/k1 12 13 Upstream boundary conditions: 14 u(0,t) = 1-sin(12*t)^4 15 16 Useful command-line parameters: 17 -ts_arkimex_fully_implicit - treats advection implicitly, only relevant with -ts_type arkimex (default) 18 -ts_type rosw - use Rosenbrock-W method (linearly implicit IMEX, amortizes assembly and preconditioner setup) 19 */ 20 21 #include <petscdm.h> 22 #include <petscdmda.h> 23 #include <petscts.h> 24 25 typedef PetscScalar Field[2]; 26 27 typedef struct _User *User; 28 struct _User { 29 PetscReal a[2]; /* Advection speeds */ 30 PetscReal k[2]; /* Reaction coefficients */ 31 PetscReal s[2]; /* Source terms */ 32 }; 33 34 static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 35 static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 36 static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 37 static PetscErrorCode FormInitialSolution(TS,Vec,void*); 38 39 int main(int argc,char **argv) 40 { 41 TS ts; /* time integrator */ 42 SNES snes; /* nonlinear solver */ 43 SNESLineSearch linesearch; /* line search */ 44 Vec X; /* solution, residual vectors */ 45 Mat J; /* Jacobian matrix */ 46 PetscInt steps,mx; 47 DM da; 48 PetscReal ftime,dt; 49 struct _User user; /* user-defined work context */ 50 TSConvergedReason reason; 51 52 PetscFunctionBeginUser; 53 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 54 55 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56 Create distributed array (DMDA) to manage parallel grid and vectors 57 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 58 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da)); 59 PetscCall(DMSetFromOptions(da)); 60 PetscCall(DMSetUp(da)); 61 62 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 63 Extract global vectors from DMDA; 64 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 65 PetscCall(DMCreateGlobalVector(da,&X)); 66 67 /* Initialize user application context */ 68 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options",""); 69 { 70 user.a[0] = 1; PetscCall(PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL)); 71 user.a[1] = 0; PetscCall(PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL)); 72 user.k[0] = 1e6; PetscCall(PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL)); 73 user.k[1] = 2*user.k[0]; PetscCall(PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL)); 74 user.s[0] = 0; PetscCall(PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL)); 75 user.s[1] = 1; PetscCall(PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL)); 76 } 77 PetscOptionsEnd(); 78 79 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80 Create timestepping solver context 81 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 82 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 83 PetscCall(TSSetDM(ts,da)); 84 PetscCall(TSSetType(ts,TSARKIMEX)); 85 PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 86 PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user)); 87 PetscCall(DMSetMatType(da,MATAIJ)); 88 PetscCall(DMCreateMatrix(da,&J)); 89 PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user)); 90 91 /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since 92 * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use 93 * SNESSetType(snes,SNESKSPONLY). */ 94 PetscCall(TSGetSNES(ts,&snes)); 95 PetscCall(SNESGetLineSearch(snes,&linesearch)); 96 PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 97 98 ftime = .1; 99 PetscCall(TSSetMaxTime(ts,ftime)); 100 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 101 102 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103 Set initial conditions 104 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 105 PetscCall(FormInitialSolution(ts,X,&user)); 106 PetscCall(TSSetSolution(ts,X)); 107 PetscCall(VecGetSize(X,&mx)); 108 dt = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */ 109 PetscCall(TSSetTimeStep(ts,dt)); 110 111 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112 Set runtime options 113 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 114 PetscCall(TSSetFromOptions(ts)); 115 116 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117 Solve nonlinear system 118 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119 PetscCall(TSSolve(ts,X)); 120 PetscCall(TSGetSolveTime(ts,&ftime)); 121 PetscCall(TSGetStepNumber(ts,&steps)); 122 PetscCall(TSGetConvergedReason(ts,&reason)); 123 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %" PetscInt_FMT " steps\n",TSConvergedReasons[reason],(double)ftime,steps)); 124 125 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126 Free work space. 127 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 128 PetscCall(MatDestroy(&J)); 129 PetscCall(VecDestroy(&X)); 130 PetscCall(TSDestroy(&ts)); 131 PetscCall(DMDestroy(&da)); 132 PetscCall(PetscFinalize()); 133 return 0; 134 } 135 136 static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 137 { 138 User user = (User)ptr; 139 DM da; 140 DMDALocalInfo info; 141 PetscInt i; 142 Field *f; 143 const Field *x,*xdot; 144 145 PetscFunctionBeginUser; 146 PetscCall(TSGetDM(ts,&da)); 147 PetscCall(DMDAGetLocalInfo(da,&info)); 148 149 /* Get pointers to vector data */ 150 PetscCall(DMDAVecGetArrayRead(da,X,(void*)&x)); 151 PetscCall(DMDAVecGetArrayRead(da,Xdot,(void*)&xdot)); 152 PetscCall(DMDAVecGetArray(da,F,&f)); 153 154 /* Compute function over the locally owned part of the grid */ 155 for (i=info.xs; i<info.xs+info.xm; i++) { 156 f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0]; 157 f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1]; 158 } 159 160 /* Restore vectors */ 161 PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&x)); 162 PetscCall(DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot)); 163 PetscCall(DMDAVecRestoreArray(da,F,&f)); 164 PetscFunctionReturn(0); 165 } 166 167 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 168 { 169 User user = (User)ptr; 170 DM da; 171 Vec Xloc; 172 DMDALocalInfo info; 173 PetscInt i,j; 174 PetscReal hx; 175 Field *f; 176 const Field *x; 177 178 PetscFunctionBeginUser; 179 PetscCall(TSGetDM(ts,&da)); 180 PetscCall(DMDAGetLocalInfo(da,&info)); 181 hx = 1.0/(PetscReal)info.mx; 182 183 /* 184 Scatter ghost points to local vector,using the 2-step process 185 DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 186 By placing code between these two statements, computations can be 187 done while messages are in transition. 188 */ 189 PetscCall(DMGetLocalVector(da,&Xloc)); 190 PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 191 PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 192 193 /* Get pointers to vector data */ 194 PetscCall(DMDAVecGetArrayRead(da,Xloc,(void*)&x)); 195 PetscCall(DMDAVecGetArray(da,F,&f)); 196 197 /* Compute function over the locally owned part of the grid */ 198 for (i=info.xs; i<info.xs+info.xm; i++) { 199 const PetscReal *a = user->a; 200 PetscReal u0t[2]; 201 u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12*t),4); 202 u0t[1] = 0.0; 203 for (j=0; j<2; j++) { 204 if (i == 0) f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]); 205 else if (i == 1) f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]); 206 else if (i == info.mx-2) f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]); 207 else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]); 208 else f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]); 209 } 210 } 211 212 /* Restore vectors */ 213 PetscCall(DMDAVecRestoreArrayRead(da,Xloc,(void*)&x)); 214 PetscCall(DMDAVecRestoreArray(da,F,&f)); 215 PetscCall(DMRestoreLocalVector(da,&Xloc)); 216 PetscFunctionReturn(0); 217 } 218 219 /* --------------------------------------------------------------------- */ 220 /* 221 IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 222 */ 223 PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 224 { 225 User user = (User)ptr; 226 DMDALocalInfo info; 227 PetscInt i; 228 DM da; 229 const Field *x,*xdot; 230 231 PetscFunctionBeginUser; 232 PetscCall(TSGetDM(ts,&da)); 233 PetscCall(DMDAGetLocalInfo(da,&info)); 234 235 /* Get pointers to vector data */ 236 PetscCall(DMDAVecGetArrayRead(da,X,(void*)&x)); 237 PetscCall(DMDAVecGetArrayRead(da,Xdot,(void*)&xdot)); 238 239 /* Compute function over the locally owned part of the grid */ 240 for (i=info.xs; i<info.xs+info.xm; i++) { 241 const PetscReal *k = user->k; 242 PetscScalar v[2][2]; 243 244 v[0][0] = a + k[0]; v[0][1] = -k[1]; 245 v[1][0] = -k[0]; v[1][1] = a+k[1]; 246 PetscCall(MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES)); 247 } 248 249 /* Restore vectors */ 250 PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&x)); 251 PetscCall(DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot)); 252 253 PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 254 PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 255 if (J != Jpre) { 256 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 257 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 258 } 259 PetscFunctionReturn(0); 260 } 261 262 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 263 { 264 User user = (User)ctx; 265 DM da; 266 PetscInt i; 267 DMDALocalInfo info; 268 Field *x; 269 PetscReal hx; 270 271 PetscFunctionBeginUser; 272 PetscCall(TSGetDM(ts,&da)); 273 PetscCall(DMDAGetLocalInfo(da,&info)); 274 hx = 1.0/(PetscReal)info.mx; 275 276 /* Get pointers to vector data */ 277 PetscCall(DMDAVecGetArray(da,X,&x)); 278 279 /* Compute function over the locally owned part of the grid */ 280 for (i=info.xs; i<info.xs+info.xm; i++) { 281 PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0; 282 x[i][0] = 1 + user->s[1]*r; 283 x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik; 284 } 285 PetscCall(DMDAVecRestoreArray(da,X,&x)); 286 PetscFunctionReturn(0); 287 } 288 289 /*TEST 290 291 test: 292 args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_arkimex_type 4 -ts_adapt_type none -ts_dt .005 -ts_max_time .1 293 requires: !single 294 295 test: 296 suffix: 2 297 args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_dt 1e-3 -ts_adapt_type none -ts_dt .005 -ts_max_time .1 298 nsize: 2 299 300 test: 301 suffix: 3 302 args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type ra34pw2 -ts_dt 5e-3 -ts_max_time .1 -ts_adapt_type none 303 nsize: 2 304 305 test: 306 suffix: 4 307 args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_dt 0.001 -ts_monitor -ts_max_steps 100 308 filter: sed "s/ITS/TIME/g" 309 nsize: 2 310 311 TEST*/ 312