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 PetscErrorCode ierr; 91c4762a1bSJed Brown PetscInt steps,m; 92c4762a1bSJed Brown PetscMPIInt size; 93c4762a1bSJed Brown PetscBool flg; 94c4762a1bSJed Brown PetscReal dt,ftime; 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 97c4762a1bSJed Brown Initialize program and set problem parameters 98c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 99c4762a1bSJed Brown 100c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 101ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 102*3c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 103c4762a1bSJed Brown 104c4762a1bSJed Brown m = 60; 105c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);CHKERRQ(ierr); 107c4762a1bSJed Brown appctx.m = m; 108c4762a1bSJed Brown appctx.h = 1.0/(m-1.0); 109c4762a1bSJed Brown appctx.norm_2 = 0.0; 110c4762a1bSJed Brown appctx.norm_max = 0.0; 111c4762a1bSJed Brown 112c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");CHKERRQ(ierr); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 115c4762a1bSJed Brown Create vector data structures 116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* 119c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 120c4762a1bSJed Brown */ 121c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,m,&u);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Set up displays to show graphs of the solution and error 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127c4762a1bSJed Brown 128c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Create timestepping solver context 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4762a1bSJed Brown 139c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_LINEAR);CHKERRQ(ierr); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Set optional user-defined monitoring routine 144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 145c4762a1bSJed Brown 146c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&appctx,NULL);CHKERRQ(ierr); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 149c4762a1bSJed Brown 150c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 152c4762a1bSJed Brown 153c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 157c4762a1bSJed Brown 158c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-time_dependent_rhs",&flg);CHKERRQ(ierr); 159c4762a1bSJed Brown if (flg) { 160c4762a1bSJed Brown /* 161c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 162c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 163c4762a1bSJed Brown as a time-dependent matrix. 164c4762a1bSJed Brown */ 165c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);CHKERRQ(ierr); 167c4762a1bSJed Brown } else { 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 170c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 171c4762a1bSJed Brown as a matrix only once, and then sets a null matrix evaluation 172c4762a1bSJed Brown routine. 173c4762a1bSJed Brown */ 174c4762a1bSJed Brown ierr = RHSMatrixHeat(ts,0.0,u,A,A,&appctx);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);CHKERRQ(ierr); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180c4762a1bSJed Brown Set solution vector and initial timestep 181c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 182c4762a1bSJed Brown 183c4762a1bSJed Brown dt = appctx.h*appctx.h/2.0; 184c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = TSSetSolution(ts,u);CHKERRQ(ierr); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188c4762a1bSJed Brown Customize timestepping solver: 189c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 190c4762a1bSJed Brown - Set timestepping duration info 191c4762a1bSJed Brown Then set runtime options, which can override these defaults. 192c4762a1bSJed Brown For example, 193c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 194c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 195c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 196c4762a1bSJed Brown 197c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,time_steps_max);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = TSSetMaxTime(ts,time_total_max);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 203c4762a1bSJed Brown Solve the problem 204c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Evaluate initial conditions 208c4762a1bSJed Brown */ 209c4762a1bSJed Brown ierr = InitialConditions(u,&appctx);CHKERRQ(ierr); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Run the timestepping solver 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown ierr = TSSolve(ts,u);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 219c4762a1bSJed Brown View timestepping solver info 220c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 221c4762a1bSJed Brown 222c4762a1bSJed Brown ierr = 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));CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = TSView(ts,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 226c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 227c4762a1bSJed Brown are no longer needed. 228c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 229c4762a1bSJed Brown 230c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = PetscViewerDestroy(&appctx.viewer1);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = PetscViewerDestroy(&appctx.viewer2);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr); 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 239c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 240c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 241c4762a1bSJed Brown options are chosen (e.g., -log_view). 242c4762a1bSJed Brown */ 243c4762a1bSJed Brown ierr = PetscFinalize(); 244c4762a1bSJed Brown return ierr; 245c4762a1bSJed Brown } 246c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 247c4762a1bSJed Brown /* 248c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 249c4762a1bSJed Brown 250c4762a1bSJed Brown Input Parameter: 251c4762a1bSJed Brown u - uninitialized solution vector (global) 252c4762a1bSJed Brown appctx - user-defined application context 253c4762a1bSJed Brown 254c4762a1bSJed Brown Output Parameter: 255c4762a1bSJed Brown u - vector with solution at initial time (global) 256c4762a1bSJed Brown */ 257c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 258c4762a1bSJed Brown { 259c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h; 260c4762a1bSJed Brown PetscInt i; 261c4762a1bSJed Brown PetscErrorCode ierr; 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* 264c4762a1bSJed Brown Get a pointer to vector data. 265c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 266c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 267c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 268c4762a1bSJed Brown the array. 269c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 270c4762a1bSJed Brown C version. See the users manual for details. 271c4762a1bSJed Brown */ 272c4762a1bSJed Brown ierr = VecGetArray(u,&u_localptr);CHKERRQ(ierr); 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* 275c4762a1bSJed Brown We initialize the solution array by simply writing the solution 276c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 277c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 278c4762a1bSJed Brown */ 279c4762a1bSJed Brown for (i=0; i<appctx->m; i++) u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h); 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* 282c4762a1bSJed Brown Restore vector 283c4762a1bSJed Brown */ 284c4762a1bSJed Brown ierr = VecRestoreArray(u,&u_localptr);CHKERRQ(ierr); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* 287c4762a1bSJed Brown Print debugging information if desired 288c4762a1bSJed Brown */ 289c4762a1bSJed Brown if (appctx->debug) { 290c4762a1bSJed Brown printf("initial guess vector\n"); 291c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 292c4762a1bSJed Brown } 293c4762a1bSJed Brown 294c4762a1bSJed Brown return 0; 295c4762a1bSJed Brown } 296c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 299c4762a1bSJed Brown 300c4762a1bSJed Brown Input Parameters: 301c4762a1bSJed Brown t - current time 302c4762a1bSJed Brown solution - vector in which exact solution will be computed 303c4762a1bSJed Brown appctx - user-defined application context 304c4762a1bSJed Brown 305c4762a1bSJed Brown Output Parameter: 306c4762a1bSJed Brown solution - vector with the newly computed exact solution 307c4762a1bSJed Brown */ 308c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 309c4762a1bSJed Brown { 310c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t; 311c4762a1bSJed Brown PetscInt i; 312c4762a1bSJed Brown PetscErrorCode ierr; 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* 315c4762a1bSJed Brown Get a pointer to vector data. 316c4762a1bSJed Brown */ 317c4762a1bSJed Brown ierr = VecGetArray(solution,&s_localptr);CHKERRQ(ierr); 318c4762a1bSJed Brown 319c4762a1bSJed Brown /* 320c4762a1bSJed Brown Simply write the solution directly into the array locations. 321c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 322c4762a1bSJed Brown */ 323c4762a1bSJed Brown ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc); 324c4762a1bSJed Brown sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 325c4762a1bSJed Brown for (i=0; i<appctx->m; i++) s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2; 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* 328c4762a1bSJed Brown Restore vector 329c4762a1bSJed Brown */ 330c4762a1bSJed Brown ierr = VecRestoreArray(solution,&s_localptr);CHKERRQ(ierr); 331c4762a1bSJed Brown return 0; 332c4762a1bSJed Brown } 333c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 334c4762a1bSJed Brown /* 335c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 336c4762a1bSJed Brown each timestep. This example plots the solution and computes the 337c4762a1bSJed Brown error in two different norms. 338c4762a1bSJed Brown 339c4762a1bSJed Brown Input Parameters: 340c4762a1bSJed Brown ts - the timestep context 341c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 342c4762a1bSJed Brown initial condition) 343c4762a1bSJed Brown time - the current time 344c4762a1bSJed Brown u - the solution at this timestep 345c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 346c4762a1bSJed Brown In this case we use the application context which contains 347c4762a1bSJed Brown information about the problem size, workspace and the exact 348c4762a1bSJed Brown solution. 349c4762a1bSJed Brown */ 350c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 351c4762a1bSJed Brown { 352c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 353c4762a1bSJed Brown PetscErrorCode ierr; 354c4762a1bSJed Brown PetscReal norm_2,norm_max; 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown View a graph of the current iterate 358c4762a1bSJed Brown */ 359c4762a1bSJed Brown ierr = VecView(u,appctx->viewer2);CHKERRQ(ierr); 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* 362c4762a1bSJed Brown Compute the exact solution 363c4762a1bSJed Brown */ 364c4762a1bSJed Brown ierr = ExactSolution(time,appctx->solution,appctx);CHKERRQ(ierr); 365c4762a1bSJed Brown 366c4762a1bSJed Brown /* 367c4762a1bSJed Brown Print debugging information if desired 368c4762a1bSJed Brown */ 369c4762a1bSJed Brown if (appctx->debug) { 370c4762a1bSJed Brown printf("Computed solution vector\n"); 371c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 372c4762a1bSJed Brown printf("Exact solution vector\n"); 373c4762a1bSJed Brown ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 374c4762a1bSJed Brown } 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* 377c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 378c4762a1bSJed Brown */ 379c4762a1bSJed Brown ierr = VecAXPY(appctx->solution,-1.0,u);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = VecNorm(appctx->solution,NORM_2,&norm_2);CHKERRQ(ierr); 381c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h)*norm_2; 382c4762a1bSJed Brown ierr = VecNorm(appctx->solution,NORM_MAX,&norm_max);CHKERRQ(ierr); 383c4762a1bSJed Brown if (norm_2 < 1e-14) norm_2 = 0; 384c4762a1bSJed Brown if (norm_max < 1e-14) norm_max = 0; 385c4762a1bSJed Brown 386c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 387c4762a1bSJed Brown appctx->norm_2 += norm_2; 388c4762a1bSJed Brown appctx->norm_max += norm_max; 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* 391c4762a1bSJed Brown View a graph of the error 392c4762a1bSJed Brown */ 393c4762a1bSJed Brown ierr = VecView(appctx->solution,appctx->viewer1);CHKERRQ(ierr); 394c4762a1bSJed Brown 395c4762a1bSJed Brown /* 396c4762a1bSJed Brown Print debugging information if desired 397c4762a1bSJed Brown */ 398c4762a1bSJed Brown if (appctx->debug) { 399c4762a1bSJed Brown printf("Error vector\n"); 400c4762a1bSJed Brown ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 401c4762a1bSJed Brown } 402c4762a1bSJed Brown 403c4762a1bSJed Brown return 0; 404c4762a1bSJed Brown } 405c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 406c4762a1bSJed Brown /* 407c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 408c4762a1bSJed Brown matrix for the heat equation. 409c4762a1bSJed Brown 410c4762a1bSJed Brown Input Parameters: 411c4762a1bSJed Brown ts - the TS context 412c4762a1bSJed Brown t - current time 413c4762a1bSJed Brown global_in - global input vector 414c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 415c4762a1bSJed Brown 416c4762a1bSJed Brown Output Parameters: 417c4762a1bSJed Brown AA - Jacobian matrix 418c4762a1bSJed Brown BB - optionally different preconditioning matrix 419c4762a1bSJed Brown str - flag indicating matrix structure 420c4762a1bSJed Brown 421c4762a1bSJed Brown Notes: 422c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 423c4762a1bSJed Brown in Fortran as well as in C. 424c4762a1bSJed Brown */ 425c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx) 426c4762a1bSJed Brown { 427c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 428c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 429c4762a1bSJed Brown PetscInt mstart = 0; 430c4762a1bSJed Brown PetscInt mend = appctx->m; 431c4762a1bSJed Brown PetscErrorCode ierr; 432c4762a1bSJed Brown PetscInt i,idx[3]; 433c4762a1bSJed Brown PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 436c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 437c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 438c4762a1bSJed Brown /* 439c4762a1bSJed Brown Set matrix rows corresponding to boundary data 440c4762a1bSJed Brown */ 441c4762a1bSJed Brown 442c4762a1bSJed Brown mstart = 0; 443c4762a1bSJed Brown v[0] = 1.0; 444c4762a1bSJed Brown ierr = MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr); 445c4762a1bSJed Brown mstart++; 446c4762a1bSJed Brown 447c4762a1bSJed Brown mend--; 448c4762a1bSJed Brown v[0] = 1.0; 449c4762a1bSJed Brown ierr = MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr); 450c4762a1bSJed Brown 451c4762a1bSJed Brown /* 452c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 453c4762a1bSJed Brown matrix one row at a time. 454c4762a1bSJed Brown */ 455c4762a1bSJed Brown v[0] = sone; v[1] = stwo; v[2] = sone; 456c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 457c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 458c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 459c4762a1bSJed Brown } 460c4762a1bSJed Brown 461c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 462c4762a1bSJed Brown Complete the matrix assembly process and set some options 463c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 464c4762a1bSJed Brown /* 465c4762a1bSJed Brown Assemble matrix, using the 2-step process: 466c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 467c4762a1bSJed Brown Computations can be done while messages are in transition 468c4762a1bSJed Brown by placing code between these two statements. 469c4762a1bSJed Brown */ 470c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 471c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 472c4762a1bSJed Brown 473c4762a1bSJed Brown /* 474c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 475c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 476c4762a1bSJed Brown */ 477c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 478c4762a1bSJed Brown 479c4762a1bSJed Brown return 0; 480c4762a1bSJed Brown } 481c4762a1bSJed Brown 482c4762a1bSJed Brown /*TEST 483c4762a1bSJed Brown 484c4762a1bSJed Brown test: 485c4762a1bSJed Brown requires: x 486c4762a1bSJed Brown 487c4762a1bSJed Brown test: 488c4762a1bSJed Brown suffix: nox 489c4762a1bSJed Brown args: -nox 490c4762a1bSJed Brown 491c4762a1bSJed Brown TEST*/ 492c4762a1bSJed Brown 493