1 2 static char help[] = "Model Equations for Advection-Diffusion\n"; 3 4 /* 5 Page 9, Section 1.2 Model Equations for Advection-Diffusion 6 7 u_t = a u_x + d u_xx 8 9 The initial conditions used here different then in the book. 10 11 */ 12 13 /* 14 Helpful runtime linear solver options: 15 -pc_type mg -da_refine 2 -snes_monitor -ksp_monitor -ts_view (geometric multigrid with three levels) 16 17 */ 18 19 /* 20 Include "petscts.h" so that we can use TS solvers. Note that this file 21 automatically includes: 22 petscsys.h - base PETSc routines petscvec.h - vectors 23 petscmat.h - matrices 24 petscis.h - index sets petscksp.h - Krylov subspace methods 25 petscviewer.h - viewers petscpc.h - preconditioners 26 petscksp.h - linear solvers petscsnes.h - nonlinear solvers 27 */ 28 29 #include <petscts.h> 30 #include <petscdm.h> 31 #include <petscdmda.h> 32 33 /* 34 User-defined application context - contains data needed by the 35 application-provided call-back routines. 36 */ 37 typedef struct { 38 PetscScalar a, d; /* advection and diffusion strength */ 39 PetscBool upwind; 40 } AppCtx; 41 42 /* 43 User-defined routines 44 */ 45 extern PetscErrorCode InitialConditions(TS, Vec, AppCtx *); 46 extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *); 47 extern PetscErrorCode Solution(TS, PetscReal, Vec, AppCtx *); 48 49 int main(int argc, char **argv) { 50 AppCtx appctx; /* user-defined application context */ 51 TS ts; /* timestepping context */ 52 Vec U; /* approximate solution vector */ 53 PetscReal dt; 54 DM da; 55 PetscInt M; 56 57 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 58 Initialize program and set problem parameters 59 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 60 61 PetscFunctionBeginUser; 62 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 63 appctx.a = 1.0; 64 appctx.d = 0.0; 65 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-a", &appctx.a, NULL)); 66 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-d", &appctx.d, NULL)); 67 appctx.upwind = PETSC_TRUE; 68 PetscCall(PetscOptionsGetBool(NULL, NULL, "-upwind", &appctx.upwind, NULL)); 69 70 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da)); 71 PetscCall(DMSetFromOptions(da)); 72 PetscCall(DMSetUp(da)); 73 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74 Create vector data structures 75 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 76 77 /* 78 Create vector data structures for approximate and exact solutions 79 */ 80 PetscCall(DMCreateGlobalVector(da, &U)); 81 82 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 83 Create timestepping solver context 84 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 85 86 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 87 PetscCall(TSSetDM(ts, da)); 88 89 /* 90 For linear problems with a time-dependent f(U,t) in the equation 91 u_t = f(u,t), the user provides the discretized right-hand-side 92 as a time-dependent matrix. 93 */ 94 PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx)); 95 PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSMatrixHeat, &appctx)); 96 PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))Solution, &appctx)); 97 98 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99 Customize timestepping solver: 100 - Set timestepping duration info 101 Then set runtime options, which can override these defaults. 102 For example, 103 -ts_max_steps <maxsteps> -ts_max_time <maxtime> 104 to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 105 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 106 107 PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 108 dt = .48 / (M * M); 109 PetscCall(TSSetTimeStep(ts, dt)); 110 PetscCall(TSSetMaxSteps(ts, 1000)); 111 PetscCall(TSSetMaxTime(ts, 100.0)); 112 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 113 PetscCall(TSSetType(ts, TSARKIMEX)); 114 PetscCall(TSSetFromOptions(ts)); 115 116 /* 117 Evaluate initial conditions 118 */ 119 PetscCall(InitialConditions(ts, U, &appctx)); 120 121 /* 122 Run the timestepping solver 123 */ 124 PetscCall(TSSolve(ts, U)); 125 126 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127 Free work space. All PETSc objects should be destroyed when they 128 are no longer needed. 129 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130 131 PetscCall(TSDestroy(&ts)); 132 PetscCall(VecDestroy(&U)); 133 PetscCall(DMDestroy(&da)); 134 135 /* 136 Always call PetscFinalize() before exiting a program. This routine 137 - finalizes the PETSc libraries as well as MPI 138 - provides summary and diagnostic information if certain runtime 139 options are chosen (e.g., -log_view). 140 */ 141 PetscCall(PetscFinalize()); 142 return 0; 143 } 144 /* --------------------------------------------------------------------- */ 145 /* 146 InitialConditions - Computes the solution at the initial time. 147 148 Input Parameter: 149 u - uninitialized solution vector (global) 150 appctx - user-defined application context 151 152 Output Parameter: 153 u - vector with solution at initial time (global) 154 */ 155 PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx) { 156 PetscScalar *u, h; 157 PetscInt i, mstart, mend, xm, M; 158 DM da; 159 160 PetscCall(TSGetDM(ts, &da)); 161 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 162 PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 163 h = 1.0 / M; 164 mend = mstart + xm; 165 /* 166 Get a pointer to vector data. 167 - For default PETSc vectors, VecGetArray() returns a pointer to 168 the data array. Otherwise, the routine is implementation dependent. 169 - You MUST call VecRestoreArray() when you no longer need access to 170 the array. 171 - Note that the Fortran interface to VecGetArray() differs from the 172 C version. See the users manual for details. 173 */ 174 PetscCall(DMDAVecGetArray(da, U, &u)); 175 176 /* 177 We initialize the solution array by simply writing the solution 178 directly into the array locations. Alternatively, we could use 179 VecSetValues() or VecSetValuesLocal(). 180 */ 181 for (i = mstart; i < mend; i++) u[i] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h); 182 183 /* 184 Restore vector 185 */ 186 PetscCall(DMDAVecRestoreArray(da, U, &u)); 187 return 0; 188 } 189 /* --------------------------------------------------------------------- */ 190 /* 191 Solution - Computes the exact solution at a given time. 192 193 Input Parameters: 194 t - current time 195 solution - vector in which exact solution will be computed 196 appctx - user-defined application context 197 198 Output Parameter: 199 solution - vector with the newly computed exact solution 200 */ 201 PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx) { 202 PetscScalar *u, ex1, ex2, sc1, sc2, h; 203 PetscInt i, mstart, mend, xm, M; 204 DM da; 205 206 PetscCall(TSGetDM(ts, &da)); 207 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 208 PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 209 h = 1.0 / M; 210 mend = mstart + xm; 211 /* 212 Get a pointer to vector data. 213 */ 214 PetscCall(DMDAVecGetArray(da, U, &u)); 215 216 /* 217 Simply write the solution directly into the array locations. 218 Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 219 */ 220 ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * appctx->d * t); 221 ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * appctx->d * t); 222 sc1 = PETSC_PI * 6. * h; 223 sc2 = PETSC_PI * 2. * h; 224 for (i = mstart; i < mend; i++) u[i] = PetscSinScalar(sc1 * (PetscReal)i + appctx->a * PETSC_PI * 6. * t) * ex1 + 3. * PetscSinScalar(sc2 * (PetscReal)i + appctx->a * PETSC_PI * 2. * t) * ex2; 225 226 /* 227 Restore vector 228 */ 229 PetscCall(DMDAVecRestoreArray(da, U, &u)); 230 return 0; 231 } 232 233 /* --------------------------------------------------------------------- */ 234 /* 235 RHSMatrixHeat - User-provided routine to compute the right-hand-side 236 matrix for the heat equation. 237 238 Input Parameters: 239 ts - the TS context 240 t - current time 241 global_in - global input vector 242 dummy - optional user-defined context, as set by TSetRHSJacobian() 243 244 Output Parameters: 245 AA - Jacobian matrix 246 BB - optionally different preconditioning matrix 247 str - flag indicating matrix structure 248 249 Notes: 250 Recall that MatSetValues() uses 0-based row and column numbers 251 in Fortran as well as in C. 252 */ 253 PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec U, Mat AA, Mat BB, void *ctx) { 254 Mat A = AA; /* Jacobian matrix */ 255 AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 256 PetscInt mstart, mend; 257 PetscInt i, idx[3], M, xm; 258 PetscScalar v[3], h; 259 DM da; 260 261 PetscCall(TSGetDM(ts, &da)); 262 PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 263 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 264 h = 1.0 / M; 265 mend = mstart + xm; 266 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 267 Compute entries for the locally owned part of the matrix 268 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 269 /* 270 Set matrix rows corresponding to boundary data 271 */ 272 273 /* diffusion */ 274 v[0] = appctx->d / (h * h); 275 v[1] = -2.0 * appctx->d / (h * h); 276 v[2] = appctx->d / (h * h); 277 if (!mstart) { 278 idx[0] = M - 1; 279 idx[1] = 0; 280 idx[2] = 1; 281 PetscCall(MatSetValues(A, 1, &mstart, 3, idx, v, INSERT_VALUES)); 282 mstart++; 283 } 284 285 if (mend == M) { 286 mend--; 287 idx[0] = M - 2; 288 idx[1] = M - 1; 289 idx[2] = 0; 290 PetscCall(MatSetValues(A, 1, &mend, 3, idx, v, INSERT_VALUES)); 291 } 292 293 /* 294 Set matrix rows corresponding to interior data. We construct the 295 matrix one row at a time. 296 */ 297 for (i = mstart; i < mend; i++) { 298 idx[0] = i - 1; 299 idx[1] = i; 300 idx[2] = i + 1; 301 PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES)); 302 } 303 PetscCall(MatAssemblyBegin(A, MAT_FLUSH_ASSEMBLY)); 304 PetscCall(MatAssemblyEnd(A, MAT_FLUSH_ASSEMBLY)); 305 306 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 307 mend = mstart + xm; 308 if (!appctx->upwind) { 309 /* advection -- centered differencing */ 310 v[0] = -.5 * appctx->a / (h); 311 v[1] = .5 * appctx->a / (h); 312 if (!mstart) { 313 idx[0] = M - 1; 314 idx[1] = 1; 315 PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES)); 316 mstart++; 317 } 318 319 if (mend == M) { 320 mend--; 321 idx[0] = M - 2; 322 idx[1] = 0; 323 PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES)); 324 } 325 326 for (i = mstart; i < mend; i++) { 327 idx[0] = i - 1; 328 idx[1] = i + 1; 329 PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES)); 330 } 331 } else { 332 /* advection -- upwinding */ 333 v[0] = -appctx->a / (h); 334 v[1] = appctx->a / (h); 335 if (!mstart) { 336 idx[0] = 0; 337 idx[1] = 1; 338 PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES)); 339 mstart++; 340 } 341 342 if (mend == M) { 343 mend--; 344 idx[0] = M - 1; 345 idx[1] = 0; 346 PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES)); 347 } 348 349 for (i = mstart; i < mend; i++) { 350 idx[0] = i; 351 idx[1] = i + 1; 352 PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES)); 353 } 354 } 355 356 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 357 Complete the matrix assembly process and set some options 358 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 359 /* 360 Assemble matrix, using the 2-step process: 361 MatAssemblyBegin(), MatAssemblyEnd() 362 Computations can be done while messages are in transition 363 by placing code between these two statements. 364 */ 365 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 366 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 367 368 /* 369 Set and option to indicate that we will never add a new nonzero location 370 to the matrix. If we do, it will generate an error. 371 */ 372 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 373 return 0; 374 } 375 376 /*TEST 377 378 test: 379 args: -pc_type mg -da_refine 2 -ts_view -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3 380 requires: double 381 filter: grep -v "total number of" 382 383 test: 384 suffix: 2 385 args: -pc_type mg -da_refine 2 -ts_view -ts_monitor_draw_solution -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3 386 requires: x 387 output_file: output/ex3_1.out 388 requires: double 389 filter: grep -v "total number of" 390 391 TEST*/ 392