1c4762a1bSJed Brown static char help[] = "Model Equations for Advection \n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Modified from ex3.c 5c4762a1bSJed Brown Page 9, Section 1.2 Model Equations for Advection-Diffusion 6c4762a1bSJed Brown 7c4762a1bSJed Brown u_t + a u_x = 0, 0<= x <= 1.0 8c4762a1bSJed Brown 9c4762a1bSJed Brown The initial conditions used here different from the book. 10c4762a1bSJed Brown 11c4762a1bSJed Brown Example: 12c4762a1bSJed Brown ./ex6 -ts_monitor -ts_view_solution -ts_max_steps 100 -ts_monitor_solution draw -draw_pause .1 13c4762a1bSJed Brown ./ex6 -ts_monitor -ts_max_steps 100 -ts_monitor_lg_error -draw_pause .1 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16c4762a1bSJed Brown #include <petscts.h> 17c4762a1bSJed Brown #include <petscdm.h> 18c4762a1bSJed Brown #include <petscdmda.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* 21c4762a1bSJed Brown User-defined application context - contains data needed by the 22c4762a1bSJed Brown application-provided call-back routines. 23c4762a1bSJed Brown */ 24c4762a1bSJed Brown typedef struct { 25c4762a1bSJed Brown PetscReal a; /* advection strength */ 26c4762a1bSJed Brown } AppCtx; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* User-defined routines */ 29c4762a1bSJed Brown extern PetscErrorCode InitialConditions(TS, Vec, AppCtx *); 30c4762a1bSJed Brown extern PetscErrorCode Solution(TS, PetscReal, Vec, AppCtx *); 31c4762a1bSJed Brown extern PetscErrorCode IFunction_LaxFriedrichs(TS, PetscReal, Vec, Vec, Vec, void *); 32c4762a1bSJed Brown extern PetscErrorCode IFunction_LaxWendroff(TS, PetscReal, Vec, Vec, Vec, void *); 33c4762a1bSJed Brown 34d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 35d71ae5a4SJacob Faibussowitsch { 36c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 37c4762a1bSJed Brown TS ts; /* timestepping context */ 38c4762a1bSJed Brown Vec U; /* approximate solution vector */ 39c4762a1bSJed Brown PetscReal dt; 40c4762a1bSJed Brown DM da; 41c4762a1bSJed Brown PetscInt M; 42c4762a1bSJed Brown PetscMPIInt rank; 43c4762a1bSJed Brown PetscBool useLaxWendroff = PETSC_TRUE; 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Initialize program and set problem parameters */ 46327415f7SBarry Smith PetscFunctionBeginUser; 47c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown appctx.a = -1.0; 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-a", &appctx.a, NULL)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da)); 549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 559566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Create vector data structures for approximate and exact solutions */ 589566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &U)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create timestepping solver context */ 619566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 629566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Function evaluation */ 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-useLaxWendroff", &useLaxWendroff, NULL)); 66c4762a1bSJed Brown if (useLaxWendroff) { 6748a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-Wendroff finite volume\n")); 689566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxWendroff, &appctx)); 69c4762a1bSJed Brown } else { 7048a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... Use Lax-LaxFriedrichs finite difference\n")); 719566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction_LaxFriedrichs, &appctx)); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Customize timestepping solver */ 759566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 76c4762a1bSJed Brown dt = 1.0 / (PetscAbsReal(appctx.a) * M); 779566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 789566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 100)); 799566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 100.0)); 809566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 819566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 829566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Evaluate initial conditions */ 859566063dSJacob Faibussowitsch PetscCall(InitialConditions(ts, U, &appctx)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */ 889566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode (*)(TS, PetscReal, Vec, void *))Solution, &appctx)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Run the timestepping solver */ 919566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Free work space */ 949566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 99b122ec5aSJacob Faibussowitsch return 0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 104c4762a1bSJed Brown 105c4762a1bSJed Brown Input Parameter: 106c4762a1bSJed Brown u - uninitialized solution vector (global) 107c4762a1bSJed Brown appctx - user-defined application context 108c4762a1bSJed Brown 109c4762a1bSJed Brown Output Parameter: 110c4762a1bSJed Brown u - vector with solution at initial time (global) 111c4762a1bSJed Brown */ 112d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx) 113d71ae5a4SJacob Faibussowitsch { 114c4762a1bSJed Brown PetscScalar *u; 115c4762a1bSJed Brown PetscInt i, mstart, mend, um, M; 116c4762a1bSJed Brown DM da; 117c4762a1bSJed Brown PetscReal h; 118c4762a1bSJed Brown 1193ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1209566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1219566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 1229566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 123c4762a1bSJed Brown h = 1.0 / M; 124c4762a1bSJed Brown mend = mstart + um; 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Get a pointer to vector data. 127c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 128c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 129c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 130c4762a1bSJed Brown the array. 131c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 132c4762a1bSJed Brown C version. See the users manual for details. 133c4762a1bSJed Brown */ 1349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* 137c4762a1bSJed Brown We initialize the solution array by simply writing the solution 138c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 139c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 140c4762a1bSJed Brown */ 141c4762a1bSJed Brown for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PETSC_PI * i * 6. * h) + 3. * PetscSinReal(PETSC_PI * i * 2. * h); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Restore vector */ 1449566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 148c4762a1bSJed Brown /* 149c4762a1bSJed Brown Solution - Computes the exact solution at a given time 150c4762a1bSJed Brown 151c4762a1bSJed Brown Input Parameters: 152c4762a1bSJed Brown t - current time 153c4762a1bSJed Brown solution - vector in which exact solution will be computed 154c4762a1bSJed Brown appctx - user-defined application context 155c4762a1bSJed Brown 156c4762a1bSJed Brown Output Parameter: 157c4762a1bSJed Brown solution - vector with the newly computed exact solution 158c4762a1bSJed Brown u(x,t) = sin(6*PI*(x - a*t)) + 3 * sin(2*PI*(x - a*t)) 159c4762a1bSJed Brown */ 160d71ae5a4SJacob Faibussowitsch PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx) 161d71ae5a4SJacob Faibussowitsch { 162c4762a1bSJed Brown PetscScalar *u; 163c4762a1bSJed Brown PetscReal a = appctx->a, h, PI6, PI2; 164c4762a1bSJed Brown PetscInt i, mstart, mend, um, M; 165c4762a1bSJed Brown DM da; 166c4762a1bSJed Brown 1673ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1689566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1699566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 1709566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 171c4762a1bSJed Brown h = 1.0 / M; 172c4762a1bSJed Brown mend = mstart + um; 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* Get a pointer to vector data. */ 1759566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */ 178c4762a1bSJed Brown PI6 = PETSC_PI * 6.; 179c4762a1bSJed Brown PI2 = PETSC_PI * 2.; 180ad540459SPierre Jolivet for (i = mstart; i < mend; i++) u[i] = PetscSinReal(PI6 * (i * h - a * t)) + 3. * PetscSinReal(PI2 * (i * h - a * t)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* Restore vector */ 1839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 188c4762a1bSJed Brown /* 189c4762a1bSJed Brown Use Lax-Friedrichs method to evaluate F(u,t) = du/dt + a * du/dx 190c4762a1bSJed Brown 191c4762a1bSJed Brown See https://en.wikipedia.org/wiki/Lax%E2%80%93Friedrichs_method 192c4762a1bSJed Brown */ 193*2a8381b2SBarry Smith PetscErrorCode IFunction_LaxFriedrichs(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx) 194d71ae5a4SJacob Faibussowitsch { 195c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 196c4762a1bSJed Brown PetscInt mstart, mend, M, i, um; 197c4762a1bSJed Brown DM da; 198c4762a1bSJed Brown Vec Uold, localUold; 199c4762a1bSJed Brown PetscScalar *uarray, *f, *uoldarray, h, uave, c; 200c4762a1bSJed Brown PetscReal dt; 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBegin; 2039566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 2049566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &Uold)); 205c4762a1bSJed Brown 2069566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2079566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 2089566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 209c4762a1bSJed Brown h = 1.0 / M; 210c4762a1bSJed Brown mend = mstart + um; 211c4762a1bSJed Brown /* printf(" mstart %d, um %d\n",mstart,um); */ 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localUold)); 2149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold)); 2159566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* Get pointers to vector data */ 2189566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, U, &uarray)); 2199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray)); 2209566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* advection */ 223c4762a1bSJed Brown c = appctx->a * dt / h; /* Courant-Friedrichs-Lewy number (CFL number) */ 224c4762a1bSJed Brown 225c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 226c4762a1bSJed Brown uave = 0.5 * (uoldarray[i - 1] + uoldarray[i + 1]); 227c4762a1bSJed Brown f[i] = uarray[i] - uave + c * 0.5 * (uoldarray[i + 1] - uoldarray[i - 1]); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* Restore vectors */ 2319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray)); 2329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray)); 2339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2349566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localUold)); 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Use Lax-Wendroff method to evaluate F(u,t) = du/dt + a * du/dx 240c4762a1bSJed Brown */ 241*2a8381b2SBarry Smith PetscErrorCode IFunction_LaxWendroff(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx) 242d71ae5a4SJacob Faibussowitsch { 243c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 244c4762a1bSJed Brown PetscInt mstart, mend, M, i, um; 245c4762a1bSJed Brown DM da; 246c4762a1bSJed Brown Vec Uold, localUold; 247c4762a1bSJed Brown PetscScalar *uarray, *f, *uoldarray, h, RFlux, LFlux, lambda; 248c4762a1bSJed Brown PetscReal dt, a; 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 2529566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &Uold)); 253c4762a1bSJed Brown 2549566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &um, 0, 0)); 257c4762a1bSJed Brown h = 1.0 / M; 258c4762a1bSJed Brown mend = mstart + um; 259c4762a1bSJed Brown /* printf(" mstart %d, um %d\n",mstart,um); */ 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localUold)); 2629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, Uold, INSERT_VALUES, localUold)); 2639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, Uold, INSERT_VALUES, localUold)); 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* Get pointers to vector data */ 2669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, U, &uarray)); 2679566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localUold, &uoldarray)); 2689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 269c4762a1bSJed Brown 270c4762a1bSJed Brown /* advection -- finite volume (appctx->a < 0 -- can be relaxed?) */ 271c4762a1bSJed Brown lambda = dt / h; 272c4762a1bSJed Brown a = appctx->a; 273c4762a1bSJed Brown 274c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 275c4762a1bSJed Brown RFlux = 0.5 * a * (uoldarray[i + 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i + 1] - uoldarray[i]); 276c4762a1bSJed Brown LFlux = 0.5 * a * (uoldarray[i - 1] + uoldarray[i]) - a * a * 0.5 * lambda * (uoldarray[i] - uoldarray[i - 1]); 277c4762a1bSJed Brown f[i] = uarray[i] - uoldarray[i] + lambda * (RFlux - LFlux); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Restore vectors */ 2819566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, U, &uarray)); 2829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localUold, &uoldarray)); 2839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2849566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localUold)); 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /*TEST 289c4762a1bSJed Brown 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor 292c4762a1bSJed Brown 293c4762a1bSJed Brown test: 294c4762a1bSJed Brown suffix: 2 295c4762a1bSJed Brown nsize: 3 296c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor 297c4762a1bSJed Brown output_file: output/ex6_1.out 298c4762a1bSJed Brown 299c4762a1bSJed Brown test: 300c4762a1bSJed Brown suffix: 3 301c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false 302c4762a1bSJed Brown 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 4 305c4762a1bSJed Brown nsize: 3 306c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false 307c4762a1bSJed Brown output_file: output/ex6_3.out 308c4762a1bSJed Brown 309c4762a1bSJed Brown TEST*/ 310