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