1 static const char help[] = "Time-dependent PDE in 1d. Simplified from ex15.c for illustrating how to solve DAEs. \n"; 2 /* 3 u_t = uxx 4 0 < x < 1; 5 At t=0: u(x) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5)) < .125 6 u(x) = 0.0 if r >= .125 7 8 Boundary conditions: 9 Dirichlet BC: 10 At x=0, x=1, u = 0.0 11 12 Neumann BC: 13 At x=0, x=1: du(x,t)/dx = 0 14 15 mpiexec -n 2 ./ex17 -da_grid_x 40 -ts_max_steps 2 -snes_monitor -ksp_monitor 16 ./ex17 -da_grid_x 40 -monitor_solution 17 ./ex17 -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5 # Midpoint is not L-stable 18 ./ex17 -jac_type fd_coloring -da_grid_x 500 -boundary 1 19 ./ex17 -da_grid_x 100 -ts_type gl -ts_adapt_type none -ts_max_steps 2 20 */ 21 22 #include <petscdm.h> 23 #include <petscdmda.h> 24 #include <petscts.h> 25 26 typedef enum {JACOBIAN_ANALYTIC,JACOBIAN_FD_COLORING,JACOBIAN_FD_FULL} JacobianType; 27 static const char *const JacobianTypes[] = {"analytic","fd_coloring","fd_full","JacobianType","fd_",0}; 28 29 /* 30 User-defined data structures and routines 31 */ 32 typedef struct { 33 PetscReal c; 34 PetscInt boundary; /* Type of boundary condition */ 35 PetscBool viewJacobian; 36 } AppCtx; 37 38 static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 39 static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 40 static PetscErrorCode FormInitialSolution(TS,Vec,void*); 41 42 int main(int argc,char **argv) 43 { 44 TS ts; /* nonlinear solver */ 45 Vec u; /* solution, residual vectors */ 46 Mat J; /* Jacobian matrix */ 47 PetscInt nsteps; 48 PetscReal vmin,vmax,norm; 49 PetscErrorCode ierr; 50 DM da; 51 PetscReal ftime,dt; 52 AppCtx user; /* user-defined work context */ 53 JacobianType jacType; 54 55 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 56 57 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 58 Create distributed array (DMDA) to manage parallel grid and vectors 59 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 60 CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,1,1,NULL,&da)); 61 CHKERRQ(DMSetFromOptions(da)); 62 CHKERRQ(DMSetUp(da)); 63 64 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65 Extract global vectors from DMDA; 66 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67 CHKERRQ(DMCreateGlobalVector(da,&u)); 68 69 /* Initialize user application context */ 70 user.c = -30.0; 71 user.boundary = 0; /* 0: Dirichlet BC; 1: Neumann BC */ 72 user.viewJacobian = PETSC_FALSE; 73 74 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL)); 75 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian)); 76 77 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78 Create timestepping solver context 79 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 80 CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 81 CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 82 CHKERRQ(TSSetType(ts,TSTHETA)); 83 CHKERRQ(TSThetaSetTheta(ts,1.0)); /* Make the Theta method behave like backward Euler */ 84 CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,&user)); 85 86 CHKERRQ(DMSetMatType(da,MATAIJ)); 87 CHKERRQ(DMCreateMatrix(da,&J)); 88 jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */ 89 90 CHKERRQ(TSSetDM(ts,da)); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */ 91 92 ftime = 1.0; 93 CHKERRQ(TSSetMaxTime(ts,ftime)); 94 CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 95 96 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 97 Set initial conditions 98 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 99 CHKERRQ(FormInitialSolution(ts,u,&user)); 100 CHKERRQ(TSSetSolution(ts,u)); 101 dt = .01; 102 CHKERRQ(TSSetTimeStep(ts,dt)); 103 104 /* Use slow fd Jacobian or fast fd Jacobian with colorings. 105 Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */ 106 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);CHKERRQ(ierr); 107 CHKERRQ(PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0)); 108 ierr = PetscOptionsEnd();CHKERRQ(ierr); 109 if (jacType == JACOBIAN_ANALYTIC) { 110 CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,&user)); 111 } else if (jacType == JACOBIAN_FD_COLORING) { 112 SNES snes; 113 CHKERRQ(TSGetSNES(ts,&snes)); 114 CHKERRQ(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0)); 115 } else if (jacType == JACOBIAN_FD_FULL) { 116 SNES snes; 117 CHKERRQ(TSGetSNES(ts,&snes)); 118 CHKERRQ(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user)); 119 } 120 121 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122 Set runtime options 123 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 124 CHKERRQ(TSSetFromOptions(ts)); 125 126 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127 Integrate ODE system 128 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 129 CHKERRQ(TSSolve(ts,u)); 130 131 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132 Compute diagnostics of the solution 133 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134 CHKERRQ(VecNorm(u,NORM_1,&norm)); 135 CHKERRQ(VecMax(u,NULL,&vmax)); 136 CHKERRQ(VecMin(u,NULL,&vmin)); 137 CHKERRQ(TSGetStepNumber(ts,&nsteps)); 138 CHKERRQ(TSGetTime(ts,&ftime)); 139 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"timestep %D: time %g, solution norm %g, max %g, min %g\n",nsteps,(double)ftime,(double)norm,(double)vmax,(double)vmin)); 140 141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142 Free work space. 143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144 CHKERRQ(MatDestroy(&J)); 145 CHKERRQ(VecDestroy(&u)); 146 CHKERRQ(TSDestroy(&ts)); 147 CHKERRQ(DMDestroy(&da)); 148 ierr = PetscFinalize(); 149 return ierr; 150 } 151 /* ------------------------------------------------------------------- */ 152 static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr) 153 { 154 AppCtx *user=(AppCtx*)ptr; 155 DM da; 156 PetscInt i,Mx,xs,xm; 157 PetscReal hx,sx; 158 PetscScalar *u,*udot,*f; 159 Vec localU; 160 161 PetscFunctionBeginUser; 162 CHKERRQ(TSGetDM(ts,&da)); 163 CHKERRQ(DMGetLocalVector(da,&localU)); 164 CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 165 166 hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); 167 168 /* 169 Scatter ghost points to local vector,using the 2-step process 170 DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 171 By placing code between these two statements, computations can be 172 done while messages are in transition. 173 */ 174 CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); 175 CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); 176 177 /* Get pointers to vector data */ 178 CHKERRQ(DMDAVecGetArrayRead(da,localU,&u)); 179 CHKERRQ(DMDAVecGetArrayRead(da,Udot,&udot)); 180 CHKERRQ(DMDAVecGetArray(da,F,&f)); 181 182 /* Get local grid boundaries */ 183 CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 184 185 /* Compute function over the locally owned part of the grid */ 186 for (i=xs; i<xs+xm; i++) { 187 if (user->boundary == 0) { /* Dirichlet BC */ 188 if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */ 189 else f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx; 190 } else { /* Neumann BC */ 191 if (i == 0) f[i] = u[0] - u[1]; 192 else if (i == Mx-1) f[i] = u[i] - u[i-1]; 193 else f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx; 194 } 195 } 196 197 /* Restore vectors */ 198 CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u)); 199 CHKERRQ(DMDAVecRestoreArrayRead(da,Udot,&udot)); 200 CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 201 CHKERRQ(DMRestoreLocalVector(da,&localU)); 202 PetscFunctionReturn(0); 203 } 204 205 /* --------------------------------------------------------------------- */ 206 /* 207 IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 208 */ 209 PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx) 210 { 211 PetscInt i,rstart,rend,Mx; 212 PetscReal hx,sx; 213 AppCtx *user = (AppCtx*)ctx; 214 DM da; 215 MatStencil col[3],row; 216 PetscInt nc; 217 PetscScalar vals[3]; 218 219 PetscFunctionBeginUser; 220 CHKERRQ(TSGetDM(ts,&da)); 221 CHKERRQ(MatGetOwnershipRange(Jpre,&rstart,&rend)); 222 CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 223 hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); 224 for (i=rstart; i<rend; i++) { 225 nc = 0; 226 row.i = i; 227 if (user->boundary == 0 && (i == 0 || i == Mx-1)) { 228 col[nc].i = i; vals[nc++] = 1.0; 229 } else if (user->boundary > 0 && i == 0) { /* Left Neumann */ 230 col[nc].i = i; vals[nc++] = 1.0; 231 col[nc].i = i+1; vals[nc++] = -1.0; 232 } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */ 233 col[nc].i = i-1; vals[nc++] = -1.0; 234 col[nc].i = i; vals[nc++] = 1.0; 235 } else { /* Interior */ 236 col[nc].i = i-1; vals[nc++] = -1.0*sx; 237 col[nc].i = i; vals[nc++] = 2.0*sx + a; 238 col[nc].i = i+1; vals[nc++] = -1.0*sx; 239 } 240 CHKERRQ(MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES)); 241 } 242 243 CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 244 CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 245 if (J != Jpre) { 246 CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 247 CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 248 } 249 if (user->viewJacobian) { 250 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n")); 251 CHKERRQ(MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD)); 252 } 253 PetscFunctionReturn(0); 254 } 255 256 /* ------------------------------------------------------------------- */ 257 PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr) 258 { 259 AppCtx *user=(AppCtx*)ptr; 260 PetscReal c =user->c; 261 DM da; 262 PetscInt i,xs,xm,Mx; 263 PetscScalar *u; 264 PetscReal hx,x,r; 265 266 PetscFunctionBeginUser; 267 CHKERRQ(TSGetDM(ts,&da)); 268 CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 269 270 hx = 1.0/(PetscReal)(Mx-1); 271 272 /* Get pointers to vector data */ 273 CHKERRQ(DMDAVecGetArray(da,U,&u)); 274 275 /* Get local grid boundaries */ 276 CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 277 278 /* Compute function over the locally owned part of the grid */ 279 for (i=xs; i<xs+xm; i++) { 280 x = i*hx; 281 r = PetscSqrtReal((x-.5)*(x-.5)); 282 if (r < .125) u[i] = PetscExpReal(c*r*r*r); 283 else u[i] = 0.0; 284 } 285 286 /* Restore vectors */ 287 CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 288 PetscFunctionReturn(0); 289 } 290 291 /*TEST 292 293 test: 294 requires: !single 295 args: -da_grid_x 40 -ts_max_steps 2 -snes_monitor_short -ksp_monitor_short -ts_monitor 296 297 test: 298 suffix: 2 299 requires: !single 300 args: -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5 301 302 TEST*/ 303