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 Concepts: TS^time-dependent linear problems 11c4762a1bSJed Brown Concepts: TS^heat equation 12c4762a1bSJed Brown Concepts: TS^diffusion equation 13c4762a1bSJed Brown Processors: 1 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* ------------------------------------------------------------------------ 17c4762a1bSJed Brown 18c4762a1bSJed Brown This program solves the one-dimensional heat equation (also called the 19c4762a1bSJed Brown diffusion equation), 20c4762a1bSJed Brown u_t = u_xx, 21c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 22c4762a1bSJed Brown u(t,0) = 1, u(t,1) = 1, 23c4762a1bSJed Brown and the initial condition 24c4762a1bSJed Brown u(0,x) = cos(6*pi*x) + 3*cos(2*pi*x). 25c4762a1bSJed Brown This is a linear, second-order, parabolic equation. 26c4762a1bSJed Brown 27c4762a1bSJed Brown We discretize the right-hand side using finite differences with 28c4762a1bSJed Brown uniform grid spacing h: 29c4762a1bSJed Brown u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) 30c4762a1bSJed Brown We then demonstrate time evolution using the various TS methods by 31c4762a1bSJed Brown running the program via 32c4762a1bSJed Brown ex3 -ts_type <timestepping solver> 33c4762a1bSJed Brown 34c4762a1bSJed Brown We compare the approximate solution with the exact solution, given by 35c4762a1bSJed Brown u_exact(x,t) = exp(-36*pi*pi*t) * cos(6*pi*x) + 36c4762a1bSJed Brown 3*exp(-4*pi*pi*t) * cos(2*pi*x) 37c4762a1bSJed Brown 38c4762a1bSJed Brown Notes: 39c4762a1bSJed Brown This code demonstrates the TS solver interface to two variants of 40c4762a1bSJed Brown linear problems, u_t = f(u,t), namely 41c4762a1bSJed Brown - time-dependent f: f(u,t) is a function of t 42c4762a1bSJed Brown - time-independent f: f(u,t) is simply just f(u) 43c4762a1bSJed Brown 44c4762a1bSJed Brown The parallel version of this code is ts/tutorials/ex4.c 45c4762a1bSJed Brown 46c4762a1bSJed Brown ------------------------------------------------------------------------- */ 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* 49c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 50c4762a1bSJed Brown automatically includes: 51c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 52c4762a1bSJed Brown petscmat.h - matrices 53c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 54c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 55c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 56c4762a1bSJed Brown */ 57c4762a1bSJed Brown #include <petscts.h> 58c4762a1bSJed Brown #include <petscdraw.h> 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* 61c4762a1bSJed Brown User-defined application context - contains data needed by the 62c4762a1bSJed Brown application-provided call-back routines. 63c4762a1bSJed Brown */ 64c4762a1bSJed Brown typedef struct { 65c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 66c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 67c4762a1bSJed Brown PetscReal h; /* mesh width h = 1/(m-1) */ 68c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 69c4762a1bSJed Brown PetscViewer viewer1,viewer2; /* viewers for the solution and error */ 70c4762a1bSJed Brown PetscReal norm_2,norm_max; /* error norms */ 71c4762a1bSJed Brown } AppCtx; 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* 74c4762a1bSJed Brown User-defined routines 75c4762a1bSJed Brown */ 76c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 77c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*); 78c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 79c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); 80c4762a1bSJed Brown 81c4762a1bSJed Brown int main(int argc,char **argv) 82c4762a1bSJed Brown { 83c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 84c4762a1bSJed Brown TS ts; /* timestepping context */ 85c4762a1bSJed Brown Mat A; /* matrix data structure */ 86c4762a1bSJed Brown Vec u; /* approximate solution vector */ 87c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 88c4762a1bSJed Brown PetscInt time_steps_max = 100; /* default max timesteps */ 89c4762a1bSJed Brown PetscDraw draw; /* drawing context */ 90c4762a1bSJed Brown PetscInt steps,m; 91c4762a1bSJed Brown PetscMPIInt size; 92c4762a1bSJed Brown PetscBool flg; 93c4762a1bSJed Brown PetscReal dt,ftime; 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Initialize program and set problem parameters 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98c4762a1bSJed Brown 99*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 1005f80ce2aSJacob Faibussowitsch CHKERRMPI(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; 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug)); 106c4762a1bSJed Brown appctx.m = m; 107c4762a1bSJed Brown appctx.h = 1.0/(m-1.0); 108c4762a1bSJed Brown appctx.norm_2 = 0.0; 109c4762a1bSJed Brown appctx.norm_max = 0.0; 110c4762a1bSJed Brown 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n")); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Create vector data structures 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* 118c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 119c4762a1bSJed Brown */ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&u)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.solution)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Set up displays to show graphs of the solution and error 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126c4762a1bSJed Brown 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Create timestepping solver context 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137c4762a1bSJed Brown 1385f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_LINEAR)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142c4762a1bSJed Brown Set optional user-defined monitoring routine 143c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144c4762a1bSJed Brown 1455f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148c4762a1bSJed Brown 149c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 150c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 151c4762a1bSJed Brown 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 156c4762a1bSJed Brown 1575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-time_dependent_rhs",&flg)); 158c4762a1bSJed Brown if (flg) { 159c4762a1bSJed Brown /* 160c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 161c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 162c4762a1bSJed Brown as a time-dependent matrix. 163c4762a1bSJed Brown */ 1645f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx)); 166c4762a1bSJed Brown } else { 167c4762a1bSJed Brown /* 168c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 169c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 170c4762a1bSJed Brown as a matrix only once, and then sets a null matrix evaluation 171c4762a1bSJed Brown routine. 172c4762a1bSJed Brown */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(RHSMatrixHeat(ts,0.0,u,A,A,&appctx)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx)); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 179c4762a1bSJed Brown Set solution vector and initial timestep 180c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 181c4762a1bSJed Brown 182c4762a1bSJed Brown dt = appctx.h*appctx.h/2.0; 1835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,u)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 187c4762a1bSJed Brown Customize timestepping solver: 188c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 189c4762a1bSJed Brown - Set timestepping duration info 190c4762a1bSJed Brown Then set runtime options, which can override these defaults. 191c4762a1bSJed Brown For example, 192c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 193c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 194c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 195c4762a1bSJed Brown 1965f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,time_steps_max)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,time_total_max)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202c4762a1bSJed Brown Solve the problem 203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* 206c4762a1bSJed Brown Evaluate initial conditions 207c4762a1bSJed Brown */ 2085f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(u,&appctx)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* 211c4762a1bSJed Brown Run the timestepping solver 212c4762a1bSJed Brown */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,u)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 218c4762a1bSJed Brown View timestepping solver info 219c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 220c4762a1bSJed Brown 2215f80ce2aSJacob 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))); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(TSView(ts,PETSC_VIEWER_STDOUT_SELF)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 226c4762a1bSJed Brown are no longer needed. 227c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 228c4762a1bSJed Brown 2295f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&appctx.viewer1)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&appctx.viewer2)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.solution)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* 237c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 238c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 239c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 240c4762a1bSJed Brown options are chosen (e.g., -log_view). 241c4762a1bSJed Brown */ 242*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 243*b122ec5aSJacob Faibussowitsch return 0; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 246c4762a1bSJed Brown /* 247c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 248c4762a1bSJed Brown 249c4762a1bSJed Brown Input Parameter: 250c4762a1bSJed Brown u - uninitialized solution vector (global) 251c4762a1bSJed Brown appctx - user-defined application context 252c4762a1bSJed Brown 253c4762a1bSJed Brown Output Parameter: 254c4762a1bSJed Brown u - vector with solution at initial time (global) 255c4762a1bSJed Brown */ 256c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 257c4762a1bSJed Brown { 258c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h; 259c4762a1bSJed Brown PetscInt i; 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* 262c4762a1bSJed Brown Get a pointer to vector data. 263c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 264c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 265c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 266c4762a1bSJed Brown the array. 267c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 268c4762a1bSJed Brown C version. See the users manual for details. 269c4762a1bSJed Brown */ 2705f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u,&u_localptr)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown We initialize the solution array by simply writing the solution 274c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 275c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 276c4762a1bSJed Brown */ 277c4762a1bSJed Brown for (i=0; i<appctx->m; i++) u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h); 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown Restore vector 281c4762a1bSJed Brown */ 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u,&u_localptr)); 283c4762a1bSJed Brown 284c4762a1bSJed Brown /* 285c4762a1bSJed Brown Print debugging information if desired 286c4762a1bSJed Brown */ 287c4762a1bSJed Brown if (appctx->debug) { 288c4762a1bSJed Brown printf("initial guess vector\n"); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown return 0; 293c4762a1bSJed Brown } 294c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 295c4762a1bSJed Brown /* 296c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 297c4762a1bSJed Brown 298c4762a1bSJed Brown Input Parameters: 299c4762a1bSJed Brown t - current time 300c4762a1bSJed Brown solution - vector in which exact solution will be computed 301c4762a1bSJed Brown appctx - user-defined application context 302c4762a1bSJed Brown 303c4762a1bSJed Brown Output Parameter: 304c4762a1bSJed Brown solution - vector with the newly computed exact solution 305c4762a1bSJed Brown */ 306c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 307c4762a1bSJed Brown { 308c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t; 309c4762a1bSJed Brown PetscInt i; 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* 312c4762a1bSJed Brown Get a pointer to vector data. 313c4762a1bSJed Brown */ 3145f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(solution,&s_localptr)); 315c4762a1bSJed Brown 316c4762a1bSJed Brown /* 317c4762a1bSJed Brown Simply write the solution directly into the array locations. 318c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 319c4762a1bSJed Brown */ 320c4762a1bSJed Brown ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc); 321c4762a1bSJed Brown sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 322c4762a1bSJed Brown for (i=0; i<appctx->m; i++) s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2; 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* 325c4762a1bSJed Brown Restore vector 326c4762a1bSJed Brown */ 3275f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(solution,&s_localptr)); 328c4762a1bSJed Brown return 0; 329c4762a1bSJed Brown } 330c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 331c4762a1bSJed Brown /* 332c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 333c4762a1bSJed Brown each timestep. This example plots the solution and computes the 334c4762a1bSJed Brown error in two different norms. 335c4762a1bSJed Brown 336c4762a1bSJed Brown Input Parameters: 337c4762a1bSJed Brown ts - the timestep context 338c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 339c4762a1bSJed Brown initial condition) 340c4762a1bSJed Brown time - the current time 341c4762a1bSJed Brown u - the solution at this timestep 342c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 343c4762a1bSJed Brown In this case we use the application context which contains 344c4762a1bSJed Brown information about the problem size, workspace and the exact 345c4762a1bSJed Brown solution. 346c4762a1bSJed Brown */ 347c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 348c4762a1bSJed Brown { 349c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 350c4762a1bSJed Brown PetscReal norm_2,norm_max; 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* 353c4762a1bSJed Brown View a graph of the current iterate 354c4762a1bSJed Brown */ 3555f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,appctx->viewer2)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown /* 358c4762a1bSJed Brown Compute the exact solution 359c4762a1bSJed Brown */ 3605f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(time,appctx->solution,appctx)); 361c4762a1bSJed Brown 362c4762a1bSJed Brown /* 363c4762a1bSJed Brown Print debugging information if desired 364c4762a1bSJed Brown */ 365c4762a1bSJed Brown if (appctx->debug) { 366c4762a1bSJed Brown printf("Computed solution vector\n"); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 368c4762a1bSJed Brown printf("Exact solution vector\n"); 3695f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* 373c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 374c4762a1bSJed Brown */ 3755f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(appctx->solution,-1.0,u)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_2,&norm_2)); 377c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h)*norm_2; 3785f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&norm_max)); 379c4762a1bSJed Brown if (norm_2 < 1e-14) norm_2 = 0; 380c4762a1bSJed Brown if (norm_max < 1e-14) norm_max = 0; 381c4762a1bSJed Brown 3825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D: time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)time,(double)norm_2,(double)norm_max)); 383c4762a1bSJed Brown appctx->norm_2 += norm_2; 384c4762a1bSJed Brown appctx->norm_max += norm_max; 385c4762a1bSJed Brown 386c4762a1bSJed Brown /* 387c4762a1bSJed Brown View a graph of the error 388c4762a1bSJed Brown */ 3895f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,appctx->viewer1)); 390c4762a1bSJed Brown 391c4762a1bSJed Brown /* 392c4762a1bSJed Brown Print debugging information if desired 393c4762a1bSJed Brown */ 394c4762a1bSJed Brown if (appctx->debug) { 395c4762a1bSJed Brown printf("Error vector\n"); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 397c4762a1bSJed Brown } 398c4762a1bSJed Brown 399c4762a1bSJed Brown return 0; 400c4762a1bSJed Brown } 401c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 402c4762a1bSJed Brown /* 403c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 404c4762a1bSJed Brown matrix for the heat equation. 405c4762a1bSJed Brown 406c4762a1bSJed Brown Input Parameters: 407c4762a1bSJed Brown ts - the TS context 408c4762a1bSJed Brown t - current time 409c4762a1bSJed Brown global_in - global input vector 410c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 411c4762a1bSJed Brown 412c4762a1bSJed Brown Output Parameters: 413c4762a1bSJed Brown AA - Jacobian matrix 414c4762a1bSJed Brown BB - optionally different preconditioning matrix 415c4762a1bSJed Brown str - flag indicating matrix structure 416c4762a1bSJed Brown 417c4762a1bSJed Brown Notes: 418c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 419c4762a1bSJed Brown in Fortran as well as in C. 420c4762a1bSJed Brown */ 421c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx) 422c4762a1bSJed Brown { 423c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 424c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 425c4762a1bSJed Brown PetscInt mstart = 0; 426c4762a1bSJed Brown PetscInt mend = appctx->m; 427c4762a1bSJed Brown PetscInt i,idx[3]; 428c4762a1bSJed Brown PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; 429c4762a1bSJed Brown 430c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 431c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 432c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 433c4762a1bSJed Brown /* 434c4762a1bSJed Brown Set matrix rows corresponding to boundary data 435c4762a1bSJed Brown */ 436c4762a1bSJed Brown 437c4762a1bSJed Brown mstart = 0; 438c4762a1bSJed Brown v[0] = 1.0; 4395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES)); 440c4762a1bSJed Brown mstart++; 441c4762a1bSJed Brown 442c4762a1bSJed Brown mend--; 443c4762a1bSJed Brown v[0] = 1.0; 4445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES)); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* 447c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 448c4762a1bSJed Brown matrix one row at a time. 449c4762a1bSJed Brown */ 450c4762a1bSJed Brown v[0] = sone; v[1] = stwo; v[2] = sone; 451c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 452c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 4535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES)); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 457c4762a1bSJed Brown Complete the matrix assembly process and set some options 458c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 459c4762a1bSJed Brown /* 460c4762a1bSJed Brown Assemble matrix, using the 2-step process: 461c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 462c4762a1bSJed Brown Computations can be done while messages are in transition 463c4762a1bSJed Brown by placing code between these two statements. 464c4762a1bSJed Brown */ 4655f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 4665f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 467c4762a1bSJed Brown 468c4762a1bSJed Brown /* 469c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 470c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 471c4762a1bSJed Brown */ 4725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 473c4762a1bSJed Brown 474c4762a1bSJed Brown return 0; 475c4762a1bSJed Brown } 476c4762a1bSJed Brown 477c4762a1bSJed Brown /*TEST 478c4762a1bSJed Brown 479c4762a1bSJed Brown test: 480c4762a1bSJed Brown requires: x 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: nox 484c4762a1bSJed Brown args: -nox 485c4762a1bSJed Brown 486c4762a1bSJed Brown TEST*/ 487