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 PetscErrorCode ierr; 97c4762a1bSJed Brown PetscInt steps,m; 98c4762a1bSJed Brown PetscMPIInt size; 99c4762a1bSJed Brown PetscReal dt; 100c4762a1bSJed Brown PetscBool flg,flg_string; 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown Initialize program and set problem parameters 104c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 105c4762a1bSJed Brown 106c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 107ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 1082c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 109c4762a1bSJed Brown 110c4762a1bSJed Brown m = 60; 111c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);CHKERRQ(ierr); 113c4762a1bSJed Brown flg_string = PETSC_FALSE; 114c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_string_viewer",&flg_string,NULL);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown appctx.m = m; 117c4762a1bSJed Brown appctx.h = 1.0/(m-1.0); 118c4762a1bSJed Brown appctx.norm_2 = 0.0; 119c4762a1bSJed Brown appctx.norm_max = 0.0; 120c4762a1bSJed Brown 121c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");CHKERRQ(ierr); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Create vector data structures 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* 128c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 129c4762a1bSJed Brown */ 130c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,m,&u);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 134c4762a1bSJed Brown Set up displays to show graphs of the solution and error 135c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 136c4762a1bSJed Brown 137c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145c4762a1bSJed Brown Create timestepping solver context 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147c4762a1bSJed Brown 148c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_LINEAR);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152c4762a1bSJed Brown Set optional user-defined monitoring routine 153c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154c4762a1bSJed Brown 155c4762a1bSJed Brown if (!flg_string) { 156c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&appctx,NULL);CHKERRQ(ierr); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown 161c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 162c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 163c4762a1bSJed Brown 164c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 168c4762a1bSJed Brown 169c4762a1bSJed Brown flg = PETSC_FALSE; 170c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-use_ifunc",&flg,NULL);CHKERRQ(ierr); 171c4762a1bSJed Brown if (!flg) { 172c4762a1bSJed Brown appctx.A = NULL; 173c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL);CHKERRQ(ierr); 174c4762a1bSJed Brown if (flg) { 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 177c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 178c4762a1bSJed Brown as a time-dependent matrix. 179c4762a1bSJed Brown */ 180c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);CHKERRQ(ierr); 182c4762a1bSJed Brown } else { 183c4762a1bSJed Brown /* 184c4762a1bSJed Brown For linear problems with a time-independent f(u) in the equation 185c4762a1bSJed Brown u_t = f(u), the user provides the discretized right-hand-side 186c4762a1bSJed Brown as a matrix only once, and then sets the special Jacobian evaluation 187c4762a1bSJed Brown routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian. 188c4762a1bSJed Brown */ 189c4762a1bSJed Brown ierr = RHSMatrixHeat(ts,0.0,u,A,A,&appctx);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);CHKERRQ(ierr); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } else { 194c4762a1bSJed Brown Mat J; 195c4762a1bSJed Brown 196c4762a1bSJed Brown ierr = RHSMatrixHeat(ts,0.0,u,A,A,&appctx);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&J);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunctionHeat,&appctx);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = TSSetIJacobian(ts,J,J,IJacobianHeat,&appctx);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 201c4762a1bSJed Brown 202c4762a1bSJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 203c4762a1bSJed Brown appctx.A = A; 204c4762a1bSJed Brown appctx.oshift = PETSC_MIN_REAL; 205c4762a1bSJed Brown } 206c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown Set solution vector and initial timestep 208c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 209c4762a1bSJed Brown 210c4762a1bSJed Brown dt = appctx.h*appctx.h/2.0; 211c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 214c4762a1bSJed Brown Customize timestepping solver: 215c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 216c4762a1bSJed Brown - Set timestepping duration info 217c4762a1bSJed Brown Then set runtime options, which can override these defaults. 218c4762a1bSJed Brown For example, 219c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 220c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 221c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 222c4762a1bSJed Brown 223c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,time_steps_max);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = TSSetMaxTime(ts,time_total_max);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 229c4762a1bSJed Brown Solve the problem 230c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* 233c4762a1bSJed Brown Evaluate initial conditions 234c4762a1bSJed Brown */ 235c4762a1bSJed Brown ierr = InitialConditions(u,&appctx);CHKERRQ(ierr); 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown Run the timestepping solver 239c4762a1bSJed Brown */ 240c4762a1bSJed Brown ierr = TSSolve(ts,u);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 244c4762a1bSJed Brown View timestepping solver info 245c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 246c4762a1bSJed Brown 247c4762a1bSJed 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); 248c4762a1bSJed Brown if (!flg_string) { 249c4762a1bSJed Brown ierr = TSView(ts,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 250c4762a1bSJed Brown } else { 251c4762a1bSJed Brown PetscViewer stringviewer; 252c4762a1bSJed Brown char string[512]; 253c4762a1bSJed Brown const char *outstring; 254c4762a1bSJed Brown 255c4762a1bSJed Brown ierr = PetscViewerStringOpen(PETSC_COMM_WORLD,string,sizeof(string),&stringviewer);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = TSView(ts,stringviewer);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = PetscViewerStringGetStringRead(stringviewer,&outstring,NULL);CHKERRQ(ierr); 2582c71b3e2SJacob Faibussowitsch PetscCheckFalse((char*)outstring != (char*)string,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"String returned from viewer does not equal original string"); 259c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Output from string viewer:%s\n",outstring);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = PetscViewerDestroy(&stringviewer);CHKERRQ(ierr); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 265c4762a1bSJed Brown are no longer needed. 266c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 267c4762a1bSJed Brown 268c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = PetscViewerDestroy(&appctx.viewer1);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = PetscViewerDestroy(&appctx.viewer2);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = MatDestroy(&appctx.A);CHKERRQ(ierr); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* 277c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 278c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 279c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 280c4762a1bSJed Brown options are chosen (e.g., -log_view). 281c4762a1bSJed Brown */ 282c4762a1bSJed Brown ierr = PetscFinalize(); 283c4762a1bSJed Brown return ierr; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 286c4762a1bSJed Brown /* 287c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 288c4762a1bSJed Brown 289c4762a1bSJed Brown Input Parameter: 290c4762a1bSJed Brown u - uninitialized solution vector (global) 291c4762a1bSJed Brown appctx - user-defined application context 292c4762a1bSJed Brown 293c4762a1bSJed Brown Output Parameter: 294c4762a1bSJed Brown u - vector with solution at initial time (global) 295c4762a1bSJed Brown */ 296c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 297c4762a1bSJed Brown { 298c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h; 299c4762a1bSJed Brown PetscErrorCode ierr; 300c4762a1bSJed Brown PetscInt i; 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* 303c4762a1bSJed Brown Get a pointer to vector data. 304c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 305c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 306c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 307c4762a1bSJed Brown the array. 308c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 309c4762a1bSJed Brown C version. See the users manual for details. 310c4762a1bSJed Brown */ 311303a5415SBarry Smith ierr = VecGetArrayWrite(u,&u_localptr);CHKERRQ(ierr); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown We initialize the solution array by simply writing the solution 315c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 316c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 317c4762a1bSJed Brown */ 318c4762a1bSJed Brown for (i=0; i<appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h); 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown Restore vector 322c4762a1bSJed Brown */ 323303a5415SBarry Smith ierr = VecRestoreArrayWrite(u,&u_localptr);CHKERRQ(ierr); 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* 326c4762a1bSJed Brown Print debugging information if desired 327c4762a1bSJed Brown */ 328c4762a1bSJed Brown if (appctx->debug) { 329c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n");CHKERRQ(ierr); 330c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown return 0; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 336c4762a1bSJed Brown /* 337c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 338c4762a1bSJed Brown 339c4762a1bSJed Brown Input Parameters: 340c4762a1bSJed Brown t - current time 341c4762a1bSJed Brown solution - vector in which exact solution will be computed 342c4762a1bSJed Brown appctx - user-defined application context 343c4762a1bSJed Brown 344c4762a1bSJed Brown Output Parameter: 345c4762a1bSJed Brown solution - vector with the newly computed exact solution 346c4762a1bSJed Brown */ 347c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 348c4762a1bSJed Brown { 349c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t; 350c4762a1bSJed Brown PetscErrorCode ierr; 351c4762a1bSJed Brown PetscInt i; 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* 354c4762a1bSJed Brown Get a pointer to vector data. 355c4762a1bSJed Brown */ 356303a5415SBarry Smith ierr = VecGetArrayWrite(solution,&s_localptr);CHKERRQ(ierr); 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* 359c4762a1bSJed Brown Simply write the solution directly into the array locations. 360c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 361c4762a1bSJed Brown */ 362c4762a1bSJed Brown ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); 363c4762a1bSJed Brown ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc); 364c4762a1bSJed Brown sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 365c4762a1bSJed Brown for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2; 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* 368c4762a1bSJed Brown Restore vector 369c4762a1bSJed Brown */ 370303a5415SBarry Smith ierr = VecRestoreArrayWrite(solution,&s_localptr);CHKERRQ(ierr); 371c4762a1bSJed Brown return 0; 372c4762a1bSJed Brown } 373c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 374c4762a1bSJed Brown /* 375c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 376c4762a1bSJed Brown each timestep. This example plots the solution and computes the 377c4762a1bSJed Brown error in two different norms. 378c4762a1bSJed Brown 379c4762a1bSJed Brown This example also demonstrates changing the timestep via TSSetTimeStep(). 380c4762a1bSJed Brown 381c4762a1bSJed Brown Input Parameters: 382c4762a1bSJed Brown ts - the timestep context 383c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 384c4762a1bSJed Brown initial condition) 385c4762a1bSJed Brown time - the current time 386c4762a1bSJed Brown u - the solution at this timestep 387c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 388c4762a1bSJed Brown In this case we use the application context which contains 389c4762a1bSJed Brown information about the problem size, workspace and the exact 390c4762a1bSJed Brown solution. 391c4762a1bSJed Brown */ 392c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 393c4762a1bSJed Brown { 394c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 395c4762a1bSJed Brown PetscErrorCode ierr; 396c4762a1bSJed Brown PetscReal norm_2,norm_max,dt,dttol; 397c4762a1bSJed Brown 398c4762a1bSJed Brown /* 399c4762a1bSJed Brown View a graph of the current iterate 400c4762a1bSJed Brown */ 401c4762a1bSJed Brown ierr = VecView(u,appctx->viewer2);CHKERRQ(ierr); 402c4762a1bSJed Brown 403c4762a1bSJed Brown /* 404c4762a1bSJed Brown Compute the exact solution 405c4762a1bSJed Brown */ 406c4762a1bSJed Brown ierr = ExactSolution(time,appctx->solution,appctx);CHKERRQ(ierr); 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* 409c4762a1bSJed Brown Print debugging information if desired 410c4762a1bSJed Brown */ 411c4762a1bSJed Brown if (appctx->debug) { 412c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n");CHKERRQ(ierr); 413c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 414c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");CHKERRQ(ierr); 415c4762a1bSJed Brown ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 418c4762a1bSJed Brown /* 419c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 420c4762a1bSJed Brown */ 421c4762a1bSJed Brown ierr = VecAXPY(appctx->solution,-1.0,u);CHKERRQ(ierr); 422c4762a1bSJed Brown ierr = VecNorm(appctx->solution,NORM_2,&norm_2);CHKERRQ(ierr); 423c4762a1bSJed Brown norm_2 = PetscSqrtReal(appctx->h)*norm_2; 424c4762a1bSJed Brown ierr = VecNorm(appctx->solution,NORM_MAX,&norm_max);CHKERRQ(ierr); 425c4762a1bSJed Brown 426c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 427c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 428c4762a1bSJed Brown 429c4762a1bSJed Brown appctx->norm_2 += norm_2; 430c4762a1bSJed Brown appctx->norm_max += norm_max; 431c4762a1bSJed Brown 432c4762a1bSJed Brown dttol = .0001; 433c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL);CHKERRQ(ierr); 434c4762a1bSJed Brown if (dt < dttol) { 435c4762a1bSJed Brown dt *= .999; 436c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* 440c4762a1bSJed Brown View a graph of the error 441c4762a1bSJed Brown */ 442c4762a1bSJed Brown ierr = VecView(appctx->solution,appctx->viewer1);CHKERRQ(ierr); 443c4762a1bSJed Brown 444c4762a1bSJed Brown /* 445c4762a1bSJed Brown Print debugging information if desired 446c4762a1bSJed Brown */ 447c4762a1bSJed Brown if (appctx->debug) { 448c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error vector\n");CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452c4762a1bSJed Brown return 0; 453c4762a1bSJed Brown } 454c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 455c4762a1bSJed Brown /* 456c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 457c4762a1bSJed Brown matrix for the heat equation. 458c4762a1bSJed Brown 459c4762a1bSJed Brown Input Parameters: 460c4762a1bSJed Brown ts - the TS context 461c4762a1bSJed Brown t - current time 462c4762a1bSJed Brown global_in - global input vector 463c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 464c4762a1bSJed Brown 465c4762a1bSJed Brown Output Parameters: 466c4762a1bSJed Brown AA - Jacobian matrix 467c4762a1bSJed Brown BB - optionally different preconditioning matrix 468c4762a1bSJed Brown str - flag indicating matrix structure 469c4762a1bSJed Brown 470c4762a1bSJed Brown Notes: 471c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 472c4762a1bSJed Brown in Fortran as well as in C. 473c4762a1bSJed Brown */ 474c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx) 475c4762a1bSJed Brown { 476c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 477c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 478c4762a1bSJed Brown PetscInt mstart = 0; 479c4762a1bSJed Brown PetscInt mend = appctx->m; 480c4762a1bSJed Brown PetscErrorCode ierr; 481c4762a1bSJed Brown PetscInt i,idx[3]; 482c4762a1bSJed Brown PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; 483c4762a1bSJed Brown 484c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 485c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 486c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 487c4762a1bSJed Brown /* 488c4762a1bSJed Brown Set matrix rows corresponding to boundary data 489c4762a1bSJed Brown */ 490c4762a1bSJed Brown 491c4762a1bSJed Brown mstart = 0; 492c4762a1bSJed Brown v[0] = 1.0; 493c4762a1bSJed Brown ierr = MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr); 494c4762a1bSJed Brown mstart++; 495c4762a1bSJed Brown 496c4762a1bSJed Brown mend--; 497c4762a1bSJed Brown v[0] = 1.0; 498c4762a1bSJed Brown ierr = MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr); 499c4762a1bSJed Brown 500c4762a1bSJed Brown /* 501c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 502c4762a1bSJed Brown matrix one row at a time. 503c4762a1bSJed Brown */ 504c4762a1bSJed Brown v[0] = sone; v[1] = stwo; v[2] = sone; 505c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 506c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 507c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 508c4762a1bSJed Brown } 509c4762a1bSJed Brown 510c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 511c4762a1bSJed Brown Complete the matrix assembly process and set some options 512c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 513c4762a1bSJed Brown /* 514c4762a1bSJed Brown Assemble matrix, using the 2-step process: 515c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 516c4762a1bSJed Brown Computations can be done while messages are in transition 517c4762a1bSJed Brown by placing code between these two statements. 518c4762a1bSJed Brown */ 519c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 520c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 521c4762a1bSJed Brown 522c4762a1bSJed Brown /* 523c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 524c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 525c4762a1bSJed Brown */ 526c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 527c4762a1bSJed Brown 528c4762a1bSJed Brown return 0; 529c4762a1bSJed Brown } 530c4762a1bSJed Brown 531c4762a1bSJed Brown PetscErrorCode IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,void *ctx) 532c4762a1bSJed Brown { 533c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 534c4762a1bSJed Brown PetscErrorCode ierr; 535c4762a1bSJed Brown 536c4762a1bSJed Brown ierr = MatMult(appctx->A,X,r);CHKERRQ(ierr); 537c4762a1bSJed Brown ierr = VecAYPX(r,-1.0,Xdot);CHKERRQ(ierr); 538c4762a1bSJed Brown return 0; 539c4762a1bSJed Brown } 540c4762a1bSJed Brown 541c4762a1bSJed Brown PetscErrorCode IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,void *ctx) 542c4762a1bSJed Brown { 543c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 544c4762a1bSJed Brown PetscErrorCode ierr; 545c4762a1bSJed Brown 546c4762a1bSJed Brown if (appctx->oshift == s) return 0; 547c4762a1bSJed Brown ierr = MatCopy(appctx->A,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 548c4762a1bSJed Brown ierr = MatScale(A,-1);CHKERRQ(ierr); 549c4762a1bSJed Brown ierr = MatShift(A,s);CHKERRQ(ierr); 550c4762a1bSJed Brown ierr = MatCopy(A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 551c4762a1bSJed Brown appctx->oshift = s; 552c4762a1bSJed Brown return 0; 553c4762a1bSJed Brown } 554c4762a1bSJed Brown 555c4762a1bSJed Brown /*TEST 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown args: -nox -ts_type ssp -ts_dt 0.0005 559c4762a1bSJed Brown 560c4762a1bSJed Brown test: 561c4762a1bSJed Brown suffix: 2 562c4762a1bSJed Brown args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1 563c4762a1bSJed Brown 564c4762a1bSJed Brown test: 565c4762a1bSJed Brown suffix: 3 566c4762a1bSJed Brown args: -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason 567c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 568c4762a1bSJed Brown requires: !single 569c4762a1bSJed Brown 570c4762a1bSJed Brown test: 571c4762a1bSJed Brown suffix: 4 572c4762a1bSJed Brown args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason 573c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 574c4762a1bSJed Brown 575c4762a1bSJed Brown test: 576c4762a1bSJed Brown suffix: 5 577c4762a1bSJed Brown args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs 578c4762a1bSJed Brown filter: sed "s/ATOL/RTOL/g" 579c4762a1bSJed Brown 580c4762a1bSJed Brown test: 581c4762a1bSJed Brown requires: !single 582c4762a1bSJed Brown suffix: pod_guess 583c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason 584c4762a1bSJed Brown 585c4762a1bSJed Brown test: 586c4762a1bSJed Brown requires: !single 587c4762a1bSJed Brown suffix: pod_guess_Ainner 588c4762a1bSJed 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 589c4762a1bSJed Brown 590c4762a1bSJed Brown test: 591c4762a1bSJed Brown requires: !single 592c4762a1bSJed Brown suffix: fischer_guess 593c4762a1bSJed Brown args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason 594c4762a1bSJed Brown 595c4762a1bSJed Brown test: 596c4762a1bSJed Brown requires: !single 597c4762a1bSJed Brown suffix: fischer_guess_2 598c4762a1bSJed 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 599c4762a1bSJed Brown 600c4762a1bSJed Brown test: 601c4762a1bSJed Brown requires: !single 602*cbb17d71SDavid Wells suffix: fischer_guess_3 603*cbb17d71SDavid 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 604*cbb17d71SDavid Wells 605*cbb17d71SDavid Wells test: 606*cbb17d71SDavid Wells requires: !single 607c4762a1bSJed Brown suffix: stringview 608c4762a1bSJed Brown args: -nox -ts_type rosw -test_string_viewer 609c4762a1bSJed Brown 610c4762a1bSJed Brown test: 611c4762a1bSJed Brown requires: !single 612c4762a1bSJed Brown suffix: stringview_euler 613c4762a1bSJed Brown args: -nox -ts_type euler -test_string_viewer 614c4762a1bSJed Brown 615c4762a1bSJed Brown TEST*/ 616