1 2 static char help[] = "Model Equations for Advection \n"; 3 4 /* 5 Modified from ex3.c 6 Page 9, Section 1.2 Model Equations for Advection-Diffusion 7 8 u_t + a u_x = 0, 0<= x <= 1.0 9 10 The initial conditions used here different from the book. 11 12 Example: 13 ./ex6 -ts_monitor -ts_view_solution -ts_max_steps 100 -ts_monitor_solution draw -draw_pause .1 14 ./ex6 -ts_monitor -ts_max_steps 100 -ts_monitor_lg_error -draw_pause .1 15 */ 16 17 #include <petscts.h> 18 #include <petscdm.h> 19 #include <petscdmda.h> 20 21 /* 22 User-defined application context - contains data needed by the 23 application-provided call-back routines. 24 */ 25 typedef struct { 26 PetscReal a; /* advection strength */ 27 } AppCtx; 28 29 /* User-defined routines */ 30 extern PetscErrorCode InitialConditions(TS, Vec, AppCtx *); 31 extern PetscErrorCode Solution(TS, PetscReal, Vec, AppCtx *); 32 extern PetscErrorCode IFunction_LaxFriedrichs(TS, PetscReal, Vec, Vec, Vec, void *); 33 extern PetscErrorCode IFunction_LaxWendroff(TS, PetscReal, Vec, Vec, Vec, void *); 34 35 int main(int argc, char **argv) 36 { 37 AppCtx appctx; /* user-defined application context */ 38 TS ts; /* timestepping context */ 39 Vec U; /* approximate solution vector */ 40 PetscReal dt; 41 DM da; 42 PetscInt M; 43 PetscMPIInt rank; 44 PetscBool useLaxWendroff = PETSC_TRUE; 45 46 /* Initialize program and set problem parameters */ 47 PetscFunctionBeginUser; 48 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 49 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 50 51 appctx.a = -1.0; 52 PetscCall(PetscOptionsGetReal(NULL, NULL, "-a", &appctx.a, NULL)); 53 54 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da)); 55 PetscCall(DMSetFromOptions(da)); 56 PetscCall(DMSetUp(da)); 57 58 /* Create vector data structures for approximate and exact solutions */ 59 PetscCall(DMCreateGlobalVector(da, &U)); 60 61 /* Create timestepping solver context */ 62 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 63 PetscCall(TSSetDM(ts, da)); 64 65 /* Function evaluation */ 66 PetscCall(PetscOptionsGetBool(NULL, NULL, "-useLaxWendroff", &useLaxWendroff, NULL)); 67 if (useLaxWendroff) { 68 if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-Wendroff finite volume\n")); 69 PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxWendroff, &appctx)); 70 } else { 71 if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-LaxFriedrichs finite difference\n")); 72 PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxFriedrichs, &appctx)); 73 } 74 75 /* Customize timestepping solver */ 76 PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 77 dt = 1.0 / (PetscAbsReal(appctx.a) * M); 78 PetscCall(TSSetTimeStep(ts, dt)); 79 PetscCall(TSSetMaxSteps(ts, 100)); 80 PetscCall(TSSetMaxTime(ts, 100.0)); 81 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 82 PetscCall(TSSetType(ts, TSBEULER)); 83 PetscCall(TSSetFromOptions(ts)); 84 85 /* Evaluate initial conditions */ 86 PetscCall(InitialConditions(ts, U, &appctx)); 87 88 /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */ 89 PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))Solution, &appctx)); 90 91 /* Run the timestepping solver */ 92 PetscCall(TSSolve(ts, U)); 93 94 /* Free work space */ 95 PetscCall(TSDestroy(&ts)); 96 PetscCall(VecDestroy(&U)); 97 PetscCall(DMDestroy(&da)); 98 99 PetscCall(PetscFinalize()); 100 return 0; 101 } 102 /* --------------------------------------------------------------------- */ 103 /* 104 InitialConditions - Computes the solution at the initial time. 105 106 Input Parameter: 107 u - uninitialized solution vector (global) 108 appctx - user-defined application context 109 110 Output Parameter: 111 u - vector with solution at initial time (global) 112 */ 113 PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx) 114 { 115 PetscScalar *u; 116 PetscInt i, mstart, mend, um, M; 117 DM da; 118 PetscReal h; 119 120 PetscCall(TSGetDM(ts, &da)); 121 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 122 PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 123 h = 1.0 / M; 124 mend = mstart + um; 125 /* 126 Get a pointer to vector data. 127 - For default PETSc vectors, VecGetArray() returns a pointer to 128 the data array. Otherwise, the routine is implementation dependent. 129 - You MUST call VecRestoreArray() when you no longer need access to 130 the array. 131 - Note that the Fortran interface to VecGetArray() differs from the 132 C version. See the users manual for details. 133 */ 134 PetscCall(DMDAVecGetArray(da, U, &u)); 135 136 /* 137 We initialize the solution array by simply writing the solution 138 directly into the array locations. Alternatively, we could use 139 VecSetValues() or VecSetValuesLocal(). 140 */ 141 for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PETSC_PI * i * 6. * h) + 3. * PetscSinReal(PETSC_PI * i * 2. * h); 142 143 /* Restore vector */ 144 PetscCall(DMDAVecRestoreArray(da, U, &u)); 145 return 0; 146 } 147 /* --------------------------------------------------------------------- */ 148 /* 149 Solution - Computes the exact solution at a given time 150 151 Input Parameters: 152 t - current time 153 solution - vector in which exact solution will be computed 154 appctx - user-defined application context 155 156 Output Parameter: 157 solution - vector with the newly computed exact solution 158 u(x,t) = sin(6*PI*(x - a*t)) + 3 * sin(2*PI*(x - a*t)) 159 */ 160 PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx) 161 { 162 PetscScalar *u; 163 PetscReal a = appctx->a, h, PI6, PI2; 164 PetscInt i, mstart, mend, um, M; 165 DM da; 166 167 PetscCall(TSGetDM(ts, &da)); 168 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 169 PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 170 h = 1.0 / M; 171 mend = mstart + um; 172 173 /* Get a pointer to vector data. */ 174 PetscCall(DMDAVecGetArray(da, U, &u)); 175 176 /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */ 177 PI6 = PETSC_PI * 6.; 178 PI2 = PETSC_PI * 2.; 179 for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PI6 * (i * h - a * t)) + 3. * PetscSinReal(PI2 * (i * h - a * t)); 180 181 /* Restore vector */ 182 PetscCall(DMDAVecRestoreArray(da, U, &u)); 183 return 0; 184 } 185 186 /* --------------------------------------------------------------------- */ 187 /* 188 Use Lax-Friedrichs method to evaluate F(u,t) = du/dt + a * du/dx 189 190 See https://en.wikipedia.org/wiki/Lax%E2%80%93Friedrichs_method 191 */ 192 PetscErrorCode IFunction_LaxFriedrichs(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 193 { 194 AppCtx *appctx = (AppCtx *)ctx; 195 PetscInt mstart, mend, M, i, um; 196 DM da; 197 Vec Uold, localUold; 198 PetscScalar *uarray, *f, *uoldarray, h, uave, c; 199 PetscReal dt; 200 201 PetscFunctionBegin; 202 PetscCall(TSGetTimeStep(ts, &dt)); 203 PetscCall(TSGetSolution(ts, &Uold)); 204 205 PetscCall(TSGetDM(ts, &da)); 206 PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 207 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 208 h = 1.0 / M; 209 mend = mstart + um; 210 /* printf(" mstart %d, um %d\n",mstart,um); */ 211 212 PetscCall(DMGetLocalVector(da, &localUold)); 213 PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold)); 214 PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold)); 215 216 /* Get pointers to vector data */ 217 PetscCall(DMDAVecGetArrayRead(da, U, &uarray)); 218 PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray)); 219 PetscCall(DMDAVecGetArray(da, F, &f)); 220 221 /* advection */ 222 c = appctx->a * dt / h; /* Courant-Friedrichs-Lewy number (CFL number) */ 223 224 for (i = mstart; i < mend; i++) { 225 uave = 0.5 * (uoldarray[i - 1] + uoldarray[i + 1]); 226 f[i] = uarray[i] - uave + c * 0.5 * (uoldarray[i + 1] - uoldarray[i - 1]); 227 } 228 229 /* Restore vectors */ 230 PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray)); 231 PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray)); 232 PetscCall(DMDAVecRestoreArray(da, F, &f)); 233 PetscCall(DMRestoreLocalVector(da, &localUold)); 234 PetscFunctionReturn(0); 235 } 236 237 /* 238 Use Lax-Wendroff method to evaluate F(u,t) = du/dt + a * du/dx 239 */ 240 PetscErrorCode IFunction_LaxWendroff(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 241 { 242 AppCtx *appctx = (AppCtx *)ctx; 243 PetscInt mstart, mend, M, i, um; 244 DM da; 245 Vec Uold, localUold; 246 PetscScalar *uarray, *f, *uoldarray, h, RFlux, LFlux, lambda; 247 PetscReal dt, a; 248 249 PetscFunctionBegin; 250 PetscCall(TSGetTimeStep(ts, &dt)); 251 PetscCall(TSGetSolution(ts, &Uold)); 252 253 PetscCall(TSGetDM(ts, &da)); 254 PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 255 PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 256 h = 1.0 / M; 257 mend = mstart + um; 258 /* printf(" mstart %d, um %d\n",mstart,um); */ 259 260 PetscCall(DMGetLocalVector(da, &localUold)); 261 PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold)); 262 PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold)); 263 264 /* Get pointers to vector data */ 265 PetscCall(DMDAVecGetArrayRead(da, U, &uarray)); 266 PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray)); 267 PetscCall(DMDAVecGetArray(da, F, &f)); 268 269 /* advection -- finite volume (appctx->a < 0 -- can be relaxed?) */ 270 lambda = dt / h; 271 a = appctx->a; 272 273 for (i = mstart; i < mend; i++) { 274 RFlux = 0.5 * a * (uoldarray[i + 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i + 1] - uoldarray[i]); 275 LFlux = 0.5 * a * (uoldarray[i - 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i] - uoldarray[i - 1]); 276 f[i] = uarray[i] - uoldarray[i] + lambda * (RFlux - LFlux); 277 } 278 279 /* Restore vectors */ 280 PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray)); 281 PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray)); 282 PetscCall(DMDAVecRestoreArray(da, F, &f)); 283 PetscCall(DMRestoreLocalVector(da, &localUold)); 284 PetscFunctionReturn(0); 285 } 286 287 /*TEST 288 289 test: 290 args: -ts_max_steps 10 -ts_monitor 291 292 test: 293 suffix: 2 294 nsize: 3 295 args: -ts_max_steps 10 -ts_monitor 296 output_file: output/ex6_1.out 297 298 test: 299 suffix: 3 300 args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false 301 302 test: 303 suffix: 4 304 nsize: 3 305 args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false 306 output_file: output/ex6_3.out 307 308 TEST*/ 309