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