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 -use_ifunc : Use IFunction/IJacobian interface\n\ 7c4762a1bSJed Brown -debug : Activate debugging printouts\n\ 8c4762a1bSJed Brown -nox : Deactivate x-window graphics\n\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* ------------------------------------------------------------------------ 11c4762a1bSJed Brown 12c4762a1bSJed Brown This program solves the one-dimensional heat equation (also called the 13c4762a1bSJed Brown diffusion equation), 14c4762a1bSJed Brown u_t = u_xx, 15c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 16c4762a1bSJed Brown u(t,0) = 0, u(t,1) = 0, 17c4762a1bSJed Brown and the initial condition 18c4762a1bSJed Brown u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x). 19c4762a1bSJed Brown This is a linear, second-order, parabolic equation. 20c4762a1bSJed Brown 21c4762a1bSJed Brown We discretize the right-hand side using finite differences with 22c4762a1bSJed Brown uniform grid spacing h: 23c4762a1bSJed Brown u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) 24c4762a1bSJed Brown We then demonstrate time evolution using the various TS methods by 25c4762a1bSJed Brown running the program via 26c4762a1bSJed Brown ex3 -ts_type <timestepping solver> 27c4762a1bSJed Brown 28c4762a1bSJed Brown We compare the approximate solution with the exact solution, given by 29c4762a1bSJed Brown u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) + 30c4762a1bSJed Brown 3*exp(-4*pi*pi*t) * sin(2*pi*x) 31c4762a1bSJed Brown 32c4762a1bSJed Brown Notes: 33c4762a1bSJed Brown This code demonstrates the TS solver interface to two variants of 34c4762a1bSJed Brown linear problems, u_t = f(u,t), namely 35c4762a1bSJed Brown - time-dependent f: f(u,t) is a function of t 36c4762a1bSJed Brown - time-independent f: f(u,t) is simply f(u) 37c4762a1bSJed Brown 38c4762a1bSJed Brown The parallel version of this code is ts/tutorials/ex4.c 39c4762a1bSJed Brown 40c4762a1bSJed Brown ------------------------------------------------------------------------- */ 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 44c4762a1bSJed Brown 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 <petscts.h> 53c4762a1bSJed Brown #include <petscdraw.h> 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 56c4762a1bSJed Brown User-defined application context - contains data needed by the 57c4762a1bSJed Brown application-provided call-back routines. 58c4762a1bSJed Brown */ 59c4762a1bSJed Brown typedef struct { 60c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 61c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 62c4762a1bSJed Brown PetscReal h; /* mesh width h = 1/(m-1) */ 63c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 64c4762a1bSJed Brown PetscViewer viewer1,viewer2; /* viewers for the solution and error */ 65c4762a1bSJed Brown PetscReal norm_2,norm_max; /* error norms */ 66c4762a1bSJed Brown Mat A; /* RHS mat, used with IFunction interface */ 67c4762a1bSJed Brown PetscReal oshift; /* old shift applied, prevent to recompute the IJacobian */ 68c4762a1bSJed Brown } AppCtx; 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* 71c4762a1bSJed Brown User-defined routines 72c4762a1bSJed Brown */ 73c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 74c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*); 75c4762a1bSJed Brown extern PetscErrorCode IFunctionHeat(TS,PetscReal,Vec,Vec,Vec,void*); 76c4762a1bSJed Brown extern PetscErrorCode IJacobianHeat(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 77c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 78c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); 79c4762a1bSJed Brown 80c4762a1bSJed Brown int main(int argc,char **argv) 81c4762a1bSJed Brown { 82c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 83c4762a1bSJed Brown TS ts; /* timestepping context */ 84c4762a1bSJed Brown Mat A; /* matrix data structure */ 85c4762a1bSJed Brown Vec u; /* approximate solution vector */ 86c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 87c4762a1bSJed Brown PetscInt time_steps_max = 100; /* default max timesteps */ 88c4762a1bSJed Brown PetscDraw draw; /* drawing context */ 89c4762a1bSJed Brown PetscInt steps,m; 90c4762a1bSJed Brown PetscMPIInt size; 91c4762a1bSJed Brown PetscReal dt; 92c4762a1bSJed Brown PetscBool flg,flg_string; 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4762a1bSJed Brown Initialize program and set problem parameters 96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 97c4762a1bSJed Brown 98*327415f7SBarry Smith PetscFunctionBeginUser; 999566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1013c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 102c4762a1bSJed Brown 103c4762a1bSJed Brown m = 60; 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug)); 106c4762a1bSJed Brown flg_string = PETSC_FALSE; 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_string_viewer",&flg_string,NULL)); 108c4762a1bSJed Brown 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 PetscCall(PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n")); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Create vector data structures 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 122c4762a1bSJed Brown */ 1239566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&u)); 1249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.solution)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Set up displays to show graphs of the solution and error 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 129c4762a1bSJed Brown 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw)); 1329566063dSJacob Faibussowitsch PetscCall(PetscDrawSetDoubleBuffer(draw)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2)); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw)); 1359566063dSJacob Faibussowitsch PetscCall(PetscDrawSetDoubleBuffer(draw)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Create timestepping solver context 139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF,&ts)); 1429566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_LINEAR)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145c4762a1bSJed Brown Set optional user-defined monitoring routine 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147c4762a1bSJed Brown 148c4762a1bSJed Brown if (!flg_string) { 1499566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&appctx,NULL)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153c4762a1bSJed Brown 154c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 155c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 1589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m)); 1599566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1609566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown flg = PETSC_FALSE; 1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_ifunc",&flg,NULL)); 164c4762a1bSJed Brown if (!flg) { 165c4762a1bSJed Brown appctx.A = NULL; 1669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL)); 167c4762a1bSJed Brown if (flg) { 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 170c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 171c4762a1bSJed Brown as a time-dependent matrix. 172c4762a1bSJed Brown */ 1739566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 1749566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx)); 175c4762a1bSJed Brown } else { 176c4762a1bSJed Brown /* 177c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 178c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 179c4762a1bSJed Brown as a matrix only once, and then sets the special Jacobian evaluation 180c4762a1bSJed Brown routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian. 181c4762a1bSJed Brown */ 1829566063dSJacob Faibussowitsch PetscCall(RHSMatrixHeat(ts,0.0,u,A,A,&appctx)); 1839566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 1849566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx)); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } else { 187c4762a1bSJed Brown Mat J; 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(RHSMatrixHeat(ts,0.0,u,A,A,&appctx)); 1909566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&J)); 1919566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,IFunctionHeat,&appctx)); 1929566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,J,J,IJacobianHeat,&appctx)); 1939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 194c4762a1bSJed Brown 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 196c4762a1bSJed Brown appctx.A = A; 197c4762a1bSJed Brown appctx.oshift = PETSC_MIN_REAL; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 200c4762a1bSJed Brown Set solution vector and initial timestep 201c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 202c4762a1bSJed Brown 203c4762a1bSJed Brown dt = appctx.h*appctx.h/2.0; 2049566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown Customize timestepping solver: 208c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 209c4762a1bSJed Brown - Set timestepping duration info 210c4762a1bSJed Brown Then set runtime options, which can override these defaults. 211c4762a1bSJed Brown For example, 212c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 213c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 214c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 215c4762a1bSJed Brown 2169566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,time_steps_max)); 2179566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,time_total_max)); 2189566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 2199566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222c4762a1bSJed Brown Solve the problem 223c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* 226c4762a1bSJed Brown Evaluate initial conditions 227c4762a1bSJed Brown */ 2289566063dSJacob Faibussowitsch PetscCall(InitialConditions(u,&appctx)); 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* 231c4762a1bSJed Brown Run the timestepping solver 232c4762a1bSJed Brown */ 2339566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,u)); 2349566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 237c4762a1bSJed Brown View timestepping solver info 238c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 239c4762a1bSJed Brown 2409566063dSJacob 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))); 241c4762a1bSJed Brown if (!flg_string) { 2429566063dSJacob Faibussowitsch PetscCall(TSView(ts,PETSC_VIEWER_STDOUT_SELF)); 243c4762a1bSJed Brown } else { 244c4762a1bSJed Brown PetscViewer stringviewer; 245c4762a1bSJed Brown char string[512]; 246c4762a1bSJed Brown const char *outstring; 247c4762a1bSJed Brown 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_WORLD,string,sizeof(string),&stringviewer)); 2499566063dSJacob Faibussowitsch PetscCall(TSView(ts,stringviewer)); 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerStringGetStringRead(stringviewer,&outstring,NULL)); 2513c633725SBarry Smith PetscCheck((char*)outstring == (char*)string,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"String returned from viewer does not equal original string"); 2529566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Output from string viewer:%s\n",outstring)); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&stringviewer)); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 257c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 258c4762a1bSJed Brown are no longer needed. 259c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&appctx.viewer1)); 2659566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&appctx.viewer2)); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.solution)); 2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.A)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* 270c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 271c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 272c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 273c4762a1bSJed Brown options are chosen (e.g., -log_view). 274c4762a1bSJed Brown */ 2759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 276b122ec5aSJacob Faibussowitsch return 0; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 281c4762a1bSJed Brown 282c4762a1bSJed Brown Input Parameter: 283c4762a1bSJed Brown u - uninitialized solution vector (global) 284c4762a1bSJed Brown appctx - user-defined application context 285c4762a1bSJed Brown 286c4762a1bSJed Brown Output Parameter: 287c4762a1bSJed Brown u - vector with solution at initial time (global) 288c4762a1bSJed Brown */ 289c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 290c4762a1bSJed Brown { 291c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h; 292c4762a1bSJed Brown PetscInt i; 293c4762a1bSJed Brown 294c4762a1bSJed Brown /* 295c4762a1bSJed Brown Get a pointer to vector data. 296c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 297c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 298c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 299c4762a1bSJed Brown the array. 300c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 301c4762a1bSJed Brown C version. See the users manual for details. 302c4762a1bSJed Brown */ 3039566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(u,&u_localptr)); 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* 306c4762a1bSJed Brown We initialize the solution array by simply writing the solution 307c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 308c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 309c4762a1bSJed Brown */ 310c4762a1bSJed Brown for (i=0; i<appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h); 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* 313c4762a1bSJed Brown Restore vector 314c4762a1bSJed Brown */ 3159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(u,&u_localptr)); 316c4762a1bSJed Brown 317c4762a1bSJed Brown /* 318c4762a1bSJed Brown Print debugging information if desired 319c4762a1bSJed Brown */ 320c4762a1bSJed Brown if (appctx->debug) { 3219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n")); 3229566063dSJacob Faibussowitsch PetscCall(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 323c4762a1bSJed Brown } 324c4762a1bSJed Brown 325c4762a1bSJed Brown return 0; 326c4762a1bSJed Brown } 327c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 328c4762a1bSJed Brown /* 329c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 330c4762a1bSJed Brown 331c4762a1bSJed Brown Input Parameters: 332c4762a1bSJed Brown t - current time 333c4762a1bSJed Brown solution - vector in which exact solution will be computed 334c4762a1bSJed Brown appctx - user-defined application context 335c4762a1bSJed Brown 336c4762a1bSJed Brown Output Parameter: 337c4762a1bSJed Brown solution - vector with the newly computed exact solution 338c4762a1bSJed Brown */ 339c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 340c4762a1bSJed Brown { 341c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t; 342c4762a1bSJed Brown PetscInt i; 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* 345c4762a1bSJed Brown Get a pointer to vector data. 346c4762a1bSJed Brown */ 3479566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(solution,&s_localptr)); 348c4762a1bSJed Brown 349c4762a1bSJed Brown /* 350c4762a1bSJed Brown Simply write the solution directly into the array locations. 351c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 352c4762a1bSJed Brown */ 353c4762a1bSJed Brown ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); 354c4762a1bSJed Brown ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc); 355c4762a1bSJed Brown sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 356c4762a1bSJed Brown for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2; 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* 359c4762a1bSJed Brown Restore vector 360c4762a1bSJed Brown */ 3619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(solution,&s_localptr)); 362c4762a1bSJed Brown return 0; 363c4762a1bSJed Brown } 364c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 365c4762a1bSJed Brown /* 366c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 367c4762a1bSJed Brown each timestep. This example plots the solution and computes the 368c4762a1bSJed Brown error in two different norms. 369c4762a1bSJed Brown 370c4762a1bSJed Brown This example also demonstrates changing the timestep via TSSetTimeStep(). 371c4762a1bSJed Brown 372c4762a1bSJed Brown Input Parameters: 373c4762a1bSJed Brown ts - the timestep context 374c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 375c4762a1bSJed Brown initial condition) 376c4762a1bSJed Brown time - the current time 377c4762a1bSJed Brown u - the solution at this timestep 378c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 379c4762a1bSJed Brown In this case we use the application context which contains 380c4762a1bSJed Brown information about the problem size, workspace and the exact 381c4762a1bSJed Brown solution. 382c4762a1bSJed Brown */ 383c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 384c4762a1bSJed Brown { 385c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 386c4762a1bSJed Brown PetscReal norm_2,norm_max,dt,dttol; 387c4762a1bSJed Brown 388c4762a1bSJed Brown /* 389c4762a1bSJed Brown View a graph of the current iterate 390c4762a1bSJed Brown */ 3919566063dSJacob Faibussowitsch PetscCall(VecView(u,appctx->viewer2)); 392c4762a1bSJed Brown 393c4762a1bSJed Brown /* 394c4762a1bSJed Brown Compute the exact solution 395c4762a1bSJed Brown */ 3969566063dSJacob Faibussowitsch PetscCall(ExactSolution(time,appctx->solution,appctx)); 397c4762a1bSJed Brown 398c4762a1bSJed Brown /* 399c4762a1bSJed Brown Print debugging information if desired 400c4762a1bSJed Brown */ 401c4762a1bSJed Brown if (appctx->debug) { 4029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n")); 4039566063dSJacob Faibussowitsch PetscCall(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 4049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n")); 4059566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 406c4762a1bSJed Brown } 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* 409c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 410c4762a1bSJed Brown */ 4119566063dSJacob Faibussowitsch PetscCall(VecAXPY(appctx->solution,-1.0,u)); 4129566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution,NORM_2,&norm_2)); 413c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h)*norm_2; 4149566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution,NORM_MAX,&norm_max)); 415c4762a1bSJed Brown 4169566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 41763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %3" PetscInt_FMT ": step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)dt,(double)time,(double)norm_2,(double)norm_max)); 418c4762a1bSJed Brown 419c4762a1bSJed Brown appctx->norm_2 += norm_2; 420c4762a1bSJed Brown appctx->norm_max += norm_max; 421c4762a1bSJed Brown 422c4762a1bSJed Brown dttol = .0001; 4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL)); 424c4762a1bSJed Brown if (dt < dttol) { 425c4762a1bSJed Brown dt *= .999; 4269566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 427c4762a1bSJed Brown } 428c4762a1bSJed Brown 429c4762a1bSJed Brown /* 430c4762a1bSJed Brown View a graph of the error 431c4762a1bSJed Brown */ 4329566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution,appctx->viewer1)); 433c4762a1bSJed Brown 434c4762a1bSJed Brown /* 435c4762a1bSJed Brown Print debugging information if desired 436c4762a1bSJed Brown */ 437c4762a1bSJed Brown if (appctx->debug) { 4389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error vector\n")); 4399566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 440c4762a1bSJed Brown } 441c4762a1bSJed Brown 442c4762a1bSJed Brown return 0; 443c4762a1bSJed Brown } 444c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 447c4762a1bSJed Brown matrix for the heat equation. 448c4762a1bSJed Brown 449c4762a1bSJed Brown Input Parameters: 450c4762a1bSJed Brown ts - the TS context 451c4762a1bSJed Brown t - current time 452c4762a1bSJed Brown global_in - global input vector 453c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 454c4762a1bSJed Brown 455c4762a1bSJed Brown Output Parameters: 456c4762a1bSJed Brown AA - Jacobian matrix 457c4762a1bSJed Brown BB - optionally different preconditioning matrix 458c4762a1bSJed Brown str - flag indicating matrix structure 459c4762a1bSJed Brown 460c4762a1bSJed Brown Notes: 461c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 462c4762a1bSJed Brown in Fortran as well as in C. 463c4762a1bSJed Brown */ 464c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx) 465c4762a1bSJed Brown { 466c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 467c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 468c4762a1bSJed Brown PetscInt mstart = 0; 469c4762a1bSJed Brown PetscInt mend = appctx->m; 470c4762a1bSJed Brown PetscInt i,idx[3]; 471c4762a1bSJed Brown PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; 472c4762a1bSJed Brown 473c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 474c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 475c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 476c4762a1bSJed Brown /* 477c4762a1bSJed Brown Set matrix rows corresponding to boundary data 478c4762a1bSJed Brown */ 479c4762a1bSJed Brown 480c4762a1bSJed Brown mstart = 0; 481c4762a1bSJed Brown v[0] = 1.0; 4829566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES)); 483c4762a1bSJed Brown mstart++; 484c4762a1bSJed Brown 485c4762a1bSJed Brown mend--; 486c4762a1bSJed Brown v[0] = 1.0; 4879566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES)); 488c4762a1bSJed Brown 489c4762a1bSJed Brown /* 490c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 491c4762a1bSJed Brown matrix one row at a time. 492c4762a1bSJed Brown */ 493c4762a1bSJed Brown v[0] = sone; v[1] = stwo; v[2] = sone; 494c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 495c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 4969566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES)); 497c4762a1bSJed Brown } 498c4762a1bSJed Brown 499c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 500c4762a1bSJed Brown Complete the matrix assembly process and set some options 501c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 502c4762a1bSJed Brown /* 503c4762a1bSJed Brown Assemble matrix, using the 2-step process: 504c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 505c4762a1bSJed Brown Computations can be done while messages are in transition 506c4762a1bSJed Brown by placing code between these two statements. 507c4762a1bSJed Brown */ 5089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 510c4762a1bSJed Brown 511c4762a1bSJed Brown /* 512c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 513c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 514c4762a1bSJed Brown */ 5159566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 516c4762a1bSJed Brown 517c4762a1bSJed Brown return 0; 518c4762a1bSJed Brown } 519c4762a1bSJed Brown 520c4762a1bSJed Brown PetscErrorCode IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,void *ctx) 521c4762a1bSJed Brown { 522c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 523c4762a1bSJed Brown 5249566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->A,X,r)); 5259566063dSJacob Faibussowitsch PetscCall(VecAYPX(r,-1.0,Xdot)); 526c4762a1bSJed Brown return 0; 527c4762a1bSJed Brown } 528c4762a1bSJed Brown 529c4762a1bSJed Brown PetscErrorCode IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,void *ctx) 530c4762a1bSJed Brown { 531c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 532c4762a1bSJed Brown 533c4762a1bSJed Brown if (appctx->oshift == s) return 0; 5349566063dSJacob Faibussowitsch PetscCall(MatCopy(appctx->A,A,SAME_NONZERO_PATTERN)); 5359566063dSJacob Faibussowitsch PetscCall(MatScale(A,-1)); 5369566063dSJacob Faibussowitsch PetscCall(MatShift(A,s)); 5379566063dSJacob Faibussowitsch PetscCall(MatCopy(A,B,SAME_NONZERO_PATTERN)); 538c4762a1bSJed Brown appctx->oshift = s; 539c4762a1bSJed Brown return 0; 540c4762a1bSJed Brown } 541c4762a1bSJed Brown 542c4762a1bSJed Brown /*TEST 543c4762a1bSJed Brown 544c4762a1bSJed Brown test: 545c4762a1bSJed Brown args: -nox -ts_type ssp -ts_dt 0.0005 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 2 549c4762a1bSJed Brown args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1 550c4762a1bSJed Brown 551c4762a1bSJed Brown test: 552c4762a1bSJed Brown suffix: 3 553c4762a1bSJed Brown args: -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason 554c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 555c4762a1bSJed Brown requires: !single 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 4 559c4762a1bSJed Brown args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason 560c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 561c4762a1bSJed Brown 562c4762a1bSJed Brown test: 563c4762a1bSJed Brown suffix: 5 564c4762a1bSJed Brown args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs 565c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 566c4762a1bSJed Brown 567c4762a1bSJed Brown test: 568c4762a1bSJed Brown requires: !single 569c4762a1bSJed Brown suffix: pod_guess 570c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason 571c4762a1bSJed Brown 572c4762a1bSJed Brown test: 573c4762a1bSJed Brown requires: !single 574c4762a1bSJed Brown suffix: pod_guess_Ainner 575c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -ksp_guess_pod_Ainner -pc_type none -ksp_converged_reason 576c4762a1bSJed Brown 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown requires: !single 579c4762a1bSJed Brown suffix: fischer_guess 580c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason 581c4762a1bSJed Brown 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown requires: !single 584c4762a1bSJed Brown suffix: fischer_guess_2 585c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 2,10 -pc_type none -ksp_converged_reason 586c4762a1bSJed Brown 587c4762a1bSJed Brown test: 588c4762a1bSJed Brown requires: !single 589cbb17d71SDavid Wells suffix: fischer_guess_3 590cbb17d71SDavid Wells args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 3,10 -pc_type none -ksp_converged_reason 591cbb17d71SDavid Wells 592cbb17d71SDavid Wells test: 593cbb17d71SDavid Wells requires: !single 594c4762a1bSJed Brown suffix: stringview 595c4762a1bSJed Brown args: -nox -ts_type rosw -test_string_viewer 596c4762a1bSJed Brown 597c4762a1bSJed Brown test: 598c4762a1bSJed Brown requires: !single 599c4762a1bSJed Brown suffix: stringview_euler 600c4762a1bSJed Brown args: -nox -ts_type euler -test_string_viewer 601c4762a1bSJed Brown 602c4762a1bSJed Brown TEST*/ 603