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