static char help[] ="Model Equations for Advection \n"; /* Modified from ex3.c Page 9, Section 1.2 Model Equations for Advection-Diffusion u_t + a u_x = 0, 0<= x <= 1.0 The initial conditions used here different from the book. Example: ./ex6 -ts_monitor -ts_view_solution -ts_max_steps 100 -ts_monitor_solution draw -draw_pause .1 ./ex6 -ts_monitor -ts_max_steps 100 -ts_monitor_lg_error -draw_pause .1 */ #include #include #include /* User-defined application context - contains data needed by the application-provided call-back routines. */ typedef struct { PetscReal a; /* advection strength */ } AppCtx; /* User-defined routines */ extern PetscErrorCode InitialConditions(TS,Vec,AppCtx*); extern PetscErrorCode Solution(TS,PetscReal,Vec,AppCtx*); extern PetscErrorCode IFunction_LaxFriedrichs(TS,PetscReal,Vec,Vec,Vec,void*); extern PetscErrorCode IFunction_LaxWendroff(TS,PetscReal,Vec,Vec,Vec,void*); int main(int argc,char **argv) { AppCtx appctx; /* user-defined application context */ TS ts; /* timestepping context */ Vec U; /* approximate solution vector */ PetscErrorCode ierr; PetscReal dt; DM da; PetscInt M; PetscMPIInt rank; PetscBool useLaxWendroff = PETSC_TRUE; /* Initialize program and set problem parameters */ ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); appctx.a = -1.0; CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-a",&appctx.a,NULL)); CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, 60, 1, 1,NULL,&da)); CHKERRQ(DMSetFromOptions(da)); CHKERRQ(DMSetUp(da)); /* Create vector data structures for approximate and exact solutions */ CHKERRQ(DMCreateGlobalVector(da,&U)); /* Create timestepping solver context */ CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); CHKERRQ(TSSetDM(ts,da)); /* Function evaluation */ CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-useLaxWendroff",&useLaxWendroff,NULL)); if (useLaxWendroff) { if (rank == 0) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"... Use Lax-Wendroff finite volume\n")); } CHKERRQ(TSSetIFunction(ts,NULL,IFunction_LaxWendroff,&appctx)); } else { if (rank == 0) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"... Use Lax-LaxFriedrichs finite difference\n")); } CHKERRQ(TSSetIFunction(ts,NULL,IFunction_LaxFriedrichs,&appctx)); } /* Customize timestepping solver */ CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0)); dt = 1.0/(PetscAbsReal(appctx.a)*M); CHKERRQ(TSSetTimeStep(ts,dt)); CHKERRQ(TSSetMaxSteps(ts,100)); CHKERRQ(TSSetMaxTime(ts,100.0)); CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); CHKERRQ(TSSetType(ts,TSBEULER)); CHKERRQ(TSSetFromOptions(ts)); /* Evaluate initial conditions */ CHKERRQ(InitialConditions(ts,U,&appctx)); /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */ CHKERRQ(TSSetSolutionFunction(ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void*))Solution,&appctx)); /* Run the timestepping solver */ CHKERRQ(TSSolve(ts,U)); /* Free work space */ CHKERRQ(TSDestroy(&ts)); CHKERRQ(VecDestroy(&U)); CHKERRQ(DMDestroy(&da)); ierr = PetscFinalize(); return ierr; } /* --------------------------------------------------------------------- */ /* InitialConditions - Computes the solution at the initial time. Input Parameter: u - uninitialized solution vector (global) appctx - user-defined application context Output Parameter: u - vector with solution at initial time (global) */ PetscErrorCode InitialConditions(TS ts,Vec U,AppCtx *appctx) { PetscScalar *u; PetscInt i,mstart,mend,um,M; DM da; PetscReal h; CHKERRQ(TSGetDM(ts,&da)); CHKERRQ(DMDAGetCorners(da,&mstart,0,0,&um,0,0)); CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0)); h = 1.0/M; mend = mstart + um; /* Get a pointer to vector data. - For default PETSc vectors, VecGetArray() returns a pointer to the data array. Otherwise, the routine is implementation dependent. - You MUST call VecRestoreArray() when you no longer need access to the array. - Note that the Fortran interface to VecGetArray() differs from the C version. See the users manual for details. */ CHKERRQ(DMDAVecGetArray(da,U,&u)); /* We initialize the solution array by simply writing the solution directly into the array locations. Alternatively, we could use VecSetValues() or VecSetValuesLocal(). */ for (i=mstart; ia,h,PI6,PI2; PetscInt i,mstart,mend,um,M; DM da; CHKERRQ(TSGetDM(ts,&da)); CHKERRQ(DMDAGetCorners(da,&mstart,0,0,&um,0,0)); CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0)); h = 1.0/M; mend = mstart + um; /* Get a pointer to vector data. */ CHKERRQ(DMDAVecGetArray(da,U,&u)); /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */ PI6 = PETSC_PI*6.; PI2 = PETSC_PI*2.; for (i=mstart; ia*dt/h; /* Courant-Friedrichs-Lewy number (CFL number) */ for (i=mstart; ia < 0 -- can be relaxed?) */ lambda = dt/h; a = appctx->a; for (i=mstart; i