static char help[] = "Nonlinear, time-dependent PDE in 2d.\n"; /* Solves the equation u_tt - \Delta u = 0 which we split into two first-order systems u_t - v = 0 so that F(u,v).u = v v_t - \Delta u = 0 so that F(u,v).v = Delta u Include "petscdmda.h" so that we can use distributed arrays (DMDAs). Include "petscts.h" so that we can use SNES solvers. Note that this file automatically includes: petscsys.h - base PETSc routines petscvec.h - vectors petscmat.h - matrices petscis.h - index sets petscksp.h - Krylov subspace methods petscviewer.h - viewers petscpc.h - preconditioners petscksp.h - linear solvers */ #include #include #include /* User-defined routines */ extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec); extern PetscErrorCode MyTSMonitor(TS,PetscInt,PetscReal,Vec,void*); extern PetscErrorCode MySNESMonitor(SNES,PetscInt,PetscReal,PetscViewerAndFormat*); int main(int argc,char **argv) { TS ts; /* nonlinear solver */ Vec x,r; /* solution, residual vectors */ PetscInt steps; /* iterations for convergence */ DM da; PetscReal ftime; SNES ts_snes; PetscBool usemonitor = PETSC_TRUE; PetscViewerAndFormat *vf; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize program - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-usemonitor",&usemonitor,NULL)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create distributed array (DMDA) to manage parallel grid and vectors - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da)); CHKERRQ(DMSetFromOptions(da)); CHKERRQ(DMSetUp(da)); CHKERRQ(DMDASetFieldName(da,0,"u")); CHKERRQ(DMDASetFieldName(da,1,"v")); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Extract global vectors from DMDA; then duplicate for remaining vectors that are the same types - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(DMCreateGlobalVector(da,&x)); CHKERRQ(VecDuplicate(x,&r)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create timestepping solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); CHKERRQ(TSSetDM(ts,da)); CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); CHKERRQ(TSSetRHSFunction(ts,NULL,FormFunction,da)); CHKERRQ(TSSetMaxTime(ts,1.0)); if (usemonitor) { CHKERRQ(TSMonitorSet(ts,MyTSMonitor,0,0)); } /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Customize nonlinear solver - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSSetType(ts,TSBEULER)); CHKERRQ(TSGetSNES(ts,&ts_snes)); if (usemonitor) { CHKERRQ(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf)); CHKERRQ(SNESMonitorSet(ts_snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void *))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy)); } /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set initial conditions - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(FormInitialSolution(da,x)); CHKERRQ(TSSetTimeStep(ts,.0001)); CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); CHKERRQ(TSSetSolution(ts,x)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSSetFromOptions(ts)); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve nonlinear system - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(TSSolve(ts,x)); CHKERRQ(TSGetSolveTime(ts,&ftime)); CHKERRQ(TSGetStepNumber(ts,&steps)); CHKERRQ(VecViewFromOptions(x,NULL,"-final_sol")); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Free work space. All PETSc objects should be destroyed when they are no longer needed. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); CHKERRQ(TSDestroy(&ts)); CHKERRQ(DMDestroy(&da)); CHKERRQ(PetscFinalize()); return 0; } /* ------------------------------------------------------------------- */ /* FormFunction - Evaluates nonlinear function, F(x). Input Parameters: . ts - the TS context . X - input vector . ptr - optional user-defined context, as set by SNESSetFunction() Output Parameter: . F - function vector */ PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr) { DM da = (DM)ptr; PetscInt i,j,Mx,My,xs,ys,xm,ym; PetscReal hx,hy,/*hxdhy,hydhx,*/ sx,sy; PetscScalar u,uxx,uyy,v,***x,***f; Vec localX; PetscFunctionBeginUser; CHKERRQ(DMGetLocalVector(da,&localX)); CHKERRQ(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); /*hxdhy = hx/hy;*/ /*hydhx = hy/hx;*/ /* 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. */ CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); /* Get pointers to vector data */ CHKERRQ(DMDAVecGetArrayDOF(da,localX,&x)); CHKERRQ(DMDAVecGetArrayDOF(da,F,&f)); /* Get local grid boundaries */ CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); /* Compute function over the locally owned part of the grid */ for (j=ys; j -1) { /* -1 is used to indicate an interpolated value */ CHKERRQ(PetscPrintf(comm,"timestep %D time %g norm %g\n",step,(double)ptime,(double)norm)); } PetscFunctionReturn(0); } /* MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES. Input Parameters: snes - the SNES context its - iteration number fnorm - 2-norm function value (may be estimated) ctx - optional user-defined context for private data for the monitor routine, as set by SNESMonitorSet() */ PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf) { PetscFunctionBeginUser; CHKERRQ(SNESMonitorDefaultShort(snes,its,fnorm,vf)); PetscFunctionReturn(0); } /*TEST test: args: -da_grid_x 20 -ts_max_time 3 -ts_dt 1e-1 -ts_theta_initial_guess_extrapolate 0 -ts_monitor -ksp_monitor_short requires: !single test: suffix: 2 args: -da_grid_x 20 -ts_max_time 0.11 -ts_dt 1e-1 -ts_type glle -ts_monitor -ksp_monitor_short requires: !single test: suffix: glvis_da_2d_vect args: -usemonitor 0 -da_grid_x 20 -ts_max_time 0.3 -ts_dt 1e-1 -ts_type glle -final_sol glvis: -viewer_glvis_dm_da_bs 2,0 requires: !single test: suffix: glvis_da_2d_vect_ll args: -usemonitor 0 -da_grid_x 20 -ts_max_time 0.3 -ts_dt 1e-1 -ts_type glle -final_sol glvis: -viewer_glvis_dm_da_bs 2,0 -viewer_glvis_dm_da_ll requires: !single TEST*/