1 static const char help[] = "Call PetscInitialize multiple times.\n"; 2 /* 3 This example is based on the Brusselator tutorial of the same name, but tests multiple calls to PetscInitialize(). 4 This is a bad "convergence study" because it only compares min and max values of the solution rather than comparing 5 norms of the errors. For convergence studies, we recommend invoking PetscInitialize() only once and comparing norms 6 of errors (perhaps estimated using an accurate reference solution). 7 8 Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and multiple solves. 9 10 u_t - alpha u_xx = A + u^2 v - (B+1) u 11 v_t - alpha v_xx = B u - u^2 v 12 0 < x < 1; 13 A = 1, B = 3, alpha = 1/10 14 15 Initial conditions: 16 u(x,0) = 1 + sin(2 pi x) 17 v(x,0) = 3 18 19 Boundary conditions: 20 u(0,t) = u(1,t) = 1 21 v(0,t) = v(1,t) = 3 22 */ 23 24 #include <petscdm.h> 25 #include <petscdmda.h> 26 #include <petscts.h> 27 28 typedef struct { 29 PetscScalar u,v; 30 } Field; 31 32 typedef struct _User *User; 33 struct _User { 34 PetscReal A,B; /* Reaction coefficients */ 35 PetscReal alpha; /* Diffusion coefficient */ 36 PetscReal uleft,uright; /* Dirichlet boundary conditions */ 37 PetscReal vleft,vright; /* Dirichlet boundary conditions */ 38 }; 39 40 static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 41 static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 42 static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 43 static PetscErrorCode FormInitialSolution(TS,Vec,void*); 44 static int Brusselator(int,char**,PetscInt); 45 46 int main(int argc,char **argv) 47 { 48 PetscInt cycle; 49 PetscErrorCode ierr; 50 51 ierr = MPI_Init(&argc,&argv);if (ierr) return ierr; 52 for (cycle=0; cycle<4; cycle++) { 53 ierr = Brusselator(argc,argv,cycle); 54 if (ierr) return 1; 55 } 56 ierr = MPI_Finalize(); 57 return ierr; 58 } 59 60 PetscErrorCode Brusselator(int argc,char **argv,PetscInt cycle) 61 { 62 TS ts; /* nonlinear solver */ 63 Vec X; /* solution, residual vectors */ 64 Mat J; /* Jacobian matrix */ 65 PetscInt steps,mx; 66 PetscErrorCode ierr; 67 DM da; 68 PetscReal ftime,hx,dt,xmax,xmin; 69 struct _User user; /* user-defined work context */ 70 TSConvergedReason reason; 71 72 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 73 74 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75 Create distributed array (DMDA) to manage parallel grid and vectors 76 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 77 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da)); 78 PetscCall(DMSetFromOptions(da)); 79 PetscCall(DMSetUp(da)); 80 81 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82 Extract global vectors from DMDA; 83 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 84 PetscCall(DMCreateGlobalVector(da,&X)); 85 86 /* Initialize user application context */ 87 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");PetscCall(ierr); 88 { 89 user.A = 1; 90 user.B = 3; 91 user.alpha = 0.1; 92 user.uleft = 1; 93 user.uright = 1; 94 user.vleft = 3; 95 user.vright = 3; 96 PetscCall(PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL)); 97 PetscCall(PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL)); 98 PetscCall(PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL)); 99 PetscCall(PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL)); 100 PetscCall(PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL)); 101 PetscCall(PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL)); 102 PetscCall(PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL)); 103 } 104 ierr = PetscOptionsEnd();PetscCall(ierr); 105 106 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107 Create timestepping solver context 108 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 110 PetscCall(TSSetDM(ts,da)); 111 PetscCall(TSSetType(ts,TSARKIMEX)); 112 PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 113 PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user)); 114 PetscCall(DMSetMatType(da,MATAIJ)); 115 PetscCall(DMCreateMatrix(da,&J)); 116 PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user)); 117 118 ftime = 1.0; 119 PetscCall(TSSetMaxTime(ts,ftime)); 120 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 121 122 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123 Set initial conditions 124 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125 PetscCall(FormInitialSolution(ts,X,&user)); 126 PetscCall(TSSetSolution(ts,X)); 127 PetscCall(VecGetSize(X,&mx)); 128 hx = 1.0/(PetscReal)(mx/2-1); 129 dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */ 130 dt *= PetscPowRealInt(0.2,cycle); /* Shrink the time step in convergence study. */ 131 PetscCall(TSSetTimeStep(ts,dt)); 132 PetscCall(TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL)); 133 134 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135 Set runtime options 136 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137 PetscCall(TSSetFromOptions(ts)); 138 139 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140 Solve nonlinear system 141 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142 PetscCall(TSSolve(ts,X)); 143 PetscCall(TSGetSolveTime(ts,&ftime)); 144 PetscCall(TSGetStepNumber(ts,&steps)); 145 PetscCall(TSGetConvergedReason(ts,&reason)); 146 PetscCall(VecMin(X,NULL,&xmin)); 147 PetscCall(VecMax(X,NULL,&xmax)); 148 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after % 3D steps. Range [%6.4f,%6.4f]\n",TSConvergedReasons[reason],(double)ftime,steps,(double)xmin,(double)xmax)); 149 150 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151 Free work space. 152 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153 PetscCall(MatDestroy(&J)); 154 PetscCall(VecDestroy(&X)); 155 PetscCall(TSDestroy(&ts)); 156 PetscCall(DMDestroy(&da)); 157 PetscCall(PetscFinalize()); 158 return 0; 159 } 160 161 static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 162 { 163 User user = (User)ptr; 164 DM da; 165 DMDALocalInfo info; 166 PetscInt i; 167 Field *x,*xdot,*f; 168 PetscReal hx; 169 Vec Xloc; 170 171 PetscFunctionBeginUser; 172 PetscCall(TSGetDM(ts,&da)); 173 PetscCall(DMDAGetLocalInfo(da,&info)); 174 hx = 1.0/(PetscReal)(info.mx-1); 175 176 /* 177 Scatter ghost points to local vector,using the 2-step process 178 DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 179 By placing code between these two statements, computations can be 180 done while messages are in transition. 181 */ 182 PetscCall(DMGetLocalVector(da,&Xloc)); 183 PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 184 PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 185 186 /* Get pointers to vector data */ 187 PetscCall(DMDAVecGetArrayRead(da,Xloc,&x)); 188 PetscCall(DMDAVecGetArrayRead(da,Xdot,&xdot)); 189 PetscCall(DMDAVecGetArray(da,F,&f)); 190 191 /* Compute function over the locally owned part of the grid */ 192 for (i=info.xs; i<info.xs+info.xm; i++) { 193 if (i == 0) { 194 f[i].u = hx * (x[i].u - user->uleft); 195 f[i].v = hx * (x[i].v - user->vleft); 196 } else if (i == info.mx-1) { 197 f[i].u = hx * (x[i].u - user->uright); 198 f[i].v = hx * (x[i].v - user->vright); 199 } else { 200 f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx; 201 f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx; 202 } 203 } 204 205 /* Restore vectors */ 206 PetscCall(DMDAVecRestoreArrayRead(da,Xloc,&x)); 207 PetscCall(DMDAVecRestoreArrayRead(da,Xdot,&xdot)); 208 PetscCall(DMDAVecRestoreArray(da,F,&f)); 209 PetscCall(DMRestoreLocalVector(da,&Xloc)); 210 PetscFunctionReturn(0); 211 } 212 213 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 214 { 215 User user = (User)ptr; 216 DM da; 217 DMDALocalInfo info; 218 PetscInt i; 219 PetscReal hx; 220 Field *x,*f; 221 222 PetscFunctionBeginUser; 223 PetscCall(TSGetDM(ts,&da)); 224 PetscCall(DMDAGetLocalInfo(da,&info)); 225 hx = 1.0/(PetscReal)(info.mx-1); 226 227 /* Get pointers to vector data */ 228 PetscCall(DMDAVecGetArrayRead(da,X,&x)); 229 PetscCall(DMDAVecGetArray(da,F,&f)); 230 231 /* Compute function over the locally owned part of the grid */ 232 for (i=info.xs; i<info.xs+info.xm; i++) { 233 PetscScalar u = x[i].u,v = x[i].v; 234 f[i].u = hx*(user->A + u*u*v - (user->B+1)*u); 235 f[i].v = hx*(user->B*u - u*u*v); 236 } 237 238 /* Restore vectors */ 239 PetscCall(DMDAVecRestoreArrayRead(da,X,&x)); 240 PetscCall(DMDAVecRestoreArray(da,F,&f)); 241 PetscFunctionReturn(0); 242 } 243 244 /* --------------------------------------------------------------------- */ 245 /* 246 IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 247 */ 248 PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 249 { 250 User user = (User)ptr; 251 DMDALocalInfo info; 252 PetscInt i; 253 PetscReal hx; 254 DM da; 255 Field *x,*xdot; 256 257 PetscFunctionBeginUser; 258 PetscCall(TSGetDM(ts,&da)); 259 PetscCall(DMDAGetLocalInfo(da,&info)); 260 hx = 1.0/(PetscReal)(info.mx-1); 261 262 /* Get pointers to vector data */ 263 PetscCall(DMDAVecGetArrayRead(da,X,&x)); 264 PetscCall(DMDAVecGetArrayRead(da,Xdot,&xdot)); 265 266 /* Compute function over the locally owned part of the grid */ 267 for (i=info.xs; i<info.xs+info.xm; i++) { 268 if (i == 0 || i == info.mx-1) { 269 const PetscInt row = i,col = i; 270 const PetscScalar vals[2][2] = {{hx,0},{0,hx}}; 271 PetscCall(MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES)); 272 } else { 273 const PetscInt row = i,col[] = {i-1,i,i+1}; 274 const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx; 275 const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}}, 276 {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}}; 277 PetscCall(MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES)); 278 } 279 } 280 281 /* Restore vectors */ 282 PetscCall(DMDAVecRestoreArrayRead(da,X,&x)); 283 PetscCall(DMDAVecRestoreArrayRead(da,Xdot,&xdot)); 284 285 PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 286 PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 287 if (J != Jpre) { 288 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 289 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 290 } 291 PetscFunctionReturn(0); 292 } 293 294 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 295 { 296 User user = (User)ctx; 297 DM da; 298 PetscInt i; 299 DMDALocalInfo info; 300 Field *x; 301 PetscReal hx; 302 303 PetscFunctionBeginUser; 304 PetscCall(TSGetDM(ts,&da)); 305 PetscCall(DMDAGetLocalInfo(da,&info)); 306 hx = 1.0/(PetscReal)(info.mx-1); 307 308 /* Get pointers to vector data */ 309 PetscCall(DMDAVecGetArray(da,X,&x)); 310 311 /* Compute function over the locally owned part of the grid */ 312 for (i=info.xs; i<info.xs+info.xm; i++) { 313 PetscReal xi = i*hx; 314 x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi); 315 x[i].v = user->vleft*(1.-xi) + user->vright*xi; 316 } 317 PetscCall(DMDAVecRestoreArray(da,X,&x)); 318 PetscFunctionReturn(0); 319 } 320 321 /*TEST 322 323 test: 324 args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3 325 326 test: 327 suffix: 2 328 args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3 329 330 TEST*/ 331