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) = 0, u(t,1) = 0, 16c4762a1bSJed Brown and the initial condition 17c4762a1bSJed Brown u(0,x) = sin(6*pi*x) + 3*sin(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 mpiexec -n <procs> 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) * sin(6*pi*x) + 29c4762a1bSJed Brown 3*exp(-4*pi*pi*t) * sin(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 f(u) 36c4762a1bSJed Brown 37c4762a1bSJed Brown The uniprocessor version of this code is ts/tutorials/ex3.c 38c4762a1bSJed Brown 39c4762a1bSJed Brown ------------------------------------------------------------------------- */ 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs) to manage 43c4762a1bSJed Brown the parallel grid. Include "petscts.h" so that we can use TS solvers. 44c4762a1bSJed Brown Note that this file automatically includes: 45c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 46c4762a1bSJed Brown petscmat.h - matrices 47c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 48c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 49c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 50c4762a1bSJed Brown */ 51c4762a1bSJed Brown 52c4762a1bSJed Brown #include <petscdm.h> 53c4762a1bSJed Brown #include <petscdmda.h> 54c4762a1bSJed Brown #include <petscts.h> 55c4762a1bSJed Brown #include <petscdraw.h> 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* 58c4762a1bSJed Brown User-defined application context - contains data needed by the 59c4762a1bSJed Brown application-provided call-back routines. 60c4762a1bSJed Brown */ 61c4762a1bSJed Brown typedef struct { 62c4762a1bSJed Brown MPI_Comm comm; /* communicator */ 63c4762a1bSJed Brown DM da; /* distributed array data structure */ 64c4762a1bSJed Brown Vec localwork; /* local ghosted work vector */ 65c4762a1bSJed Brown Vec u_local; /* local ghosted approximate solution vector */ 66c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 67c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 68c4762a1bSJed Brown PetscReal h; /* mesh width h = 1/(m-1) */ 69c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 70c4762a1bSJed Brown PetscViewer viewer1, viewer2; /* viewers for the solution and error */ 71c4762a1bSJed Brown PetscReal norm_2, norm_max; /* error norms */ 72c4762a1bSJed Brown } AppCtx; 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* 75c4762a1bSJed Brown User-defined routines 76c4762a1bSJed Brown */ 77c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *); 78c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS, PetscReal, Vec, Mat, Mat, void *); 79c4762a1bSJed Brown extern PetscErrorCode RHSFunctionHeat(TS, PetscReal, Vec, Vec, void *); 80c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *); 81c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *); 82c4762a1bSJed Brown 83d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 84d71ae5a4SJacob Faibussowitsch { 85c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 86c4762a1bSJed Brown TS ts; /* timestepping context */ 87c4762a1bSJed Brown Mat A; /* matrix data structure */ 88c4762a1bSJed Brown Vec u; /* approximate solution vector */ 89c4762a1bSJed Brown PetscReal time_total_max = 1.0; /* default max total time */ 90c4762a1bSJed Brown PetscInt time_steps_max = 100; /* default max timesteps */ 91c4762a1bSJed Brown PetscDraw draw; /* drawing context */ 92c4762a1bSJed Brown PetscInt steps, m; 93c4762a1bSJed Brown PetscMPIInt size; 94c4762a1bSJed Brown PetscReal dt, ftime; 95c4762a1bSJed Brown PetscBool flg; 96c4762a1bSJed Brown TSProblemType tsproblem = TS_LINEAR; 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99c4762a1bSJed Brown Initialize program and set problem parameters 100c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 101c4762a1bSJed Brown 102327415f7SBarry Smith PetscFunctionBeginUser; 1039566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 104c4762a1bSJed Brown appctx.comm = PETSC_COMM_WORLD; 105c4762a1bSJed Brown 106c4762a1bSJed Brown m = 60; 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug)); 109c4762a1bSJed Brown appctx.m = m; 110c4762a1bSJed Brown appctx.h = 1.0 / (m - 1.0); 111c4762a1bSJed Brown appctx.norm_2 = 0.0; 112c4762a1bSJed Brown appctx.norm_max = 0.0; 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solving a linear TS problem, number of processors = %d\n", size)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Create vector data structures 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 122c4762a1bSJed Brown and to set up the ghost point communication pattern. There are M 123c4762a1bSJed Brown total grid values spread equally among all the processors. 124c4762a1bSJed Brown */ 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, m, 1, 1, NULL, &appctx.da)); 1279566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(appctx.da)); 1289566063dSJacob Faibussowitsch PetscCall(DMSetUp(appctx.da)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* 131c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 132c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 133c4762a1bSJed Brown have the same types. 134c4762a1bSJed Brown */ 1359566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(appctx.da, &u)); 1369566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(appctx.da, &appctx.u_local)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown Create local work vector for use in evaluating right-hand-side function; 140c4762a1bSJed Brown create global work vector for storing exact solution. 141c4762a1bSJed Brown */ 1429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.u_local, &appctx.localwork)); 1439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.solution)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Set up displays to show graphs of the solution and error 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148c4762a1bSJed Brown 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "", 80, 380, 400, 160, &appctx.viewer1)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(appctx.viewer1, 0, &draw)); 1519566063dSJacob Faibussowitsch PetscCall(PetscDrawSetDoubleBuffer(draw)); 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, "", 80, 0, 400, 160, &appctx.viewer2)); 1539566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(appctx.viewer2, 0, &draw)); 1549566063dSJacob Faibussowitsch PetscCall(PetscDrawSetDoubleBuffer(draw)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157c4762a1bSJed Brown Create timestepping solver context 158c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 159c4762a1bSJed Brown 1609566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown flg = PETSC_FALSE; 1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonlinear", &flg, NULL)); 1649566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, flg ? TS_NONLINEAR : TS_LINEAR)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 167c4762a1bSJed Brown Set optional user-defined monitoring routine 168c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1699566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL)); 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172c4762a1bSJed Brown 173c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1779566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, m)); 1789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1799566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown flg = PETSC_FALSE; 1829566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-time_dependent_rhs", &flg, NULL)); 183c4762a1bSJed Brown if (flg) { 184c4762a1bSJed Brown /* 185c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 186c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 187c4762a1bSJed Brown as a time-dependent matrix. 188c4762a1bSJed Brown */ 1899566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx)); 1909566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSMatrixHeat, &appctx)); 191c4762a1bSJed Brown } else { 192c4762a1bSJed Brown /* 193c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 194c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 195c4762a1bSJed Brown as a matrix only once, and then sets a null matrix evaluation 196c4762a1bSJed Brown routine. 197c4762a1bSJed Brown */ 1989566063dSJacob Faibussowitsch PetscCall(RHSMatrixHeat(ts, 0.0, u, A, A, &appctx)); 1999566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx)); 2009566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, &appctx)); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown if (tsproblem == TS_NONLINEAR) { 204c4762a1bSJed Brown SNES snes; 2059566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionHeat, &appctx)); 2069566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2079566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESComputeJacobianDefault, NULL)); 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211c4762a1bSJed Brown Set solution vector and initial timestep 212c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 213c4762a1bSJed Brown 214c4762a1bSJed Brown dt = appctx.h * appctx.h / 2.0; 2159566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 2169566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, u)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 219c4762a1bSJed Brown Customize timestepping solver: 220c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 221c4762a1bSJed Brown - Set timestepping duration info 222c4762a1bSJed Brown Then set runtime options, which can override these defaults. 223c4762a1bSJed Brown For example, 224c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 225c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 226c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 227c4762a1bSJed Brown 2289566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, time_steps_max)); 2299566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, time_total_max)); 2309566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 2319566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 234c4762a1bSJed Brown Solve the problem 235c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown Evaluate initial conditions 239c4762a1bSJed Brown */ 2409566063dSJacob Faibussowitsch PetscCall(InitialConditions(u, &appctx)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* 243c4762a1bSJed Brown Run the timestepping solver 244c4762a1bSJed Brown */ 2459566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 2469566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2479566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 250c4762a1bSJed Brown View timestepping solver info 251c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 25263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total timesteps %" PetscInt_FMT ", Final time %g\n", steps, (double)ftime)); 2539566063dSJacob Faibussowitsch PetscCall(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))); 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 256c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 257c4762a1bSJed Brown are no longer needed. 258c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2639566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&appctx.viewer1)); 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&appctx.viewer2)); 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.localwork)); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.solution)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.u_local)); 2689566063dSJacob Faibussowitsch PetscCall(DMDestroy(&appctx.da)); 269c4762a1bSJed Brown 270c4762a1bSJed Brown /* 271c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 272c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 273c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 274c4762a1bSJed Brown options are chosen (e.g., -log_view). 275c4762a1bSJed Brown */ 2769566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 277b122ec5aSJacob Faibussowitsch return 0; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 280c4762a1bSJed Brown /* 281c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 282c4762a1bSJed Brown 283c4762a1bSJed Brown Input Parameter: 284c4762a1bSJed Brown u - uninitialized solution vector (global) 285c4762a1bSJed Brown appctx - user-defined application context 286c4762a1bSJed Brown 287c4762a1bSJed Brown Output Parameter: 288c4762a1bSJed Brown u - vector with solution at initial time (global) 289c4762a1bSJed Brown */ 290d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx) 291d71ae5a4SJacob Faibussowitsch { 292c4762a1bSJed Brown PetscScalar *u_localptr, h = appctx->h; 293c4762a1bSJed Brown PetscInt i, mybase, myend; 294c4762a1bSJed Brown 295*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 296c4762a1bSJed Brown /* 297c4762a1bSJed Brown Determine starting point of each processor's range of 298c4762a1bSJed Brown grid values. 299c4762a1bSJed Brown */ 3009566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(u, &mybase, &myend)); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* 303c4762a1bSJed Brown Get a pointer to vector data. 304c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 305c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 306c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 307c4762a1bSJed Brown the array. 308c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 309c4762a1bSJed Brown C version. See the users manual for details. 310c4762a1bSJed Brown */ 3119566063dSJacob Faibussowitsch PetscCall(VecGetArray(u, &u_localptr)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown We initialize the solution array by simply writing the solution 315c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 316c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 317c4762a1bSJed Brown */ 318c4762a1bSJed Brown for (i = mybase; i < myend; i++) u_localptr[i - mybase] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h); 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown Restore vector 322c4762a1bSJed Brown */ 3239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(u, &u_localptr)); 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* 326c4762a1bSJed Brown Print debugging information if desired 327c4762a1bSJed Brown */ 328c4762a1bSJed Brown if (appctx->debug) { 3299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(appctx->comm, "initial guess vector\n")); 3309566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 336c4762a1bSJed Brown /* 337c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 338c4762a1bSJed Brown 339c4762a1bSJed Brown Input Parameters: 340c4762a1bSJed Brown t - current time 341c4762a1bSJed Brown solution - vector in which exact solution will be computed 342c4762a1bSJed Brown appctx - user-defined application context 343c4762a1bSJed Brown 344c4762a1bSJed Brown Output Parameter: 345c4762a1bSJed Brown solution - vector with the newly computed exact solution 346c4762a1bSJed Brown */ 347d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx) 348d71ae5a4SJacob Faibussowitsch { 349c4762a1bSJed Brown PetscScalar *s_localptr, h = appctx->h, ex1, ex2, sc1, sc2; 350c4762a1bSJed Brown PetscInt i, mybase, myend; 351c4762a1bSJed Brown 352*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 353c4762a1bSJed Brown /* 354c4762a1bSJed Brown Determine starting and ending points of each processor's 355c4762a1bSJed Brown range of grid values 356c4762a1bSJed Brown */ 3579566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(solution, &mybase, &myend)); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* 360c4762a1bSJed Brown Get a pointer to vector data. 361c4762a1bSJed Brown */ 3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(solution, &s_localptr)); 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* 365c4762a1bSJed Brown Simply write the solution directly into the array locations. 366c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 367c4762a1bSJed Brown */ 3689371c9d4SSatish Balay ex1 = PetscExpReal(-36. * PETSC_PI * PETSC_PI * t); 3699371c9d4SSatish Balay ex2 = PetscExpReal(-4. * PETSC_PI * PETSC_PI * t); 3709371c9d4SSatish Balay sc1 = PETSC_PI * 6. * h; 3719371c9d4SSatish Balay sc2 = PETSC_PI * 2. * h; 372c4762a1bSJed Brown for (i = mybase; i < myend; i++) s_localptr[i - mybase] = PetscSinScalar(sc1 * (PetscReal)i) * ex1 + 3. * PetscSinScalar(sc2 * (PetscReal)i) * ex2; 373c4762a1bSJed Brown 374c4762a1bSJed Brown /* 375c4762a1bSJed Brown Restore vector 376c4762a1bSJed Brown */ 3779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(solution, &s_localptr)); 378*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 379c4762a1bSJed Brown } 380c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 381c4762a1bSJed Brown /* 382c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 383c4762a1bSJed Brown each timestep. This example plots the solution and computes the 384c4762a1bSJed Brown error in two different norms. 385c4762a1bSJed Brown 386c4762a1bSJed Brown Input Parameters: 387c4762a1bSJed Brown ts - the timestep context 388c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 389c4762a1bSJed Brown initial condition) 390c4762a1bSJed Brown time - the current time 391c4762a1bSJed Brown u - the solution at this timestep 392c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 393c4762a1bSJed Brown In this case we use the application context which contains 394c4762a1bSJed Brown information about the problem size, workspace and the exact 395c4762a1bSJed Brown solution. 396c4762a1bSJed Brown */ 397d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) 398d71ae5a4SJacob Faibussowitsch { 399c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 400c4762a1bSJed Brown PetscReal norm_2, norm_max; 401c4762a1bSJed Brown 402*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 403c4762a1bSJed Brown /* 404c4762a1bSJed Brown View a graph of the current iterate 405c4762a1bSJed Brown */ 4069566063dSJacob Faibussowitsch PetscCall(VecView(u, appctx->viewer2)); 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* 409c4762a1bSJed Brown Compute the exact solution 410c4762a1bSJed Brown */ 4119566063dSJacob Faibussowitsch PetscCall(ExactSolution(time, appctx->solution, appctx)); 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* 414c4762a1bSJed Brown Print debugging information if desired 415c4762a1bSJed Brown */ 416c4762a1bSJed Brown if (appctx->debug) { 4179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(appctx->comm, "Computed solution vector\n")); 4189566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 4199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(appctx->comm, "Exact solution vector\n")); 4209566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD)); 421c4762a1bSJed Brown } 422c4762a1bSJed Brown 423c4762a1bSJed Brown /* 424c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 425c4762a1bSJed Brown */ 4269566063dSJacob Faibussowitsch PetscCall(VecAXPY(appctx->solution, -1.0, u)); 4279566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2)); 428c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h) * norm_2; 4299566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max)); 430c4762a1bSJed Brown if (norm_2 < 1e-14) norm_2 = 0; 431c4762a1bSJed Brown if (norm_max < 1e-14) norm_max = 0; 432c4762a1bSJed Brown 433c4762a1bSJed Brown /* 434c4762a1bSJed Brown PetscPrintf() causes only the first processor in this 435c4762a1bSJed Brown communicator to print the timestep information. 436c4762a1bSJed Brown */ 43763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(appctx->comm, "Timestep %" PetscInt_FMT ": time = %g 2-norm error = %g max norm error = %g\n", step, (double)time, (double)norm_2, (double)norm_max)); 438c4762a1bSJed Brown appctx->norm_2 += norm_2; 439c4762a1bSJed Brown appctx->norm_max += norm_max; 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* 442c4762a1bSJed Brown View a graph of the error 443c4762a1bSJed Brown */ 4449566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution, appctx->viewer1)); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* 447c4762a1bSJed Brown Print debugging information if desired 448c4762a1bSJed Brown */ 449c4762a1bSJed Brown if (appctx->debug) { 4509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(appctx->comm, "Error vector\n")); 4519566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD)); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 457c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 458c4762a1bSJed Brown /* 459c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 460c4762a1bSJed Brown matrix for the heat equation. 461c4762a1bSJed Brown 462c4762a1bSJed Brown Input Parameters: 463c4762a1bSJed Brown ts - the TS context 464c4762a1bSJed Brown t - current time 465c4762a1bSJed Brown global_in - global input vector 466c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 467c4762a1bSJed Brown 468c4762a1bSJed Brown Output Parameters: 469c4762a1bSJed Brown AA - Jacobian matrix 470c4762a1bSJed Brown BB - optionally different preconditioning matrix 471c4762a1bSJed Brown str - flag indicating matrix structure 472c4762a1bSJed Brown 473c4762a1bSJed Brown Notes: 474c4762a1bSJed Brown RHSMatrixHeat computes entries for the locally owned part of the system. 475c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 476c4762a1bSJed Brown contiguous chunks of rows across the processors. 477c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 478c4762a1bSJed Brown locally (but any non-local elements will be sent to the 479c4762a1bSJed Brown appropriate processor during matrix assembly). 480c4762a1bSJed Brown - Always specify global row and columns of matrix entries when 481c4762a1bSJed Brown using MatSetValues(); we could alternatively use MatSetValuesLocal(). 482c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 483c4762a1bSJed Brown - Note that MatSetValues() uses 0-based row and column numbers 484c4762a1bSJed Brown in Fortran as well as in C. 485c4762a1bSJed Brown */ 486d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec X, Mat AA, Mat BB, void *ctx) 487d71ae5a4SJacob Faibussowitsch { 488c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 489c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 490c4762a1bSJed Brown PetscInt i, mstart, mend, idx[3]; 491c4762a1bSJed Brown PetscScalar v[3], stwo = -2. / (appctx->h * appctx->h), sone = -.5 * stwo; 492c4762a1bSJed Brown 493*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 494c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 495c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 496c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 497c4762a1bSJed Brown 4989566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &mstart, &mend)); 499c4762a1bSJed Brown 500c4762a1bSJed Brown /* 501c4762a1bSJed Brown Set matrix rows corresponding to boundary data 502c4762a1bSJed Brown */ 503c4762a1bSJed Brown 504c4762a1bSJed Brown if (mstart == 0) { /* first processor only */ 505c4762a1bSJed Brown v[0] = 1.0; 5069566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mstart, 1, &mstart, v, INSERT_VALUES)); 507c4762a1bSJed Brown mstart++; 508c4762a1bSJed Brown } 509c4762a1bSJed Brown 510c4762a1bSJed Brown if (mend == appctx->m) { /* last processor only */ 511c4762a1bSJed Brown mend--; 512c4762a1bSJed Brown v[0] = 1.0; 5139566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mend, 1, &mend, v, INSERT_VALUES)); 514c4762a1bSJed Brown } 515c4762a1bSJed Brown 516c4762a1bSJed Brown /* 517c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 518c4762a1bSJed Brown matrix one row at a time. 519c4762a1bSJed Brown */ 5209371c9d4SSatish Balay v[0] = sone; 5219371c9d4SSatish Balay v[1] = stwo; 5229371c9d4SSatish Balay v[2] = sone; 523c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 5249371c9d4SSatish Balay idx[0] = i - 1; 5259371c9d4SSatish Balay idx[1] = i; 5269371c9d4SSatish Balay idx[2] = i + 1; 5279566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES)); 528c4762a1bSJed Brown } 529c4762a1bSJed Brown 530c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 531c4762a1bSJed Brown Complete the matrix assembly process and set some options 532c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 533c4762a1bSJed Brown /* 534c4762a1bSJed Brown Assemble matrix, using the 2-step process: 535c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 536c4762a1bSJed Brown Computations can be done while messages are in transition 537c4762a1bSJed Brown by placing code between these two statements. 538c4762a1bSJed Brown */ 5399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 541c4762a1bSJed Brown 542c4762a1bSJed Brown /* 543c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 544c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 545c4762a1bSJed Brown */ 5469566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 547c4762a1bSJed Brown 548*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 549c4762a1bSJed Brown } 550c4762a1bSJed Brown 551d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunctionHeat(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) 552d71ae5a4SJacob Faibussowitsch { 553c4762a1bSJed Brown Mat A; 554c4762a1bSJed Brown 555c4762a1bSJed Brown PetscFunctionBeginUser; 5569566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &A, NULL, NULL, &ctx)); 5579566063dSJacob Faibussowitsch PetscCall(RHSMatrixHeat(ts, t, globalin, A, NULL, ctx)); 5589566063dSJacob Faibussowitsch /* PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); */ 5599566063dSJacob Faibussowitsch PetscCall(MatMult(A, globalin, globalout)); 560*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 561c4762a1bSJed Brown } 562c4762a1bSJed Brown 563c4762a1bSJed Brown /*TEST 564c4762a1bSJed Brown 565c4762a1bSJed Brown test: 566c4762a1bSJed Brown args: -ts_view -nox 567c4762a1bSJed Brown 568c4762a1bSJed Brown test: 569c4762a1bSJed Brown suffix: 2 570c4762a1bSJed Brown args: -ts_view -nox 571c4762a1bSJed Brown nsize: 3 572c4762a1bSJed Brown 573c4762a1bSJed Brown test: 574c4762a1bSJed Brown suffix: 3 575c4762a1bSJed Brown args: -ts_view -nox -nonlinear 576c4762a1bSJed Brown 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown suffix: 4 579c4762a1bSJed Brown args: -ts_view -nox -nonlinear 580c4762a1bSJed Brown nsize: 3 581c4762a1bSJed Brown timeoutfactor: 3 582c4762a1bSJed Brown 583c4762a1bSJed Brown test: 584c4762a1bSJed Brown suffix: sundials 585e808b789SPatrick Sanan requires: sundials2 586c4762a1bSJed Brown args: -nox -ts_type sundials -ts_max_steps 5 -nonlinear 587c4762a1bSJed Brown nsize: 4 588c4762a1bSJed Brown 5897324063eSPatrick Sanan test: 5907324063eSPatrick Sanan suffix: sundials_dense 5917324063eSPatrick Sanan requires: sundials2 5927324063eSPatrick Sanan args: -nox -ts_type sundials -ts_sundials_use_dense -ts_max_steps 5 -nonlinear 5937324063eSPatrick Sanan nsize: 1 5947324063eSPatrick Sanan 595c4762a1bSJed Brown TEST*/ 596