1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d. Simplified from ex7.c for illustrating how to use TS on a structured domain. \n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown u_t = uxx + uyy 5c4762a1bSJed Brown 0 < x < 1, 0 < y < 1; 6c4762a1bSJed Brown At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125 7c4762a1bSJed Brown u(x,y) = 0.0 if r >= .125 8c4762a1bSJed Brown 9c4762a1bSJed Brown mpiexec -n 2 ./ex13 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor 10c4762a1bSJed Brown mpiexec -n 1 ./ex13 -snes_fd_color -ts_monitor_draw_solution 11c4762a1bSJed Brown mpiexec -n 2 ./ex13 -ts_type sundials -ts_monitor 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include <petscdm.h> 15c4762a1bSJed Brown #include <petscdmda.h> 16c4762a1bSJed Brown #include <petscts.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* 19c4762a1bSJed Brown User-defined data structures and routines 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscReal c; 23c4762a1bSJed Brown } AppCtx; 24c4762a1bSJed Brown 25c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 26c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 27c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(DM,Vec,void*); 28c4762a1bSJed Brown 29c4762a1bSJed Brown int main(int argc,char **argv) 30c4762a1bSJed Brown { 31c4762a1bSJed Brown TS ts; /* nonlinear solver */ 32c4762a1bSJed Brown Vec u,r; /* solution, residual vector */ 33c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 34c4762a1bSJed Brown PetscInt steps; /* iterations for convergence */ 35c4762a1bSJed Brown DM da; 36c4762a1bSJed Brown PetscReal ftime,dt; 37c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 38c4762a1bSJed Brown 39*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 40c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 41c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 42c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 48c4762a1bSJed Brown Extract global vectors from DMDA; 49c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 505f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&u)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&r)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Initialize user application context */ 54c4762a1bSJed Brown user.c = -30.0; 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 57c4762a1bSJed Brown Create timestepping solver context 58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 595f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSBEULER)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,r,RHSFunction,&user)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Set Jacobian */ 655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATAIJ)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&J)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,J,J,RHSJacobian,NULL)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown ftime = 1.0; 705f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,ftime)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74c4762a1bSJed Brown Set initial conditions 75c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 765f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialSolution(da,u,&user)); 77c4762a1bSJed Brown dt = .01; 785f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81c4762a1bSJed Brown Set runtime options 82c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 86c4762a1bSJed Brown Solve nonlinear system 87c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 885f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,u)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93c4762a1bSJed Brown Free work space. 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 100c4762a1bSJed Brown 101*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 102*b122ec5aSJacob Faibussowitsch return 0; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown RHSFunction - Evaluates nonlinear function, F(u). 107c4762a1bSJed Brown 108c4762a1bSJed Brown Input Parameters: 109c4762a1bSJed Brown . ts - the TS context 110c4762a1bSJed Brown . U - input vector 111c4762a1bSJed Brown . ptr - optional user-defined context, as set by TSSetFunction() 112c4762a1bSJed Brown 113c4762a1bSJed Brown Output Parameter: 114c4762a1bSJed Brown . F - function vector 115c4762a1bSJed Brown */ 116c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr) 117c4762a1bSJed Brown { 118c4762a1bSJed Brown /* PETSC_UNUSED AppCtx *user=(AppCtx*)ptr; */ 119c4762a1bSJed Brown DM da; 120c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 121c4762a1bSJed Brown PetscReal two = 2.0,hx,hy,sx,sy; 122c4762a1bSJed Brown PetscScalar u,uxx,uyy,**uarray,**f; 123c4762a1bSJed Brown Vec localU; 124c4762a1bSJed Brown 125c4762a1bSJed Brown PetscFunctionBeginUser; 1265f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localU)); 1285f80ce2aSJacob Faibussowitsch 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)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); 131c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* 134c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 135c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 136c4762a1bSJed Brown By placing code between these two statements, computations can be 137c4762a1bSJed Brown done while messages are in transition. 138c4762a1bSJed Brown */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Get pointers to vector data */ 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localU,&uarray)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* Get local grid boundaries */ 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 150c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 151c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 152c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 153c4762a1bSJed Brown f[j][i] = uarray[j][i]; 154c4762a1bSJed Brown continue; 155c4762a1bSJed Brown } 156c4762a1bSJed Brown u = uarray[j][i]; 157c4762a1bSJed Brown uxx = (-two*u + uarray[j][i-1] + uarray[j][i+1])*sx; 158c4762a1bSJed Brown uyy = (-two*u + uarray[j-1][i] + uarray[j+1][i])*sy; 159c4762a1bSJed Brown f[j][i] = uxx + uyy; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* Restore vectors */ 1645f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&uarray)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localU)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(11.0*ym*xm)); 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown RHSJacobian - User-provided routine to compute the Jacobian of 174c4762a1bSJed Brown the nonlinear right-hand-side function of the ODE. 175c4762a1bSJed Brown 176c4762a1bSJed Brown Input Parameters: 177c4762a1bSJed Brown ts - the TS context 178c4762a1bSJed Brown t - current time 179c4762a1bSJed Brown U - global input vector 180c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 181c4762a1bSJed Brown 182c4762a1bSJed Brown Output Parameters: 183c4762a1bSJed Brown J - Jacobian matrix 184c4762a1bSJed Brown Jpre - optionally different preconditioning matrix 185c4762a1bSJed Brown str - flag indicating matrix structure 186c4762a1bSJed Brown */ 187c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat Jpre,void *ctx) 188c4762a1bSJed Brown { 189c4762a1bSJed Brown DM da; 190c4762a1bSJed Brown DMDALocalInfo info; 191c4762a1bSJed Brown PetscInt i,j; 192c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 193c4762a1bSJed Brown 194c4762a1bSJed Brown PetscFunctionBeginUser; 1955f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 197c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); sx = 1.0/(hx*hx); 198c4762a1bSJed Brown hy = 1.0/(PetscReal)(info.my-1); sy = 1.0/(hy*hy); 199c4762a1bSJed Brown for (j=info.ys; j<info.ys+info.ym; j++) { 200c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 201c4762a1bSJed Brown PetscInt nc = 0; 202c4762a1bSJed Brown MatStencil row,col[5]; 203c4762a1bSJed Brown PetscScalar val[5]; 204c4762a1bSJed Brown row.i = i; row.j = j; 205c4762a1bSJed Brown if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) { 206c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = 1.0; 207c4762a1bSJed Brown } else { 208c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = sx; 209c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = sx; 210c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = sy; 211c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = sy; 212c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -2*sx - 2*sy; 213c4762a1bSJed Brown } 2145f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES)); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown } 2175f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 219c4762a1bSJed Brown if (J != Jpre) { 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown PetscFunctionReturn(0); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 227c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U,void* ptr) 228c4762a1bSJed Brown { 229c4762a1bSJed Brown AppCtx *user=(AppCtx*)ptr; 230c4762a1bSJed Brown PetscReal c=user->c; 231c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 232c4762a1bSJed Brown PetscScalar **u; 233c4762a1bSJed Brown PetscReal hx,hy,x,y,r; 234c4762a1bSJed Brown 235c4762a1bSJed Brown PetscFunctionBeginUser; 2365f80ce2aSJacob Faibussowitsch 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)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 239c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* Get pointers to vector data */ 2425f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Get local grid boundaries */ 2455f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 248c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 249c4762a1bSJed Brown y = j*hy; 250c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 251c4762a1bSJed Brown x = i*hx; 252c4762a1bSJed Brown r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)); 253c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(c*r*r*r); 254c4762a1bSJed Brown else u[j][i] = 0.0; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* Restore vectors */ 2595f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 260c4762a1bSJed Brown PetscFunctionReturn(0); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /*TEST 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown args: -ts_max_steps 5 -ts_monitor 267c4762a1bSJed Brown 268c4762a1bSJed Brown test: 269c4762a1bSJed Brown suffix: 2 270c4762a1bSJed Brown args: -ts_max_steps 5 -ts_monitor 271c4762a1bSJed Brown 272c4762a1bSJed Brown test: 273c4762a1bSJed Brown suffix: 3 274c4762a1bSJed Brown args: -ts_max_steps 5 -snes_fd_color -ts_monitor 275c4762a1bSJed Brown 276c4762a1bSJed Brown TEST*/ 277