/* Program usage: ex3 [-help] [all PETSc options] */ static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\ Input parameters include:\n\ -m , where = number of grid points\n\ -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\ -debug : Activate debugging printouts\n\ -nox : Deactivate x-window graphics\n\n"; /* Concepts: TS^time-dependent linear problems Concepts: TS^heat equation Concepts: TS^diffusion equation Processors: 1 */ /* ------------------------------------------------------------------------ This program solves the one-dimensional heat equation (also called the diffusion equation), u_t = u_xx, on the domain 0 <= x <= 1, with the boundary conditions u(t,0) = 1, u(t,1) = 1, and the initial condition u(0,x) = cos(6*pi*x) + 3*cos(2*pi*x). This is a linear, second-order, parabolic equation. We discretize the right-hand side using finite differences with uniform grid spacing h: u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2) We then demonstrate time evolution using the various TS methods by running the program via ex3 -ts_type We compare the approximate solution with the exact solution, given by u_exact(x,t) = exp(-36*pi*pi*t) * cos(6*pi*x) + 3*exp(-4*pi*pi*t) * cos(2*pi*x) Notes: This code demonstrates the TS solver interface to two variants of linear problems, u_t = f(u,t), namely - time-dependent f: f(u,t) is a function of t - time-independent f: f(u,t) is simply just f(u) The parallel version of this code is ts/examples/tutorials/ex4.c ------------------------------------------------------------------------- */ /* Include "petscts.h" so that we can use TS solvers. Note that this file automatically includes: petsc.h - base PETSc routines petscvec.h - vectors petscsys.h - system routines petscmat.h - matrices petscis.h - index sets petscksp.h - Krylov subspace methods petscviewer.h - viewers petscpc.h - preconditioners petscksp.h - linear solvers petscsnes.h - nonlinear solvers */ #include "petscts.h" /* User-defined application context - contains data needed by the application-provided call-back routines. */ typedef struct { Vec solution; /* global exact solution vector */ PetscInt m; /* total number of grid points */ PetscReal h; /* mesh width h = 1/(m-1) */ PetscTruth debug; /* flag (1 indicates activation of debugging printouts) */ PetscViewer viewer1,viewer2; /* viewers for the solution and error */ PetscReal norm_2,norm_max; /* error norms */ } AppCtx; /* User-defined routines */ extern PetscErrorCode InitialConditions(Vec,AppCtx*); extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Mat*,Mat*,MatStructure*,void*); extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **argv) { AppCtx appctx; /* user-defined application context */ TS ts; /* timestepping context */ Mat A; /* matrix data structure */ Vec u; /* approximate solution vector */ PetscReal time_total_max = 100.0; /* default max total time */ PetscInt time_steps_max = 100; /* default max timesteps */ PetscDraw draw; /* drawing context */ PetscErrorCode ierr; PetscInt steps,m; PetscMPIInt size; PetscTruth flg; PetscReal dt,ftime; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize program and set problem parameters - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscInitialize(&argc,&argv,(char*)0,help);CHKERRQ(ierr); ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); if (size != 1) SETERRQ(1,"This is a uniprocessor example only!"); m = 60; ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsHasName(PETSC_NULL,"-debug",&appctx.debug);CHKERRQ(ierr); appctx.m = m; appctx.h = 1.0/(m-1.0); appctx.norm_2 = 0.0; appctx.norm_max = 0.0; ierr = PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n");CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create vector data structures - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Create vector data structures for approximate and exact solutions */ ierr = VecCreateSeq(PETSC_COMM_SELF,m,&u);CHKERRQ(ierr); ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set up displays to show graphs of the solution and error - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1);CHKERRQ(ierr); ierr = PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);CHKERRQ(ierr); ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2);CHKERRQ(ierr); ierr = PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);CHKERRQ(ierr); ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create timestepping solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr); ierr = TSSetProblemType(ts,TS_LINEAR);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set optional user-defined monitoring routine - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSetMonitor(ts,Monitor,&appctx,PETSC_NULL);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create matrix data structure; set matrix evaluation routine. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = PetscOptionsHasName(PETSC_NULL,"-time_dependent_rhs",&flg);CHKERRQ(ierr); if (flg) { /* For linear problems with a time-dependent f(u,t) in the equation u_t = f(u,t), the user provides the discretized right-hand-side as a time-dependent matrix. */ ierr = TSSetRHSMatrix(ts,A,A,RHSMatrixHeat,&appctx);CHKERRQ(ierr); } else { /* For linear problems with a time-independent f(u) in the equation u_t = f(u), the user provides the discretized right-hand-side as a matrix only once, and then sets a null matrix evaluation routine. */ MatStructure A_structure; ierr = RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);CHKERRQ(ierr); ierr = TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);CHKERRQ(ierr); } /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set solution vector and initial timestep - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ dt = appctx.h*appctx.h/2.0; ierr = TSSetInitialTimeStep(ts,0.0,dt);CHKERRQ(ierr); ierr = TSSetSolution(ts,u);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Customize timestepping solver: - Set the solution method to be the Backward Euler method. - Set timestepping duration info Then set runtime options, which can override these defaults. For example, -ts_max_steps -ts_max_time to override the defaults set by TSSetDuration(). - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSetDuration(ts,time_steps_max,time_total_max);CHKERRQ(ierr); ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve the problem - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Evaluate initial conditions */ ierr = InitialConditions(u,&appctx);CHKERRQ(ierr); /* Run the timestepping solver */ ierr = TSStep(ts,&steps,&ftime);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - View timestepping solver info - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %G, avg. error (max norm) = %G\n", appctx.norm_2/steps,appctx.norm_max/steps);CHKERRQ(ierr); ierr = TSView(ts,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Free work space. All PETSc objects should be destroyed when they are no longer needed. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSDestroy(ts);CHKERRQ(ierr); ierr = MatDestroy(A);CHKERRQ(ierr); ierr = VecDestroy(u);CHKERRQ(ierr); ierr = PetscViewerDestroy(appctx.viewer1);CHKERRQ(ierr); ierr = PetscViewerDestroy(appctx.viewer2);CHKERRQ(ierr); ierr = VecDestroy(appctx.solution);CHKERRQ(ierr); /* Always call PetscFinalize() before exiting a program. This routine - finalizes the PETSc libraries as well as MPI - provides summary and diagnostic information if certain runtime options are chosen (e.g., -log_summary). */ ierr = PetscFinalize();CHKERRQ(ierr); return 0; } /* --------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "InitialConditions" /* InitialConditions - Computes the solution at the initial time. Input Parameter: u - uninitialized solution vector (global) appctx - user-defined application context Output Parameter: u - vector with solution at initial time (global) */ PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) { PetscScalar *u_localptr,h = appctx->h; PetscInt i; PetscErrorCode ierr; /* Get a pointer to vector data. - For default PETSc vectors, VecGetArray() returns a pointer to the data array. Otherwise, the routine is implementation dependent. - You MUST call VecRestoreArray() when you no longer need access to the array. - Note that the Fortran interface to VecGetArray() differs from the C version. See the users manual for details. */ ierr = VecGetArray(u,&u_localptr);CHKERRQ(ierr); /* We initialize the solution array by simply writing the solution directly into the array locations. Alternatively, we could use VecSetValues() or VecSetValuesLocal(). */ for (i=0; im; i++) { u_localptr[i] = PetscCosScalar(PETSC_PI*i*6.*h) + 3.*PetscCosScalar(PETSC_PI*i*2.*h); } /* Restore vector */ ierr = VecRestoreArray(u,&u_localptr);CHKERRQ(ierr); /* Print debugging information if desired */ if (appctx->debug) { printf("initial guess vector\n"); ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); } return 0; } /* --------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "ExactSolution" /* ExactSolution - Computes the exact solution at a given time. Input Parameters: t - current time solution - vector in which exact solution will be computed appctx - user-defined application context Output Parameter: solution - vector with the newly computed exact solution */ PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) { PetscScalar *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t; PetscInt i; PetscErrorCode ierr; /* Get a pointer to vector data. */ ierr = VecGetArray(solution,&s_localptr);CHKERRQ(ierr); /* Simply write the solution directly into the array locations. Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). */ ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc); ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc); sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; for (i=0; im; i++) { s_localptr[i] = PetscCosScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscCosScalar(sc2*(PetscReal)i)*ex2; } /* Restore vector */ ierr = VecRestoreArray(solution,&s_localptr);CHKERRQ(ierr); return 0; } /* --------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "Monitor" /* Monitor - User-provided routine to monitor the solution computed at each timestep. This example plots the solution and computes the error in two different norms. Input Parameters: ts - the timestep context step - the count of the current step (with 0 meaning the initial condition) time - the current time u - the solution at this timestep ctx - the user-provided context for this monitoring routine. In this case we use the application context which contains information about the problem size, workspace and the exact solution. */ PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) { AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ PetscErrorCode ierr; PetscReal norm_2,norm_max; /* View a graph of the current iterate */ ierr = VecView(u,appctx->viewer2);CHKERRQ(ierr); /* Compute the exact solution */ ierr = ExactSolution(time,appctx->solution,appctx);CHKERRQ(ierr); /* Print debugging information if desired */ if (appctx->debug) { printf("Computed solution vector\n"); ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); printf("Exact solution vector\n"); ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); } /* Compute the 2-norm and max-norm of the error */ ierr = VecAXPY(appctx->solution,-1.0,u);CHKERRQ(ierr); ierr = VecNorm(appctx->solution,NORM_2,&norm_2);CHKERRQ(ierr); norm_2 = sqrt(appctx->h)*norm_2; ierr = VecNorm(appctx->solution,NORM_MAX,&norm_max);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Timestep %D: time = %G, 2-norm error = %G, max norm error = %G\n", step,time,norm_2,norm_max);CHKERRQ(ierr); appctx->norm_2 += norm_2; appctx->norm_max += norm_max; /* View a graph of the error */ ierr = VecView(appctx->solution,appctx->viewer1);CHKERRQ(ierr); /* Print debugging information if desired */ if (appctx->debug) { printf("Error vector\n"); ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); } return 0; } /* --------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "RHSMatrixHeat" /* RHSMatrixHeat - User-provided routine to compute the right-hand-side matrix for the heat equation. Input Parameters: ts - the TS context t - current time global_in - global input vector dummy - optional user-defined context, as set by TSetRHSJacobian() Output Parameters: AA - Jacobian matrix BB - optionally different preconditioning matrix str - flag indicating matrix structure Notes: Recall that MatSetValues() uses 0-based row and column numbers in Fortran as well as in C. */ PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx) { Mat A = *AA; /* Jacobian matrix */ AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ PetscInt mstart = 0; PetscInt mend = appctx->m; PetscErrorCode ierr; PetscInt i,idx[3]; PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Compute entries for the locally owned part of the matrix - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Set matrix rows corresponding to boundary data */ mstart = 0; v[0] = 1.0; ierr = MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr); mstart++; mend--; v[0] = 1.0; ierr = MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr); /* Set matrix rows corresponding to interior data. We construct the matrix one row at a time. */ v[0] = sone; v[1] = stwo; v[2] = sone; for (i=mstart; i