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