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 Concepts: TS^time-dependent linear problems 12c4762a1bSJed Brown Concepts: TS^heat equation 13c4762a1bSJed Brown Concepts: TS^diffusion equation 14c4762a1bSJed Brown Processors: 1 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* ------------------------------------------------------------------------ 18c4762a1bSJed Brown 19c4762a1bSJed Brown This program solves the one-dimensional heat equation (also called the 20c4762a1bSJed Brown diffusion equation), 21c4762a1bSJed Brown u_t = u_xx, 22c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 23c4762a1bSJed Brown u(t,0) = 0, u(t,1) = 0, 24c4762a1bSJed Brown and the initial condition 25c4762a1bSJed Brown u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x). 26c4762a1bSJed Brown This is a linear, second-order, parabolic equation. 27c4762a1bSJed Brown 28c4762a1bSJed Brown We discretize the right-hand side using finite differences with 29c4762a1bSJed Brown uniform grid spacing h: 30c4762a1bSJed Brown u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) 31c4762a1bSJed Brown We then demonstrate time evolution using the various TS methods by 32c4762a1bSJed Brown running the program via 33c4762a1bSJed Brown ex3 -ts_type <timestepping solver> 34c4762a1bSJed Brown 35c4762a1bSJed Brown We compare the approximate solution with the exact solution, given by 36c4762a1bSJed Brown u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) + 37c4762a1bSJed Brown 3*exp(-4*pi*pi*t) * sin(2*pi*x) 38c4762a1bSJed Brown 39c4762a1bSJed Brown Notes: 40c4762a1bSJed Brown This code demonstrates the TS solver interface to two variants of 41c4762a1bSJed Brown linear problems, u_t = f(u,t), namely 42c4762a1bSJed Brown - time-dependent f: f(u,t) is a function of t 43c4762a1bSJed Brown - time-independent f: f(u,t) is simply f(u) 44c4762a1bSJed Brown 45c4762a1bSJed Brown The parallel version of this code is ts/tutorials/ex4.c 46c4762a1bSJed Brown 47c4762a1bSJed Brown ------------------------------------------------------------------------- */ 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 51c4762a1bSJed Brown automatically includes: 52c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 53c4762a1bSJed Brown petscmat.h - matrices 54c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 55c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 56c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 57c4762a1bSJed Brown */ 58c4762a1bSJed Brown 59c4762a1bSJed Brown #include <petscts.h> 60c4762a1bSJed Brown #include <petscdraw.h> 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown User-defined application context - contains data needed by the 64c4762a1bSJed Brown application-provided call-back routines. 65c4762a1bSJed Brown */ 66c4762a1bSJed Brown typedef struct { 67c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 68c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 69c4762a1bSJed Brown PetscReal h; /* mesh width h = 1/(m-1) */ 70c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 71c4762a1bSJed Brown PetscViewer viewer1,viewer2; /* viewers for the solution and error */ 72c4762a1bSJed Brown PetscReal norm_2,norm_max; /* error norms */ 73c4762a1bSJed Brown Mat A; /* RHS mat, used with IFunction interface */ 74c4762a1bSJed Brown PetscReal oshift; /* old shift applied, prevent to recompute the IJacobian */ 75c4762a1bSJed Brown } AppCtx; 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* 78c4762a1bSJed Brown User-defined routines 79c4762a1bSJed Brown */ 80c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 81c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*); 82c4762a1bSJed Brown extern PetscErrorCode IFunctionHeat(TS,PetscReal,Vec,Vec,Vec,void*); 83c4762a1bSJed Brown extern PetscErrorCode IJacobianHeat(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 84c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 85c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); 86c4762a1bSJed Brown 87c4762a1bSJed Brown int main(int argc,char **argv) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 90c4762a1bSJed Brown TS ts; /* timestepping context */ 91c4762a1bSJed Brown Mat A; /* matrix data structure */ 92c4762a1bSJed Brown Vec u; /* approximate solution vector */ 93c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 94c4762a1bSJed Brown PetscInt time_steps_max = 100; /* default max timesteps */ 95c4762a1bSJed Brown PetscDraw draw; /* drawing context */ 96c4762a1bSJed Brown PetscInt steps,m; 97c4762a1bSJed Brown PetscMPIInt size; 98c4762a1bSJed Brown PetscReal dt; 99c4762a1bSJed Brown PetscBool flg,flg_string; 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102c4762a1bSJed Brown Initialize program and set problem parameters 103c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 104c4762a1bSJed Brown 105*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 1065f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1073c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 108c4762a1bSJed Brown 109c4762a1bSJed Brown m = 60; 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug)); 112c4762a1bSJed Brown flg_string = PETSC_FALSE; 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_string_viewer",&flg_string,NULL)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown appctx.m = m; 116c4762a1bSJed Brown appctx.h = 1.0/(m-1.0); 117c4762a1bSJed Brown appctx.norm_2 = 0.0; 118c4762a1bSJed Brown appctx.norm_max = 0.0; 119c4762a1bSJed Brown 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n")); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Create vector data structures 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* 127c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 128c4762a1bSJed Brown */ 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&u)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.solution)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 133c4762a1bSJed Brown Set up displays to show graphs of the solution and error 134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 135c4762a1bSJed Brown 1365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Create timestepping solver context 145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 146c4762a1bSJed Brown 1475f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_LINEAR)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Set optional user-defined monitoring routine 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153c4762a1bSJed Brown 154c4762a1bSJed Brown if (!flg_string) { 1555f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL)); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159c4762a1bSJed Brown 160c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 161c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 162c4762a1bSJed Brown 1635f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown flg = PETSC_FALSE; 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_ifunc",&flg,NULL)); 170c4762a1bSJed Brown if (!flg) { 171c4762a1bSJed Brown appctx.A = NULL; 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL)); 173c4762a1bSJed Brown if (flg) { 174c4762a1bSJed Brown /* 175c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 176c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 177c4762a1bSJed Brown as a time-dependent matrix. 178c4762a1bSJed Brown */ 1795f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx)); 181c4762a1bSJed Brown } else { 182c4762a1bSJed Brown /* 183c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 184c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 185c4762a1bSJed Brown as a matrix only once, and then sets the special Jacobian evaluation 186c4762a1bSJed Brown routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian. 187c4762a1bSJed Brown */ 1885f80ce2aSJacob Faibussowitsch CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx)); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown } else { 193c4762a1bSJed Brown Mat J; 194c4762a1bSJed Brown 1955f80ce2aSJacob Faibussowitsch CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&J)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunctionHeat,&appctx)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,IJacobianHeat,&appctx)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 200c4762a1bSJed Brown 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)A)); 202c4762a1bSJed Brown appctx.A = A; 203c4762a1bSJed Brown appctx.oshift = PETSC_MIN_REAL; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206c4762a1bSJed Brown Set solution vector and initial timestep 207c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 208c4762a1bSJed Brown 209c4762a1bSJed Brown dt = appctx.h*appctx.h/2.0; 2105f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 213c4762a1bSJed Brown Customize timestepping solver: 214c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 215c4762a1bSJed Brown - Set timestepping duration info 216c4762a1bSJed Brown Then set runtime options, which can override these defaults. 217c4762a1bSJed Brown For example, 218c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 219c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 220c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 221c4762a1bSJed Brown 2225f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,time_steps_max)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,time_total_max)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 228c4762a1bSJed Brown Solve the problem 229c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* 232c4762a1bSJed Brown Evaluate initial conditions 233c4762a1bSJed Brown */ 2345f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(u,&appctx)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* 237c4762a1bSJed Brown Run the timestepping solver 238c4762a1bSJed Brown */ 2395f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,u)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 243c4762a1bSJed Brown View timestepping solver info 244c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 245c4762a1bSJed Brown 2465f80ce2aSJacob Faibussowitsch CHKERRQ(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))); 247c4762a1bSJed Brown if (!flg_string) { 2485f80ce2aSJacob Faibussowitsch CHKERRQ(TSView(ts,PETSC_VIEWER_STDOUT_SELF)); 249c4762a1bSJed Brown } else { 250c4762a1bSJed Brown PetscViewer stringviewer; 251c4762a1bSJed Brown char string[512]; 252c4762a1bSJed Brown const char *outstring; 253c4762a1bSJed Brown 2545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringOpen(PETSC_COMM_WORLD,string,sizeof(string),&stringviewer)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(TSView(ts,stringviewer)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringGetStringRead(stringviewer,&outstring,NULL)); 2573c633725SBarry Smith PetscCheck((char*)outstring == (char*)string,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"String returned from viewer does not equal original string"); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Output from string viewer:%s\n",outstring)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&stringviewer)); 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 263c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 264c4762a1bSJed Brown are no longer needed. 265c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 266c4762a1bSJed Brown 2675f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&appctx.viewer1)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&appctx.viewer2)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.solution)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.A)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* 276c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 277c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 278c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 279c4762a1bSJed Brown options are chosen (e.g., -log_view). 280c4762a1bSJed Brown */ 281*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 282*b122ec5aSJacob Faibussowitsch return 0; 283c4762a1bSJed Brown } 284c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 285c4762a1bSJed Brown /* 286c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 287c4762a1bSJed Brown 288c4762a1bSJed Brown Input Parameter: 289c4762a1bSJed Brown u - uninitialized solution vector (global) 290c4762a1bSJed Brown appctx - user-defined application context 291c4762a1bSJed Brown 292c4762a1bSJed Brown Output Parameter: 293c4762a1bSJed Brown u - vector with solution at initial time (global) 294c4762a1bSJed Brown */ 295c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 296c4762a1bSJed Brown { 297c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h; 298c4762a1bSJed Brown PetscInt i; 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* 301c4762a1bSJed Brown Get a pointer to vector data. 302c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 303c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 304c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 305c4762a1bSJed Brown the array. 306c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 307c4762a1bSJed Brown C version. See the users manual for details. 308c4762a1bSJed Brown */ 3095f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(u,&u_localptr)); 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* 312c4762a1bSJed Brown We initialize the solution array by simply writing the solution 313c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 314c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 315c4762a1bSJed Brown */ 316c4762a1bSJed Brown for (i=0; i<appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h); 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* 319c4762a1bSJed Brown Restore vector 320c4762a1bSJed Brown */ 3215f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(u,&u_localptr)); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* 324c4762a1bSJed Brown Print debugging information if desired 325c4762a1bSJed Brown */ 326c4762a1bSJed Brown if (appctx->debug) { 3275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n")); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown return 0; 332c4762a1bSJed Brown } 333c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 334c4762a1bSJed Brown /* 335c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 336c4762a1bSJed Brown 337c4762a1bSJed Brown Input Parameters: 338c4762a1bSJed Brown t - current time 339c4762a1bSJed Brown solution - vector in which exact solution will be computed 340c4762a1bSJed Brown appctx - user-defined application context 341c4762a1bSJed Brown 342c4762a1bSJed Brown Output Parameter: 343c4762a1bSJed Brown solution - vector with the newly computed exact solution 344c4762a1bSJed Brown */ 345c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 346c4762a1bSJed Brown { 347c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t; 348c4762a1bSJed Brown PetscInt i; 349c4762a1bSJed Brown 350c4762a1bSJed Brown /* 351c4762a1bSJed Brown Get a pointer to vector data. 352c4762a1bSJed Brown */ 3535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(solution,&s_localptr)); 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* 356c4762a1bSJed Brown Simply write the solution directly into the array locations. 357c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 358c4762a1bSJed Brown */ 359c4762a1bSJed Brown ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); 360c4762a1bSJed Brown ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc); 361c4762a1bSJed Brown sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 362c4762a1bSJed Brown for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2; 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* 365c4762a1bSJed Brown Restore vector 366c4762a1bSJed Brown */ 3675f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(solution,&s_localptr)); 368c4762a1bSJed Brown return 0; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 371c4762a1bSJed Brown /* 372c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 373c4762a1bSJed Brown each timestep. This example plots the solution and computes the 374c4762a1bSJed Brown error in two different norms. 375c4762a1bSJed Brown 376c4762a1bSJed Brown This example also demonstrates changing the timestep via TSSetTimeStep(). 377c4762a1bSJed Brown 378c4762a1bSJed Brown Input Parameters: 379c4762a1bSJed Brown ts - the timestep context 380c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 381c4762a1bSJed Brown initial condition) 382c4762a1bSJed Brown time - the current time 383c4762a1bSJed Brown u - the solution at this timestep 384c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 385c4762a1bSJed Brown In this case we use the application context which contains 386c4762a1bSJed Brown information about the problem size, workspace and the exact 387c4762a1bSJed Brown solution. 388c4762a1bSJed Brown */ 389c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 390c4762a1bSJed Brown { 391c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 392c4762a1bSJed Brown PetscReal norm_2,norm_max,dt,dttol; 393c4762a1bSJed Brown 394c4762a1bSJed Brown /* 395c4762a1bSJed Brown View a graph of the current iterate 396c4762a1bSJed Brown */ 3975f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,appctx->viewer2)); 398c4762a1bSJed Brown 399c4762a1bSJed Brown /* 400c4762a1bSJed Brown Compute the exact solution 401c4762a1bSJed Brown */ 4025f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(time,appctx->solution,appctx)); 403c4762a1bSJed Brown 404c4762a1bSJed Brown /* 405c4762a1bSJed Brown Print debugging information if desired 406c4762a1bSJed Brown */ 407c4762a1bSJed Brown if (appctx->debug) { 4085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n")); 4095f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 4105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n")); 4115f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* 415c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 416c4762a1bSJed Brown */ 4175f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(appctx->solution,-1.0,u)); 4185f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_2,&norm_2)); 419c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h)*norm_2; 4205f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&norm_max)); 421c4762a1bSJed Brown 4225f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Timestep %3D: 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)); 424c4762a1bSJed Brown 425c4762a1bSJed Brown appctx->norm_2 += norm_2; 426c4762a1bSJed Brown appctx->norm_max += norm_max; 427c4762a1bSJed Brown 428c4762a1bSJed Brown dttol = .0001; 4295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL)); 430c4762a1bSJed Brown if (dt < dttol) { 431c4762a1bSJed Brown dt *= .999; 4325f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* 436c4762a1bSJed Brown View a graph of the error 437c4762a1bSJed Brown */ 4385f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,appctx->viewer1)); 439c4762a1bSJed Brown 440c4762a1bSJed Brown /* 441c4762a1bSJed Brown Print debugging information if desired 442c4762a1bSJed Brown */ 443c4762a1bSJed Brown if (appctx->debug) { 4445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error vector\n")); 4455f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 446c4762a1bSJed Brown } 447c4762a1bSJed Brown 448c4762a1bSJed Brown return 0; 449c4762a1bSJed Brown } 450c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 451c4762a1bSJed Brown /* 452c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 453c4762a1bSJed Brown matrix for the heat equation. 454c4762a1bSJed Brown 455c4762a1bSJed Brown Input Parameters: 456c4762a1bSJed Brown ts - the TS context 457c4762a1bSJed Brown t - current time 458c4762a1bSJed Brown global_in - global input vector 459c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 460c4762a1bSJed Brown 461c4762a1bSJed Brown Output Parameters: 462c4762a1bSJed Brown AA - Jacobian matrix 463c4762a1bSJed Brown BB - optionally different preconditioning matrix 464c4762a1bSJed Brown str - flag indicating matrix structure 465c4762a1bSJed Brown 466c4762a1bSJed Brown Notes: 467c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 468c4762a1bSJed Brown in Fortran as well as in C. 469c4762a1bSJed Brown */ 470c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx) 471c4762a1bSJed Brown { 472c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 473c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 474c4762a1bSJed Brown PetscInt mstart = 0; 475c4762a1bSJed Brown PetscInt mend = appctx->m; 476c4762a1bSJed Brown PetscInt i,idx[3]; 477c4762a1bSJed Brown PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; 478c4762a1bSJed Brown 479c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 480c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 481c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 482c4762a1bSJed Brown /* 483c4762a1bSJed Brown Set matrix rows corresponding to boundary data 484c4762a1bSJed Brown */ 485c4762a1bSJed Brown 486c4762a1bSJed Brown mstart = 0; 487c4762a1bSJed Brown v[0] = 1.0; 4885f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES)); 489c4762a1bSJed Brown mstart++; 490c4762a1bSJed Brown 491c4762a1bSJed Brown mend--; 492c4762a1bSJed Brown v[0] = 1.0; 4935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES)); 494c4762a1bSJed Brown 495c4762a1bSJed Brown /* 496c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 497c4762a1bSJed Brown matrix one row at a time. 498c4762a1bSJed Brown */ 499c4762a1bSJed Brown v[0] = sone; v[1] = stwo; v[2] = sone; 500c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 501c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 5025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES)); 503c4762a1bSJed Brown } 504c4762a1bSJed Brown 505c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 506c4762a1bSJed Brown Complete the matrix assembly process and set some options 507c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 508c4762a1bSJed Brown /* 509c4762a1bSJed Brown Assemble matrix, using the 2-step process: 510c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 511c4762a1bSJed Brown Computations can be done while messages are in transition 512c4762a1bSJed Brown by placing code between these two statements. 513c4762a1bSJed Brown */ 5145f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5155f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 516c4762a1bSJed Brown 517c4762a1bSJed Brown /* 518c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 519c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 520c4762a1bSJed Brown */ 5215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 522c4762a1bSJed Brown 523c4762a1bSJed Brown return 0; 524c4762a1bSJed Brown } 525c4762a1bSJed Brown 526c4762a1bSJed Brown PetscErrorCode IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,void *ctx) 527c4762a1bSJed Brown { 528c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 529c4762a1bSJed Brown 5305f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(appctx->A,X,r)); 5315f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(r,-1.0,Xdot)); 532c4762a1bSJed Brown return 0; 533c4762a1bSJed Brown } 534c4762a1bSJed Brown 535c4762a1bSJed Brown PetscErrorCode IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,void *ctx) 536c4762a1bSJed Brown { 537c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 538c4762a1bSJed Brown 539c4762a1bSJed Brown if (appctx->oshift == s) return 0; 5405f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(appctx->A,A,SAME_NONZERO_PATTERN)); 5415f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,-1)); 5425f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,s)); 5435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(A,B,SAME_NONZERO_PATTERN)); 544c4762a1bSJed Brown appctx->oshift = s; 545c4762a1bSJed Brown return 0; 546c4762a1bSJed Brown } 547c4762a1bSJed Brown 548c4762a1bSJed Brown /*TEST 549c4762a1bSJed Brown 550c4762a1bSJed Brown test: 551c4762a1bSJed Brown args: -nox -ts_type ssp -ts_dt 0.0005 552c4762a1bSJed Brown 553c4762a1bSJed Brown test: 554c4762a1bSJed Brown suffix: 2 555c4762a1bSJed Brown args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 3 559c4762a1bSJed Brown args: -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason 560c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 561c4762a1bSJed Brown requires: !single 562c4762a1bSJed Brown 563c4762a1bSJed Brown test: 564c4762a1bSJed Brown suffix: 4 565c4762a1bSJed Brown args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason 566c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 567c4762a1bSJed Brown 568c4762a1bSJed Brown test: 569c4762a1bSJed Brown suffix: 5 570c4762a1bSJed Brown args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs 571c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 572c4762a1bSJed Brown 573c4762a1bSJed Brown test: 574c4762a1bSJed Brown requires: !single 575c4762a1bSJed Brown suffix: pod_guess 576c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason 577c4762a1bSJed Brown 578c4762a1bSJed Brown test: 579c4762a1bSJed Brown requires: !single 580c4762a1bSJed Brown suffix: pod_guess_Ainner 581c4762a1bSJed 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 582c4762a1bSJed Brown 583c4762a1bSJed Brown test: 584c4762a1bSJed Brown requires: !single 585c4762a1bSJed Brown suffix: fischer_guess 586c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason 587c4762a1bSJed Brown 588c4762a1bSJed Brown test: 589c4762a1bSJed Brown requires: !single 590c4762a1bSJed Brown suffix: fischer_guess_2 591c4762a1bSJed 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 592c4762a1bSJed Brown 593c4762a1bSJed Brown test: 594c4762a1bSJed Brown requires: !single 595cbb17d71SDavid Wells suffix: fischer_guess_3 596cbb17d71SDavid 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 597cbb17d71SDavid Wells 598cbb17d71SDavid Wells test: 599cbb17d71SDavid Wells requires: !single 600c4762a1bSJed Brown suffix: stringview 601c4762a1bSJed Brown args: -nox -ts_type rosw -test_string_viewer 602c4762a1bSJed Brown 603c4762a1bSJed Brown test: 604c4762a1bSJed Brown requires: !single 605c4762a1bSJed Brown suffix: stringview_euler 606c4762a1bSJed Brown args: -nox -ts_type euler -test_string_viewer 607c4762a1bSJed Brown 608c4762a1bSJed Brown TEST*/ 609