static char help[] = "Time-dependent PDE in 2d. Simplified from ex7.c for illustrating how to use TS on a structured domain. \n"; /* u_t = uxx + uyy 0 < x < 1, 0 < y < 1; At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125 u(x,y) = 0.0 if r >= .125 mpiexec -n 2 ./ex13 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor mpiexec -n 1 ./ex13 -snes_fd_color -ts_monitor_draw_solution mpiexec -n 2 ./ex13 -ts_type sundials -ts_monitor */ #include #include #include /* User-defined data structures and routines */ typedef struct { PetscReal c; } AppCtx; extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); extern PetscErrorCode FormInitialSolution(DM,Vec,void*); int main(int argc,char **argv) { TS ts; /* nonlinear solver */ Vec u,r; /* solution, residual vector */ Mat J; /* Jacobian matrix */ PetscInt steps; /* iterations for convergence */ DM da; PetscReal ftime,dt; AppCtx user; /* user-defined work context */ PetscFunctionBeginUser; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create distributed array (DMDA) to manage parallel grid and vectors - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); PetscCall(DMSetFromOptions(da)); PetscCall(DMSetUp(da)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Extract global vectors from DMDA; - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(DMCreateGlobalVector(da,&u)); PetscCall(VecDuplicate(u,&r)); /* Initialize user application context */ user.c = -30.0; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create timestepping solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); PetscCall(TSSetDM(ts,da)); PetscCall(TSSetType(ts,TSBEULER)); PetscCall(TSSetRHSFunction(ts,r,RHSFunction,&user)); /* Set Jacobian */ PetscCall(DMSetMatType(da,MATAIJ)); PetscCall(DMCreateMatrix(da,&J)); PetscCall(TSSetRHSJacobian(ts,J,J,RHSJacobian,NULL)); ftime = 1.0; PetscCall(TSSetMaxTime(ts,ftime)); PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set initial conditions - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(FormInitialSolution(da,u,&user)); dt = .01; PetscCall(TSSetTimeStep(ts,dt)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(TSSetFromOptions(ts)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve nonlinear system - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(TSSolve(ts,u)); PetscCall(TSGetSolveTime(ts,&ftime)); PetscCall(TSGetStepNumber(ts,&steps)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Free work space. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ PetscCall(MatDestroy(&J)); PetscCall(VecDestroy(&u)); PetscCall(VecDestroy(&r)); PetscCall(TSDestroy(&ts)); PetscCall(DMDestroy(&da)); PetscCall(PetscFinalize()); return 0; } /* ------------------------------------------------------------------- */ /* RHSFunction - Evaluates nonlinear function, F(u). Input Parameters: . ts - the TS context . U - input vector . ptr - optional user-defined context, as set by TSSetFunction() Output Parameter: . F - function vector */ PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr) { /* PETSC_UNUSED AppCtx *user=(AppCtx*)ptr; */ DM da; PetscInt i,j,Mx,My,xs,ys,xm,ym; PetscReal two = 2.0,hx,hy,sx,sy; PetscScalar u,uxx,uyy,**uarray,**f; Vec localU; PetscFunctionBeginUser; PetscCall(TSGetDM(ts,&da)); PetscCall(DMGetLocalVector(da,&localU)); PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy); /* Scatter ghost points to local vector,using the 2-step process DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). By placing code between these two statements, computations can be done while messages are in transition. */ PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); /* Get pointers to vector data */ PetscCall(DMDAVecGetArrayRead(da,localU,&uarray)); PetscCall(DMDAVecGetArray(da,F,&f)); /* Get local grid boundaries */ PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); /* Compute function over the locally owned part of the grid */ for (j=ys; jc; PetscInt i,j,xs,ys,xm,ym,Mx,My; PetscScalar **u; PetscReal hx,hy,x,y,r; PetscFunctionBeginUser; PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); hx = 1.0/(PetscReal)(Mx-1); hy = 1.0/(PetscReal)(My-1); /* Get pointers to vector data */ PetscCall(DMDAVecGetArray(da,U,&u)); /* Get local grid boundaries */ PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); /* Compute function over the locally owned part of the grid */ for (j=ys; j