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 PetscReal dt; 41c4762a1bSJed Brown DM da; 42c4762a1bSJed Brown PetscInt M; 43c4762a1bSJed Brown PetscMPIInt rank; 44c4762a1bSJed Brown PetscBool useLaxWendroff = PETSC_TRUE; 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Initialize program and set problem parameters */ 47*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 485f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown appctx.a = -1.0; 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-a",&appctx.a,NULL)); 52c4762a1bSJed Brown 535f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, 60, 1, 1,NULL,&da)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Create vector data structures for approximate and exact solutions */ 585f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&U)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create timestepping solver context */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Function evaluation */ 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-useLaxWendroff",&useLaxWendroff,NULL)); 66c4762a1bSJed Brown if (useLaxWendroff) { 67dd400576SPatrick Sanan if (rank == 0) { 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"... Use Lax-Wendroff finite volume\n")); 69c4762a1bSJed Brown } 705f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunction_LaxWendroff,&appctx)); 71c4762a1bSJed Brown } else { 72dd400576SPatrick Sanan if (rank == 0) { 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"... Use Lax-LaxFriedrichs finite difference\n")); 74c4762a1bSJed Brown } 755f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunction_LaxFriedrichs,&appctx)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Customize timestepping solver */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0)); 80c4762a1bSJed Brown dt = 1.0/(PetscAbsReal(appctx.a)*M); 815f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,100)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,100.0)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSBEULER)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Evaluate initial conditions */ 895f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(ts,U,&appctx)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolutionFunction(ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void*))Solution,&appctx)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Run the timestepping solver */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Free work space */ 985f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 101c4762a1bSJed Brown 102*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 103*b122ec5aSJacob Faibussowitsch return 0; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 108c4762a1bSJed Brown 109c4762a1bSJed Brown Input Parameter: 110c4762a1bSJed Brown u - uninitialized solution vector (global) 111c4762a1bSJed Brown appctx - user-defined application context 112c4762a1bSJed Brown 113c4762a1bSJed Brown Output Parameter: 114c4762a1bSJed Brown u - vector with solution at initial time (global) 115c4762a1bSJed Brown */ 116c4762a1bSJed Brown PetscErrorCode InitialConditions(TS ts,Vec U,AppCtx *appctx) 117c4762a1bSJed Brown { 118c4762a1bSJed Brown PetscScalar *u; 119c4762a1bSJed Brown PetscInt i,mstart,mend,um,M; 120c4762a1bSJed Brown DM da; 121c4762a1bSJed Brown PetscReal h; 122c4762a1bSJed Brown 1235f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&mstart,0,0,&um,0,0)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0)); 126c4762a1bSJed Brown h = 1.0/M; 127c4762a1bSJed Brown mend = mstart + um; 128c4762a1bSJed Brown /* 129c4762a1bSJed Brown Get a pointer to vector data. 130c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 131c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 132c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 133c4762a1bSJed Brown the array. 134c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 135c4762a1bSJed Brown C version. See the users manual for details. 136c4762a1bSJed Brown */ 1375f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* 140c4762a1bSJed Brown We initialize the solution array by simply writing the solution 141c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 142c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 143c4762a1bSJed Brown */ 144c4762a1bSJed Brown for (i=mstart; i<mend; i++) u[i] = PetscSinReal(PETSC_PI*i*6.*h) + 3.*PetscSinReal(PETSC_PI*i*2.*h); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* Restore vector */ 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 148c4762a1bSJed Brown return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 151c4762a1bSJed Brown /* 152c4762a1bSJed Brown Solution - Computes the exact solution at a given time 153c4762a1bSJed Brown 154c4762a1bSJed Brown Input Parameters: 155c4762a1bSJed Brown t - current time 156c4762a1bSJed Brown solution - vector in which exact solution will be computed 157c4762a1bSJed Brown appctx - user-defined application context 158c4762a1bSJed Brown 159c4762a1bSJed Brown Output Parameter: 160c4762a1bSJed Brown solution - vector with the newly computed exact solution 161c4762a1bSJed Brown u(x,t) = sin(6*PI*(x - a*t)) + 3 * sin(2*PI*(x - a*t)) 162c4762a1bSJed Brown */ 163c4762a1bSJed Brown PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *appctx) 164c4762a1bSJed Brown { 165c4762a1bSJed Brown PetscScalar *u; 166c4762a1bSJed Brown PetscReal a=appctx->a,h,PI6,PI2; 167c4762a1bSJed Brown PetscInt i,mstart,mend,um,M; 168c4762a1bSJed Brown DM da; 169c4762a1bSJed Brown 1705f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&mstart,0,0,&um,0,0)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0)); 173c4762a1bSJed Brown h = 1.0/M; 174c4762a1bSJed Brown mend = mstart + um; 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Get a pointer to vector data. */ 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */ 180c4762a1bSJed Brown PI6 = PETSC_PI*6.; 181c4762a1bSJed Brown PI2 = PETSC_PI*2.; 182c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 183c4762a1bSJed Brown u[i] = PetscSinReal(PI6*(i*h - a*t)) + 3.*PetscSinReal(PI2*(i*h - a*t)); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Restore vector */ 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 188c4762a1bSJed Brown return 0; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 192c4762a1bSJed Brown /* 193c4762a1bSJed Brown Use Lax-Friedrichs method to evaluate F(u,t) = du/dt + a * du/dx 194c4762a1bSJed Brown 195c4762a1bSJed Brown See https://en.wikipedia.org/wiki/Lax%E2%80%93Friedrichs_method 196c4762a1bSJed Brown */ 197c4762a1bSJed Brown PetscErrorCode IFunction_LaxFriedrichs(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void* ctx) 198c4762a1bSJed Brown { 199c4762a1bSJed Brown AppCtx *appctx=(AppCtx*)ctx; 200c4762a1bSJed Brown PetscInt mstart,mend,M,i,um; 201c4762a1bSJed Brown DM da; 202c4762a1bSJed Brown Vec Uold,localUold; 203c4762a1bSJed Brown PetscScalar *uarray,*f,*uoldarray,h,uave,c; 204c4762a1bSJed Brown PetscReal dt; 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscFunctionBegin; 2075f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts,&Uold)); 209c4762a1bSJed Brown 2105f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&mstart,0,0,&um,0,0)); 213c4762a1bSJed Brown h = 1.0/M; 214c4762a1bSJed Brown mend = mstart + um; 215c4762a1bSJed Brown /* printf(" mstart %d, um %d\n",mstart,um); */ 216c4762a1bSJed Brown 2175f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localUold)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,Uold,INSERT_VALUES,localUold)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,Uold,INSERT_VALUES,localUold)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* Get pointers to vector data */ 2225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,U,&uarray)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localUold,&uoldarray)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* advection */ 227c4762a1bSJed Brown c = appctx->a*dt/h; /* Courant-Friedrichs-Lewy number (CFL number) */ 228c4762a1bSJed Brown 229c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 230c4762a1bSJed Brown uave = 0.5*(uoldarray[i-1] + uoldarray[i+1]); 231c4762a1bSJed Brown f[i] = uarray[i] - uave + c*0.5*(uoldarray[i+1] - uoldarray[i-1]); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* Restore vectors */ 2355f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,U,&uarray)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localUold,&uoldarray)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localUold)); 239c4762a1bSJed Brown PetscFunctionReturn(0); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* 243c4762a1bSJed Brown Use Lax-Wendroff method to evaluate F(u,t) = du/dt + a * du/dx 244c4762a1bSJed Brown */ 245c4762a1bSJed Brown PetscErrorCode IFunction_LaxWendroff(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void* ctx) 246c4762a1bSJed Brown { 247c4762a1bSJed Brown AppCtx *appctx=(AppCtx*)ctx; 248c4762a1bSJed Brown PetscInt mstart,mend,M,i,um; 249c4762a1bSJed Brown DM da; 250c4762a1bSJed Brown Vec Uold,localUold; 251c4762a1bSJed Brown PetscScalar *uarray,*f,*uoldarray,h,RFlux,LFlux,lambda; 252c4762a1bSJed Brown PetscReal dt,a; 253c4762a1bSJed Brown 254c4762a1bSJed Brown PetscFunctionBegin; 2555f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts,&Uold)); 257c4762a1bSJed Brown 2585f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&mstart,0,0,&um,0,0)); 261c4762a1bSJed Brown h = 1.0/M; 262c4762a1bSJed Brown mend = mstart + um; 263c4762a1bSJed Brown /* printf(" mstart %d, um %d\n",mstart,um); */ 264c4762a1bSJed Brown 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localUold)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,Uold,INSERT_VALUES,localUold)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,Uold,INSERT_VALUES,localUold)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* Get pointers to vector data */ 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,U,&uarray)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localUold,&uoldarray)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* advection -- finite volume (appctx->a < 0 -- can be relaxed?) */ 275c4762a1bSJed Brown lambda = dt/h; 276c4762a1bSJed Brown a = appctx->a; 277c4762a1bSJed Brown 278c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 279c4762a1bSJed Brown RFlux = 0.5 * a * (uoldarray[i+1] + uoldarray[i]) - a*a*0.5*lambda * (uoldarray[i+1] - uoldarray[i]); 280c4762a1bSJed Brown LFlux = 0.5 * a * (uoldarray[i-1] + uoldarray[i]) - a*a*0.5*lambda * (uoldarray[i] - uoldarray[i-1]); 281c4762a1bSJed Brown f[i] = uarray[i] - uoldarray[i] + lambda * (RFlux - LFlux); 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown /* Restore vectors */ 2855f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,U,&uarray)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localUold,&uoldarray)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localUold)); 289c4762a1bSJed Brown PetscFunctionReturn(0); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown /*TEST 293c4762a1bSJed Brown 294c4762a1bSJed Brown test: 295c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor 296c4762a1bSJed Brown 297c4762a1bSJed Brown test: 298c4762a1bSJed Brown suffix: 2 299c4762a1bSJed Brown nsize: 3 300c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor 301c4762a1bSJed Brown output_file: output/ex6_1.out 302c4762a1bSJed Brown 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 3 305c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false 306c4762a1bSJed Brown 307c4762a1bSJed Brown test: 308c4762a1bSJed Brown suffix: 4 309c4762a1bSJed Brown nsize: 3 310c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false 311c4762a1bSJed Brown output_file: output/ex6_3.out 312c4762a1bSJed Brown 313c4762a1bSJed Brown TEST*/ 314