1 2 static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\ 3 Input parameters include:\n\ 4 -m <points>, where <points> = number of grid points\n\ 5 -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\ 6 -debug : Activate debugging printouts\n\ 7 -nox : Deactivate x-window graphics\n\n"; 8 9 /* 10 Concepts: TS^time-dependent linear problems 11 Concepts: TS^heat equation 12 Concepts: TS^diffusion equation 13 Processors: n 14 */ 15 16 /* ------------------------------------------------------------------------ 17 18 This program solves the one-dimensional heat equation (also called the 19 diffusion equation), 20 u_t = u_xx, 21 on the domain 0 <= x <= 1, with the boundary conditions 22 u(t,0) = 0, u(t,1) = 0, 23 and the initial condition 24 u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x). 25 This is a linear, second-order, parabolic equation. 26 27 We discretize the right-hand side using finite differences with 28 uniform grid spacing h: 29 u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) 30 We then demonstrate time evolution using the various TS methods by 31 running the program via 32 mpiexec -n <procs> ex3 -ts_type <timestepping solver> 33 34 We compare the approximate solution with the exact solution, given by 35 u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) + 36 3*exp(-4*pi*pi*t) * sin(2*pi*x) 37 38 Notes: 39 This code demonstrates the TS solver interface to two variants of 40 linear problems, u_t = f(u,t), namely 41 - time-dependent f: f(u,t) is a function of t 42 - time-independent f: f(u,t) is simply f(u) 43 44 The uniprocessor version of this code is ts/tutorials/ex3.c 45 46 ------------------------------------------------------------------------- */ 47 48 /* 49 Include "petscdmda.h" so that we can use distributed arrays (DMDAs) to manage 50 the parallel grid. Include "petscts.h" so that we can use TS solvers. 51 Note that this file automatically includes: 52 petscsys.h - base PETSc routines petscvec.h - vectors 53 petscmat.h - matrices 54 petscis.h - index sets petscksp.h - Krylov subspace methods 55 petscviewer.h - viewers petscpc.h - preconditioners 56 petscksp.h - linear solvers petscsnes.h - nonlinear solvers 57 */ 58 59 #include <petscdm.h> 60 #include <petscdmda.h> 61 #include <petscts.h> 62 #include <petscdraw.h> 63 64 /* 65 User-defined application context - contains data needed by the 66 application-provided call-back routines. 67 */ 68 typedef struct { 69 MPI_Comm comm; /* communicator */ 70 DM da; /* distributed array data structure */ 71 Vec localwork; /* local ghosted work vector */ 72 Vec u_local; /* local ghosted approximate solution vector */ 73 Vec solution; /* global exact solution vector */ 74 PetscInt m; /* total number of grid points */ 75 PetscReal h; /* mesh width h = 1/(m-1) */ 76 PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 77 PetscViewer viewer1,viewer2; /* viewers for the solution and error */ 78 PetscReal norm_2,norm_max; /* error norms */ 79 } AppCtx; 80 81 /* 82 User-defined routines 83 */ 84 extern PetscErrorCode InitialConditions(Vec,AppCtx*); 85 extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*); 86 extern PetscErrorCode RHSFunctionHeat(TS,PetscReal,Vec,Vec,void*); 87 extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 88 extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); 89 90 int main(int argc,char **argv) 91 { 92 AppCtx appctx; /* user-defined application context */ 93 TS ts; /* timestepping context */ 94 Mat A; /* matrix data structure */ 95 Vec u; /* approximate solution vector */ 96 PetscReal time_total_max = 1.0; /* default max total time */ 97 PetscInt time_steps_max = 100; /* default max timesteps */ 98 PetscDraw draw; /* drawing context */ 99 PetscErrorCode ierr; 100 PetscInt steps,m; 101 PetscMPIInt size; 102 PetscReal dt,ftime; 103 PetscBool flg; 104 TSProblemType tsproblem = TS_LINEAR; 105 106 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107 Initialize program and set problem parameters 108 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109 110 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 111 appctx.comm = PETSC_COMM_WORLD; 112 113 m = 60; 114 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 115 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug)); 116 appctx.m = m; 117 appctx.h = 1.0/(m-1.0); 118 appctx.norm_2 = 0.0; 119 appctx.norm_max = 0.0; 120 121 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 122 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Solving a linear TS problem, number of processors = %d\n",size)); 123 124 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125 Create vector data structures 126 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127 /* 128 Create distributed array (DMDA) to manage parallel grid and vectors 129 and to set up the ghost point communication pattern. There are M 130 total grid values spread equally among all the processors. 131 */ 132 133 CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,m,1,1,NULL,&appctx.da)); 134 CHKERRQ(DMSetFromOptions(appctx.da)); 135 CHKERRQ(DMSetUp(appctx.da)); 136 137 /* 138 Extract global and local vectors from DMDA; we use these to store the 139 approximate solution. Then duplicate these for remaining vectors that 140 have the same types. 141 */ 142 CHKERRQ(DMCreateGlobalVector(appctx.da,&u)); 143 CHKERRQ(DMCreateLocalVector(appctx.da,&appctx.u_local)); 144 145 /* 146 Create local work vector for use in evaluating right-hand-side function; 147 create global work vector for storing exact solution. 148 */ 149 CHKERRQ(VecDuplicate(appctx.u_local,&appctx.localwork)); 150 CHKERRQ(VecDuplicate(u,&appctx.solution)); 151 152 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153 Set up displays to show graphs of the solution and error 154 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 155 156 CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,380,400,160,&appctx.viewer1)); 157 CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw)); 158 CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 159 CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,0,400,160,&appctx.viewer2)); 160 CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw)); 161 CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 162 163 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 164 Create timestepping solver context 165 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 166 167 CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 168 169 flg = PETSC_FALSE; 170 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-nonlinear",&flg,NULL)); 171 CHKERRQ(TSSetProblemType(ts,flg ? TS_NONLINEAR : TS_LINEAR)); 172 173 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 174 Set optional user-defined monitoring routine 175 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 176 CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL)); 177 178 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 179 180 Create matrix data structure; set matrix evaluation routine. 181 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 182 183 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 184 CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m)); 185 CHKERRQ(MatSetFromOptions(A)); 186 CHKERRQ(MatSetUp(A)); 187 188 flg = PETSC_FALSE; 189 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL)); 190 if (flg) { 191 /* 192 For linear problems with a time-dependent f(u,t) in the equation 193 u_t = f(u,t), the user provides the discretized right-hand-side 194 as a time-dependent matrix. 195 */ 196 CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 197 CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx)); 198 } else { 199 /* 200 For linear problems with a time-independent f(u) in the equation 201 u_t = f(u), the user provides the discretized right-hand-side 202 as a matrix only once, and then sets a null matrix evaluation 203 routine. 204 */ 205 CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx)); 206 CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 207 CHKERRQ(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx)); 208 } 209 210 if (tsproblem == TS_NONLINEAR) { 211 SNES snes; 212 CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunctionHeat,&appctx)); 213 CHKERRQ(TSGetSNES(ts,&snes)); 214 CHKERRQ(SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefault,NULL)); 215 } 216 217 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 218 Set solution vector and initial timestep 219 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 220 221 dt = appctx.h*appctx.h/2.0; 222 CHKERRQ(TSSetTimeStep(ts,dt)); 223 CHKERRQ(TSSetSolution(ts,u)); 224 225 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 226 Customize timestepping solver: 227 - Set the solution method to be the Backward Euler method. 228 - Set timestepping duration info 229 Then set runtime options, which can override these defaults. 230 For example, 231 -ts_max_steps <maxsteps> -ts_max_time <maxtime> 232 to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 233 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 234 235 CHKERRQ(TSSetMaxSteps(ts,time_steps_max)); 236 CHKERRQ(TSSetMaxTime(ts,time_total_max)); 237 CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 238 CHKERRQ(TSSetFromOptions(ts)); 239 240 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241 Solve the problem 242 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 243 244 /* 245 Evaluate initial conditions 246 */ 247 CHKERRQ(InitialConditions(u,&appctx)); 248 249 /* 250 Run the timestepping solver 251 */ 252 CHKERRQ(TSSolve(ts,u)); 253 CHKERRQ(TSGetSolveTime(ts,&ftime)); 254 CHKERRQ(TSGetStepNumber(ts,&steps)); 255 256 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 257 View timestepping solver info 258 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 259 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Total timesteps %D, Final time %g\n",steps,(double)ftime)); 260 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Avg. error (2 norm) = %g Avg. error (max norm) = %g\n",(double)(appctx.norm_2/steps),(double)(appctx.norm_max/steps))); 261 262 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 263 Free work space. All PETSc objects should be destroyed when they 264 are no longer needed. 265 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 266 267 CHKERRQ(TSDestroy(&ts)); 268 CHKERRQ(MatDestroy(&A)); 269 CHKERRQ(VecDestroy(&u)); 270 CHKERRQ(PetscViewerDestroy(&appctx.viewer1)); 271 CHKERRQ(PetscViewerDestroy(&appctx.viewer2)); 272 CHKERRQ(VecDestroy(&appctx.localwork)); 273 CHKERRQ(VecDestroy(&appctx.solution)); 274 CHKERRQ(VecDestroy(&appctx.u_local)); 275 CHKERRQ(DMDestroy(&appctx.da)); 276 277 /* 278 Always call PetscFinalize() before exiting a program. This routine 279 - finalizes the PETSc libraries as well as MPI 280 - provides summary and diagnostic information if certain runtime 281 options are chosen (e.g., -log_view). 282 */ 283 ierr = PetscFinalize(); 284 return ierr; 285 } 286 /* --------------------------------------------------------------------- */ 287 /* 288 InitialConditions - Computes the solution at the initial time. 289 290 Input Parameter: 291 u - uninitialized solution vector (global) 292 appctx - user-defined application context 293 294 Output Parameter: 295 u - vector with solution at initial time (global) 296 */ 297 PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 298 { 299 PetscScalar *u_localptr,h = appctx->h; 300 PetscInt i,mybase,myend; 301 302 /* 303 Determine starting point of each processor's range of 304 grid values. 305 */ 306 CHKERRQ(VecGetOwnershipRange(u,&mybase,&myend)); 307 308 /* 309 Get a pointer to vector data. 310 - For default PETSc vectors, VecGetArray() returns a pointer to 311 the data array. Otherwise, the routine is implementation dependent. 312 - You MUST call VecRestoreArray() when you no longer need access to 313 the array. 314 - Note that the Fortran interface to VecGetArray() differs from the 315 C version. See the users manual for details. 316 */ 317 CHKERRQ(VecGetArray(u,&u_localptr)); 318 319 /* 320 We initialize the solution array by simply writing the solution 321 directly into the array locations. Alternatively, we could use 322 VecSetValues() or VecSetValuesLocal(). 323 */ 324 for (i=mybase; i<myend; i++) u_localptr[i-mybase] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h); 325 326 /* 327 Restore vector 328 */ 329 CHKERRQ(VecRestoreArray(u,&u_localptr)); 330 331 /* 332 Print debugging information if desired 333 */ 334 if (appctx->debug) { 335 CHKERRQ(PetscPrintf(appctx->comm,"initial guess vector\n")); 336 CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 337 } 338 339 return 0; 340 } 341 /* --------------------------------------------------------------------- */ 342 /* 343 ExactSolution - Computes the exact solution at a given time. 344 345 Input Parameters: 346 t - current time 347 solution - vector in which exact solution will be computed 348 appctx - user-defined application context 349 350 Output Parameter: 351 solution - vector with the newly computed exact solution 352 */ 353 PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 354 { 355 PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2; 356 PetscInt i,mybase,myend; 357 358 /* 359 Determine starting and ending points of each processor's 360 range of grid values 361 */ 362 CHKERRQ(VecGetOwnershipRange(solution,&mybase,&myend)); 363 364 /* 365 Get a pointer to vector data. 366 */ 367 CHKERRQ(VecGetArray(solution,&s_localptr)); 368 369 /* 370 Simply write the solution directly into the array locations. 371 Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 372 */ 373 ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t); ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t); 374 sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 375 for (i=mybase; i<myend; i++) s_localptr[i-mybase] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2; 376 377 /* 378 Restore vector 379 */ 380 CHKERRQ(VecRestoreArray(solution,&s_localptr)); 381 return 0; 382 } 383 /* --------------------------------------------------------------------- */ 384 /* 385 Monitor - User-provided routine to monitor the solution computed at 386 each timestep. This example plots the solution and computes the 387 error in two different norms. 388 389 Input Parameters: 390 ts - the timestep context 391 step - the count of the current step (with 0 meaning the 392 initial condition) 393 time - the current time 394 u - the solution at this timestep 395 ctx - the user-provided context for this monitoring routine. 396 In this case we use the application context which contains 397 information about the problem size, workspace and the exact 398 solution. 399 */ 400 PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 401 { 402 AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 403 PetscReal norm_2,norm_max; 404 405 /* 406 View a graph of the current iterate 407 */ 408 CHKERRQ(VecView(u,appctx->viewer2)); 409 410 /* 411 Compute the exact solution 412 */ 413 CHKERRQ(ExactSolution(time,appctx->solution,appctx)); 414 415 /* 416 Print debugging information if desired 417 */ 418 if (appctx->debug) { 419 CHKERRQ(PetscPrintf(appctx->comm,"Computed solution vector\n")); 420 CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 421 CHKERRQ(PetscPrintf(appctx->comm,"Exact solution vector\n")); 422 CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD)); 423 } 424 425 /* 426 Compute the 2-norm and max-norm of the error 427 */ 428 CHKERRQ(VecAXPY(appctx->solution,-1.0,u)); 429 CHKERRQ(VecNorm(appctx->solution,NORM_2,&norm_2)); 430 norm_2 = PetscSqrtReal(appctx->h)*norm_2; 431 CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&norm_max)); 432 if (norm_2 < 1e-14) norm_2 = 0; 433 if (norm_max < 1e-14) norm_max = 0; 434 435 /* 436 PetscPrintf() causes only the first processor in this 437 communicator to print the timestep information. 438 */ 439 CHKERRQ(PetscPrintf(appctx->comm,"Timestep %D: time = %g 2-norm error = %g max norm error = %g\n",step,(double)time,(double)norm_2,(double)norm_max)); 440 appctx->norm_2 += norm_2; 441 appctx->norm_max += norm_max; 442 443 /* 444 View a graph of the error 445 */ 446 CHKERRQ(VecView(appctx->solution,appctx->viewer1)); 447 448 /* 449 Print debugging information if desired 450 */ 451 if (appctx->debug) { 452 CHKERRQ(PetscPrintf(appctx->comm,"Error vector\n")); 453 CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD)); 454 } 455 456 return 0; 457 } 458 459 /* --------------------------------------------------------------------- */ 460 /* 461 RHSMatrixHeat - User-provided routine to compute the right-hand-side 462 matrix for the heat equation. 463 464 Input Parameters: 465 ts - the TS context 466 t - current time 467 global_in - global input vector 468 dummy - optional user-defined context, as set by TSetRHSJacobian() 469 470 Output Parameters: 471 AA - Jacobian matrix 472 BB - optionally different preconditioning matrix 473 str - flag indicating matrix structure 474 475 Notes: 476 RHSMatrixHeat computes entries for the locally owned part of the system. 477 - Currently, all PETSc parallel matrix formats are partitioned by 478 contiguous chunks of rows across the processors. 479 - Each processor needs to insert only elements that it owns 480 locally (but any non-local elements will be sent to the 481 appropriate processor during matrix assembly). 482 - Always specify global row and columns of matrix entries when 483 using MatSetValues(); we could alternatively use MatSetValuesLocal(). 484 - Here, we set all entries for a particular row at once. 485 - Note that MatSetValues() uses 0-based row and column numbers 486 in Fortran as well as in C. 487 */ 488 PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx) 489 { 490 Mat A = AA; /* Jacobian matrix */ 491 AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 492 PetscInt i,mstart,mend,idx[3]; 493 PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; 494 495 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 496 Compute entries for the locally owned part of the matrix 497 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 498 499 CHKERRQ(MatGetOwnershipRange(A,&mstart,&mend)); 500 501 /* 502 Set matrix rows corresponding to boundary data 503 */ 504 505 if (mstart == 0) { /* first processor only */ 506 v[0] = 1.0; 507 CHKERRQ(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES)); 508 mstart++; 509 } 510 511 if (mend == appctx->m) { /* last processor only */ 512 mend--; 513 v[0] = 1.0; 514 CHKERRQ(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES)); 515 } 516 517 /* 518 Set matrix rows corresponding to interior data. We construct the 519 matrix one row at a time. 520 */ 521 v[0] = sone; v[1] = stwo; v[2] = sone; 522 for (i=mstart; i<mend; i++) { 523 idx[0] = i-1; idx[1] = i; idx[2] = i+1; 524 CHKERRQ(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES)); 525 } 526 527 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 528 Complete the matrix assembly process and set some options 529 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 530 /* 531 Assemble matrix, using the 2-step process: 532 MatAssemblyBegin(), MatAssemblyEnd() 533 Computations can be done while messages are in transition 534 by placing code between these two statements. 535 */ 536 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 537 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 538 539 /* 540 Set and option to indicate that we will never add a new nonzero location 541 to the matrix. If we do, it will generate an error. 542 */ 543 CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 544 545 return 0; 546 } 547 548 PetscErrorCode RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 549 { 550 Mat A; 551 552 PetscFunctionBeginUser; 553 CHKERRQ(TSGetRHSJacobian(ts,&A,NULL,NULL,&ctx)); 554 CHKERRQ(RHSMatrixHeat(ts,t,globalin,A,NULL,ctx)); 555 /* CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); */ 556 CHKERRQ(MatMult(A,globalin,globalout)); 557 PetscFunctionReturn(0); 558 } 559 560 /*TEST 561 562 test: 563 args: -ts_view -nox 564 565 test: 566 suffix: 2 567 args: -ts_view -nox 568 nsize: 3 569 570 test: 571 suffix: 3 572 args: -ts_view -nox -nonlinear 573 574 test: 575 suffix: 4 576 args: -ts_view -nox -nonlinear 577 nsize: 3 578 timeoutfactor: 3 579 580 test: 581 suffix: sundials 582 requires: sundials2 583 args: -nox -ts_type sundials -ts_max_steps 5 -nonlinear 584 nsize: 4 585 586 test: 587 suffix: sundials_dense 588 requires: sundials2 589 args: -nox -ts_type sundials -ts_sundials_use_dense -ts_max_steps 5 -nonlinear 590 nsize: 1 591 592 TEST*/ 593