1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves a simple time-dependent linear PDE (the heat equation).\n\ 3c4762a1bSJed Brown Input parameters include:\n\ 4c4762a1bSJed Brown -m <points>, where <points> = number of grid points\n\ 5c4762a1bSJed Brown -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\ 6c4762a1bSJed Brown -debug : Activate debugging printouts\n\ 7c4762a1bSJed Brown -nox : Deactivate x-window graphics\n\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* ------------------------------------------------------------------------ 10c4762a1bSJed Brown 11c4762a1bSJed Brown This program solves the one-dimensional heat equation (also called the 12c4762a1bSJed Brown diffusion equation), 13c4762a1bSJed Brown u_t = u_xx, 14c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 15c4762a1bSJed Brown u(t,0) = 1, u(t,1) = 1, 16c4762a1bSJed Brown and the initial condition 17c4762a1bSJed Brown u(0,x) = cos(6*pi*x) + 3*cos(2*pi*x). 18c4762a1bSJed Brown This is a linear, second-order, parabolic equation. 19c4762a1bSJed Brown 20c4762a1bSJed Brown We discretize the right-hand side using finite differences with 21c4762a1bSJed Brown uniform grid spacing h: 22c4762a1bSJed Brown u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) 23c4762a1bSJed Brown We then demonstrate time evolution using the various TS methods by 24c4762a1bSJed Brown running the program via 25c4762a1bSJed Brown ex3 -ts_type <timestepping solver> 26c4762a1bSJed Brown 27c4762a1bSJed Brown We compare the approximate solution with the exact solution, given by 28c4762a1bSJed Brown u_exact(x,t) = exp(-36*pi*pi*t) * cos(6*pi*x) + 29c4762a1bSJed Brown 3*exp(-4*pi*pi*t) * cos(2*pi*x) 30c4762a1bSJed Brown 31c4762a1bSJed Brown Notes: 32c4762a1bSJed Brown This code demonstrates the TS solver interface to two variants of 33c4762a1bSJed Brown linear problems, u_t = f(u,t), namely 34c4762a1bSJed Brown - time-dependent f: f(u,t) is a function of t 35c4762a1bSJed Brown - time-independent f: f(u,t) is simply just f(u) 36c4762a1bSJed Brown 37c4762a1bSJed Brown The parallel version of this code is ts/tutorials/ex4.c 38c4762a1bSJed Brown 39c4762a1bSJed Brown ------------------------------------------------------------------------- */ 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 43c4762a1bSJed Brown automatically includes: 44c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 45c4762a1bSJed Brown petscmat.h - matrices 46c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 47c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 48c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown #include <petscts.h> 51c4762a1bSJed Brown #include <petscdraw.h> 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* 54c4762a1bSJed Brown User-defined application context - contains data needed by the 55c4762a1bSJed Brown application-provided call-back routines. 56c4762a1bSJed Brown */ 57c4762a1bSJed Brown typedef struct { 58c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 59c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 60c4762a1bSJed Brown PetscReal h; /* mesh width h = 1/(m-1) */ 61c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 62c4762a1bSJed Brown PetscViewer viewer1, viewer2; /* viewers for the solution and error */ 63c4762a1bSJed Brown PetscReal norm_2, norm_max; /* error norms */ 64c4762a1bSJed Brown } AppCtx; 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown User-defined routines 68c4762a1bSJed Brown */ 69c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *); 70c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *); 71c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *); 72c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *); 73c4762a1bSJed Brown 74*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 75*d71ae5a4SJacob Faibussowitsch { 76c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 77c4762a1bSJed Brown TS ts; /* timestepping context */ 78c4762a1bSJed Brown Mat A; /* matrix data structure */ 79c4762a1bSJed Brown Vec u; /* approximate solution vector */ 80c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 81c4762a1bSJed Brown PetscInt time_steps_max = 100; /* default max timesteps */ 82c4762a1bSJed Brown PetscDraw draw; /* drawing context */ 83c4762a1bSJed Brown PetscInt steps, m; 84c4762a1bSJed Brown PetscMPIInt size; 85c4762a1bSJed Brown PetscBool flg; 86c4762a1bSJed Brown PetscReal dt, ftime; 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89c4762a1bSJed Brown Initialize program and set problem parameters 90c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 91c4762a1bSJed Brown 92327415f7SBarry Smith PetscFunctionBeginUser; 939566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 953c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 96c4762a1bSJed Brown 97c4762a1bSJed Brown m = 60; 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug)); 100c4762a1bSJed Brown appctx.m = m; 101c4762a1bSJed Brown appctx.h = 1.0 / (m - 1.0); 102c4762a1bSJed Brown appctx.norm_2 = 0.0; 103c4762a1bSJed Brown appctx.norm_max = 0.0; 104c4762a1bSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving a linear TS problem on 1 processor\n")); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Create vector data structures 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* 112c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 113c4762a1bSJed Brown */ 1149566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &u)); 1159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.solution)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Set up displays to show graphs of the solution and error 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 380, 400, 160, &appctx.viewer1)); 1229566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw)); 1239566063dSJacob Faibussowitsch PetscCall(PetscDrawSetDoubleBuffer(draw)); 1249566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, "", 80, 0, 400, 160, &appctx.viewer2)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw)); 1269566063dSJacob Faibussowitsch PetscCall(PetscDrawSetDoubleBuffer(draw)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Create timestepping solver context 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131c4762a1bSJed Brown 1329566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF, &ts)); 1339566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_LINEAR)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Set optional user-defined monitoring routine 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142c4762a1bSJed Brown 143c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 1479566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m)); 1489566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1499566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-time_dependent_rhs", &flg)); 152c4762a1bSJed Brown if (flg) { 153c4762a1bSJed Brown /* 154c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 155c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 156c4762a1bSJed Brown as a time-dependent matrix. 157c4762a1bSJed Brown */ 1589566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx)); 1599566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx)); 160c4762a1bSJed Brown } else { 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 163c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 164c4762a1bSJed Brown as a matrix only once, and then sets a null matrix evaluation 165c4762a1bSJed Brown routine. 166c4762a1bSJed Brown */ 1679566063dSJacob Faibussowitsch PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx)); 1689566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx)); 1699566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx)); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173c4762a1bSJed Brown Set solution vector and initial timestep 174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 175c4762a1bSJed Brown 176c4762a1bSJed Brown dt = appctx.h * appctx.h / 2.0; 1779566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 1789566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, u)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Customize timestepping solver: 182c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 183c4762a1bSJed Brown - Set timestepping duration info 184c4762a1bSJed Brown Then set runtime options, which can override these defaults. 185c4762a1bSJed Brown For example, 186c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 187c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, time_steps_max)); 1919566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, time_total_max)); 1929566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1939566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196c4762a1bSJed Brown Solve the problem 197c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* 200c4762a1bSJed Brown Evaluate initial conditions 201c4762a1bSJed Brown */ 2029566063dSJacob Faibussowitsch PetscCall(InitialConditions(u, &appctx)); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* 205c4762a1bSJed Brown Run the timestepping solver 206c4762a1bSJed Brown */ 2079566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 2089566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2099566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212c4762a1bSJed Brown View timestepping solver info 213c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 214c4762a1bSJed Brown 2159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "avg. error (2 norm) = %g, avg. error (max norm) = %g\n", (double)(appctx.norm_2 / steps), (double)(appctx.norm_max / steps))); 2169566063dSJacob Faibussowitsch PetscCall(TSView(ts, PETSC_VIEWER_STDOUT_SELF)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 219c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 220c4762a1bSJed Brown are no longer needed. 221c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2269566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&appctx.viewer1)); 2279566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&appctx.viewer2)); 2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.solution)); 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* 231c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 232c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 233c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 234c4762a1bSJed Brown options are chosen (e.g., -log_view). 235c4762a1bSJed Brown */ 2369566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 237b122ec5aSJacob Faibussowitsch return 0; 238c4762a1bSJed Brown } 239c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 240c4762a1bSJed Brown /* 241c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 242c4762a1bSJed Brown 243c4762a1bSJed Brown Input Parameter: 244c4762a1bSJed Brown u - uninitialized solution vector (global) 245c4762a1bSJed Brown appctx - user-defined application context 246c4762a1bSJed Brown 247c4762a1bSJed Brown Output Parameter: 248c4762a1bSJed Brown u - vector with solution at initial time (global) 249c4762a1bSJed Brown */ 250*d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx) 251*d71ae5a4SJacob Faibussowitsch { 252c4762a1bSJed Brown PetscScalar *u_localptr, h = appctx->h; 253c4762a1bSJed Brown PetscInt i; 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* 256c4762a1bSJed Brown Get a pointer to vector data. 257c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 258c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 259c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 260c4762a1bSJed Brown the array. 261c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 262c4762a1bSJed Brown C version. See the users manual for details. 263c4762a1bSJed Brown */ 2649566063dSJacob Faibussowitsch PetscCall(VecGetArray(u, &u_localptr)); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* 267c4762a1bSJed Brown We initialize the solution array by simply writing the solution 268c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 269c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 270c4762a1bSJed Brown */ 271c4762a1bSJed Brown for (i = 0; i < appctx->m; i++) u_localptr[i] = PetscCosScalar(PETSC_PI * i * 6. * h) + 3. * PetscCosScalar(PETSC_PI * i * 2. * h); 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* 274c4762a1bSJed Brown Restore vector 275c4762a1bSJed Brown */ 2769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(u, &u_localptr)); 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* 279c4762a1bSJed Brown Print debugging information if desired 280c4762a1bSJed Brown */ 281c4762a1bSJed Brown if (appctx->debug) { 282c4762a1bSJed Brown printf("initial guess vector\n"); 2839566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF)); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown return 0; 287c4762a1bSJed Brown } 288c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 289c4762a1bSJed Brown /* 290c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 291c4762a1bSJed Brown 292c4762a1bSJed Brown Input Parameters: 293c4762a1bSJed Brown t - current time 294c4762a1bSJed Brown solution - vector in which exact solution will be computed 295c4762a1bSJed Brown appctx - user-defined application context 296c4762a1bSJed Brown 297c4762a1bSJed Brown Output Parameter: 298c4762a1bSJed Brown solution - vector with the newly computed exact solution 299c4762a1bSJed Brown */ 300*d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx) 301*d71ae5a4SJacob Faibussowitsch { 302c4762a1bSJed Brown PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2, tc = t; 303c4762a1bSJed Brown PetscInt i; 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* 306c4762a1bSJed Brown Get a pointer to vector data. 307c4762a1bSJed Brown */ 3089566063dSJacob Faibussowitsch PetscCall(VecGetArray(solution, &s_localptr)); 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* 311c4762a1bSJed Brown Simply write the solution directly into the array locations. 312c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 313c4762a1bSJed Brown */ 3149371c9d4SSatish Balay ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * tc); 3159371c9d4SSatish Balay ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * tc); 3169371c9d4SSatish Balay sc1 = PETSC_PI * 6. * h; 3179371c9d4SSatish Balay sc2 = PETSC_PI * 2. * h; 318c4762a1bSJed Brown for (i = 0; i < appctx->m; i++) s_localptr[i] = PetscCosScalar(sc1 * (PetscReal)i) * ex1 + 3. * PetscCosScalar(sc2 * (PetscReal)i) * ex2; 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown Restore vector 322c4762a1bSJed Brown */ 3239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(solution, &s_localptr)); 324c4762a1bSJed Brown return 0; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 327c4762a1bSJed Brown /* 328c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 329c4762a1bSJed Brown each timestep. This example plots the solution and computes the 330c4762a1bSJed Brown error in two different norms. 331c4762a1bSJed Brown 332c4762a1bSJed Brown Input Parameters: 333c4762a1bSJed Brown ts - the timestep context 334c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 335c4762a1bSJed Brown initial condition) 336c4762a1bSJed Brown time - the current time 337c4762a1bSJed Brown u - the solution at this timestep 338c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 339c4762a1bSJed Brown In this case we use the application context which contains 340c4762a1bSJed Brown information about the problem size, workspace and the exact 341c4762a1bSJed Brown solution. 342c4762a1bSJed Brown */ 343*d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) 344*d71ae5a4SJacob Faibussowitsch { 345c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 346c4762a1bSJed Brown PetscReal norm_2, norm_max; 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* 349c4762a1bSJed Brown View a graph of the current iterate 350c4762a1bSJed Brown */ 3519566063dSJacob Faibussowitsch PetscCall(VecView(u, appctx->viewer2)); 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* 354c4762a1bSJed Brown Compute the exact solution 355c4762a1bSJed Brown */ 3569566063dSJacob Faibussowitsch PetscCall(ExactSolution(time, appctx->solution, appctx)); 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* 359c4762a1bSJed Brown Print debugging information if desired 360c4762a1bSJed Brown */ 361c4762a1bSJed Brown if (appctx->debug) { 362c4762a1bSJed Brown printf("Computed solution vector\n"); 3639566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF)); 364c4762a1bSJed Brown printf("Exact solution vector\n"); 3659566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF)); 366c4762a1bSJed Brown } 367c4762a1bSJed Brown 368c4762a1bSJed Brown /* 369c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 370c4762a1bSJed Brown */ 3719566063dSJacob Faibussowitsch PetscCall(VecAXPY(appctx->solution, -1.0, u)); 3729566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2)); 373c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h) * norm_2; 3749566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max)); 375c4762a1bSJed Brown if (norm_2 < 1e-14) norm_2 = 0; 376c4762a1bSJed Brown if (norm_max < 1e-14) norm_max = 0; 377c4762a1bSJed Brown 37863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep %" PetscInt_FMT ": time = %g, 2-norm error = %g, max norm error = %g\n", step, (double)time, (double)norm_2, (double)norm_max)); 379c4762a1bSJed Brown appctx->norm_2 += norm_2; 380c4762a1bSJed Brown appctx->norm_max += norm_max; 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* 383c4762a1bSJed Brown View a graph of the error 384c4762a1bSJed Brown */ 3859566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution, appctx->viewer1)); 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* 388c4762a1bSJed Brown Print debugging information if desired 389c4762a1bSJed Brown */ 390c4762a1bSJed Brown if (appctx->debug) { 391c4762a1bSJed Brown printf("Error vector\n"); 3929566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF)); 393c4762a1bSJed Brown } 394c4762a1bSJed Brown 395c4762a1bSJed Brown return 0; 396c4762a1bSJed Brown } 397c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 398c4762a1bSJed Brown /* 399c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 400c4762a1bSJed Brown matrix for the heat equation. 401c4762a1bSJed Brown 402c4762a1bSJed Brown Input Parameters: 403c4762a1bSJed Brown ts - the TS context 404c4762a1bSJed Brown t - current time 405c4762a1bSJed Brown global_in - global input vector 406c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 407c4762a1bSJed Brown 408c4762a1bSJed Brown Output Parameters: 409c4762a1bSJed Brown AA - Jacobian matrix 410c4762a1bSJed Brown BB - optionally different preconditioning matrix 411c4762a1bSJed Brown str - flag indicating matrix structure 412c4762a1bSJed Brown 413c4762a1bSJed Brown Notes: 414c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 415c4762a1bSJed Brown in Fortran as well as in C. 416c4762a1bSJed Brown */ 417*d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, void *ctx) 418*d71ae5a4SJacob Faibussowitsch { 419c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 420c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 421c4762a1bSJed Brown PetscInt mstart = 0; 422c4762a1bSJed Brown PetscInt mend = appctx->m; 423c4762a1bSJed Brown PetscInt i, idx[3]; 424c4762a1bSJed Brown PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo; 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 427c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 428c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 429c4762a1bSJed Brown /* 430c4762a1bSJed Brown Set matrix rows corresponding to boundary data 431c4762a1bSJed Brown */ 432c4762a1bSJed Brown 433c4762a1bSJed Brown mstart = 0; 434c4762a1bSJed Brown v[0] = 1.0; 4359566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES)); 436c4762a1bSJed Brown mstart++; 437c4762a1bSJed Brown 438c4762a1bSJed Brown mend--; 439c4762a1bSJed Brown v[0] = 1.0; 4409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES)); 441c4762a1bSJed Brown 442c4762a1bSJed Brown /* 443c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 444c4762a1bSJed Brown matrix one row at a time. 445c4762a1bSJed Brown */ 4469371c9d4SSatish Balay v[0] = sone; 4479371c9d4SSatish Balay v[1] = stwo; 4489371c9d4SSatish Balay v[2] = sone; 449c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 4509371c9d4SSatish Balay idx[0] = i - 1; 4519371c9d4SSatish Balay idx[1] = i; 4529371c9d4SSatish Balay idx[2] = i + 1; 4539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES)); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 457c4762a1bSJed Brown Complete the matrix assembly process and set some options 458c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 459c4762a1bSJed Brown /* 460c4762a1bSJed Brown Assemble matrix, using the 2-step process: 461c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 462c4762a1bSJed Brown Computations can be done while messages are in transition 463c4762a1bSJed Brown by placing code between these two statements. 464c4762a1bSJed Brown */ 4659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 4669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 467c4762a1bSJed Brown 468c4762a1bSJed Brown /* 469c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 470c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 471c4762a1bSJed Brown */ 4729566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 473c4762a1bSJed Brown 474c4762a1bSJed Brown return 0; 475c4762a1bSJed Brown } 476c4762a1bSJed Brown 477c4762a1bSJed Brown /*TEST 478c4762a1bSJed Brown 479c4762a1bSJed Brown test: 480c4762a1bSJed Brown requires: x 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: nox 484c4762a1bSJed Brown args: -nox 485c4762a1bSJed Brown 486c4762a1bSJed Brown TEST*/ 487