1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Solves a time-dependent nonlinear PDE. Uses implicit\n\ 3c4762a1bSJed Brown timestepping. Runtime options include:\n\ 4c4762a1bSJed Brown -M <xg>, where <xg> = number of grid points\n\ 5c4762a1bSJed Brown -debug : Activate debugging printouts\n\ 6c4762a1bSJed Brown -nox : Deactivate x-window graphics\n\n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown /* 9c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 10c4762a1bSJed Brown Processors: n 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* ------------------------------------------------------------------------ 14c4762a1bSJed Brown 15c4762a1bSJed Brown This program solves the PDE 16c4762a1bSJed Brown 17c4762a1bSJed Brown u * u_xx 18c4762a1bSJed Brown u_t = --------- 19c4762a1bSJed Brown 2*(t+1)^2 20c4762a1bSJed Brown 21c4762a1bSJed Brown on the domain 0 <= x <= 1, with boundary conditions 22c4762a1bSJed Brown u(t,0) = t + 1, u(t,1) = 2*t + 2, 23c4762a1bSJed Brown and initial condition 24c4762a1bSJed Brown u(0,x) = 1 + x*x. 25c4762a1bSJed Brown 26c4762a1bSJed Brown The exact solution is: 27c4762a1bSJed Brown u(t,x) = (1 + x*x) * (1 + t) 28c4762a1bSJed Brown 29c4762a1bSJed Brown Note that since the solution is linear in time and quadratic in x, 30c4762a1bSJed Brown the finite difference scheme actually computes the "exact" solution. 31c4762a1bSJed Brown 32c4762a1bSJed Brown We use by default the backward Euler method. 33c4762a1bSJed Brown 34c4762a1bSJed Brown ------------------------------------------------------------------------- */ 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown Include "petscts.h" to use the PETSc timestepping routines. Note that 38c4762a1bSJed Brown this file automatically includes "petscsys.h" and other lower-level 39c4762a1bSJed Brown PETSc include files. 40c4762a1bSJed Brown 41c4762a1bSJed Brown Include the "petscdmda.h" to allow us to use the distributed array data 42c4762a1bSJed Brown structures to manage the parallel grid. 43c4762a1bSJed Brown */ 44c4762a1bSJed Brown #include <petscts.h> 45c4762a1bSJed Brown #include <petscdm.h> 46c4762a1bSJed Brown #include <petscdmda.h> 47c4762a1bSJed Brown #include <petscdraw.h> 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown User-defined application context - contains data needed by the 51c4762a1bSJed Brown application-provided callback routines. 52c4762a1bSJed Brown */ 53c4762a1bSJed Brown typedef struct { 54c4762a1bSJed Brown MPI_Comm comm; /* communicator */ 55c4762a1bSJed Brown DM da; /* distributed array data structure */ 56c4762a1bSJed Brown Vec localwork; /* local ghosted work vector */ 57c4762a1bSJed Brown Vec u_local; /* local ghosted approximate solution vector */ 58c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 59c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 60c4762a1bSJed Brown PetscReal h; /* mesh width: h = 1/(m-1) */ 61c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 62c4762a1bSJed Brown } AppCtx; 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown User-defined routines, provided below. 66c4762a1bSJed Brown */ 67c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 68c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 69c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 70c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 71c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); 72c4762a1bSJed Brown 73c4762a1bSJed Brown int main(int argc,char **argv) 74c4762a1bSJed Brown { 75c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 76c4762a1bSJed Brown TS ts; /* timestepping context */ 77c4762a1bSJed Brown Mat A; /* Jacobian matrix data structure */ 78c4762a1bSJed Brown Vec u; /* approximate solution vector */ 79c4762a1bSJed Brown PetscInt time_steps_max = 100; /* default max timesteps */ 80c4762a1bSJed Brown PetscErrorCode ierr; 81c4762a1bSJed Brown PetscReal dt; 82c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 83c4762a1bSJed Brown PetscBool mymonitor = PETSC_FALSE; 84c4762a1bSJed Brown PetscReal bounds[] = {1.0, 3.3}; 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 87c4762a1bSJed Brown Initialize program and set problem parameters 88c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 89c4762a1bSJed Brown 90c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 91c4762a1bSJed Brown ierr = PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds);CHKERRQ(ierr); 92c4762a1bSJed Brown 93c4762a1bSJed Brown appctx.comm = PETSC_COMM_WORLD; 94c4762a1bSJed Brown appctx.m = 60; 95c4762a1bSJed Brown 96c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&appctx.m,NULL);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);CHKERRQ(ierr); 99c4762a1bSJed Brown 100c4762a1bSJed Brown appctx.h = 1.0/(appctx.m-1.0); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown Create vector data structures 104c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 108c4762a1bSJed Brown and to set up the ghost point communication pattern. There are M 109c4762a1bSJed Brown total grid values spread equally among all the processors. 110c4762a1bSJed Brown */ 111c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = DMSetFromOptions(appctx.da);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = DMSetUp(appctx.da);CHKERRQ(ierr); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* 116c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 117c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 118c4762a1bSJed Brown have the same types. 119c4762a1bSJed Brown */ 120c4762a1bSJed Brown ierr = DMCreateGlobalVector(appctx.da,&u);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = DMCreateLocalVector(appctx.da,&appctx.u_local);CHKERRQ(ierr); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* 124c4762a1bSJed Brown Create local work vector for use in evaluating right-hand-side function; 125c4762a1bSJed Brown create global work vector for storing exact solution. 126c4762a1bSJed Brown */ 127c4762a1bSJed Brown ierr = VecDuplicate(appctx.u_local,&appctx.localwork);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131c4762a1bSJed Brown Create timestepping solver context; set callback routine for 132c4762a1bSJed Brown right-hand-side function evaluation. 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134c4762a1bSJed Brown 135c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140c4762a1bSJed Brown Set optional user-defined monitoring routine 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142c4762a1bSJed Brown 143c4762a1bSJed Brown if (mymonitor) { 144c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&appctx,NULL);CHKERRQ(ierr); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148c4762a1bSJed Brown For nonlinear problems, the user can provide a Jacobian evaluation 149c4762a1bSJed Brown routine (or use a finite differencing approximation). 150c4762a1bSJed Brown 151c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine. 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153c4762a1bSJed Brown 154c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx);CHKERRQ(ierr); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161c4762a1bSJed Brown Set solution vector and initial timestep 162c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 163c4762a1bSJed Brown 164c4762a1bSJed Brown dt = appctx.h/2.0; 165c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168c4762a1bSJed Brown Customize timestepping solver: 169c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 170c4762a1bSJed Brown - Set timestepping duration info 171c4762a1bSJed Brown Then set runtime options, which can override these defaults. 172c4762a1bSJed Brown For example, 173c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 174c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 175c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 176c4762a1bSJed Brown 177c4762a1bSJed Brown ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,time_steps_max);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = TSSetMaxTime(ts,time_total_max);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184c4762a1bSJed Brown Solve the problem 185c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* 188c4762a1bSJed Brown Evaluate initial conditions 189c4762a1bSJed Brown */ 190c4762a1bSJed Brown ierr = InitialConditions(u,&appctx);CHKERRQ(ierr); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* 193c4762a1bSJed Brown Run the timestepping solver 194c4762a1bSJed Brown */ 195c4762a1bSJed Brown ierr = TSSolve(ts,u);CHKERRQ(ierr); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 199c4762a1bSJed Brown are no longer needed. 200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 201c4762a1bSJed Brown 202c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = DMDestroy(&appctx.da);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = VecDestroy(&appctx.localwork);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = VecDestroy(&appctx.u_local);CHKERRQ(ierr); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* 211c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 212c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 213c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 214c4762a1bSJed Brown options are chosen (e.g., -log_view). 215c4762a1bSJed Brown */ 216c4762a1bSJed Brown ierr = PetscFinalize(); 217c4762a1bSJed Brown return ierr; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 220c4762a1bSJed Brown /* 221c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 222c4762a1bSJed Brown 223c4762a1bSJed Brown Input Parameters: 224c4762a1bSJed Brown u - uninitialized solution vector (global) 225c4762a1bSJed Brown appctx - user-defined application context 226c4762a1bSJed Brown 227c4762a1bSJed Brown Output Parameter: 228c4762a1bSJed Brown u - vector with solution at initial time (global) 229c4762a1bSJed Brown */ 230c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 231c4762a1bSJed Brown { 232c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h,x; 233c4762a1bSJed Brown PetscInt i,mybase,myend; 234c4762a1bSJed Brown PetscErrorCode ierr; 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* 237c4762a1bSJed Brown Determine starting point of each processor's range of 238c4762a1bSJed Brown grid values. 239c4762a1bSJed Brown */ 240c4762a1bSJed Brown ierr = VecGetOwnershipRange(u,&mybase,&myend);CHKERRQ(ierr); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* 243c4762a1bSJed Brown Get a pointer to vector data. 244c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 245c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 246c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 247c4762a1bSJed Brown the array. 248c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 249c4762a1bSJed Brown C version. See the users manual for details. 250c4762a1bSJed Brown */ 251c4762a1bSJed Brown ierr = VecGetArray(u,&u_localptr);CHKERRQ(ierr); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown We initialize the solution array by simply writing the solution 255c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 256c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 257c4762a1bSJed Brown */ 258c4762a1bSJed Brown for (i=mybase; i<myend; i++) { 259c4762a1bSJed Brown x = h*(PetscReal)i; /* current location in global grid */ 260c4762a1bSJed Brown u_localptr[i-mybase] = 1.0 + x*x; 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* 264c4762a1bSJed Brown Restore vector 265c4762a1bSJed Brown */ 266c4762a1bSJed Brown ierr = VecRestoreArray(u,&u_localptr);CHKERRQ(ierr); 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* 269c4762a1bSJed Brown Print debugging information if desired 270c4762a1bSJed Brown */ 271c4762a1bSJed Brown if (appctx->debug) { 272c4762a1bSJed Brown ierr = PetscPrintf(appctx->comm,"initial guess vector\n");CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276c4762a1bSJed Brown return 0; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 281c4762a1bSJed Brown 282c4762a1bSJed Brown Input Parameters: 283c4762a1bSJed Brown t - current time 284c4762a1bSJed Brown solution - vector in which exact solution will be computed 285c4762a1bSJed Brown appctx - user-defined application context 286c4762a1bSJed Brown 287c4762a1bSJed Brown Output Parameter: 288c4762a1bSJed Brown solution - vector with the newly computed exact solution 289c4762a1bSJed Brown */ 290c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 291c4762a1bSJed Brown { 292c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,x; 293c4762a1bSJed Brown PetscInt i,mybase,myend; 294c4762a1bSJed Brown PetscErrorCode ierr; 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* 297c4762a1bSJed Brown Determine starting and ending points of each processor's 298c4762a1bSJed Brown range of grid values 299c4762a1bSJed Brown */ 300c4762a1bSJed Brown ierr = VecGetOwnershipRange(solution,&mybase,&myend);CHKERRQ(ierr); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* 303c4762a1bSJed Brown Get a pointer to vector data. 304c4762a1bSJed Brown */ 305c4762a1bSJed Brown ierr = VecGetArray(solution,&s_localptr);CHKERRQ(ierr); 306c4762a1bSJed Brown 307c4762a1bSJed Brown /* 308c4762a1bSJed Brown Simply write the solution directly into the array locations. 309c4762a1bSJed Brown Alternatively, we could use VecSetValues() or VecSetValuesLocal(). 310c4762a1bSJed Brown */ 311c4762a1bSJed Brown for (i=mybase; i<myend; i++) { 312c4762a1bSJed Brown x = h*(PetscReal)i; 313c4762a1bSJed Brown s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x); 314c4762a1bSJed Brown } 315c4762a1bSJed Brown 316c4762a1bSJed Brown /* 317c4762a1bSJed Brown Restore vector 318c4762a1bSJed Brown */ 319c4762a1bSJed Brown ierr = VecRestoreArray(solution,&s_localptr);CHKERRQ(ierr); 320c4762a1bSJed Brown return 0; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 323c4762a1bSJed Brown /* 324c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 325c4762a1bSJed Brown each timestep. This example plots the solution and computes the 326c4762a1bSJed Brown error in two different norms. 327c4762a1bSJed Brown 328c4762a1bSJed Brown Input Parameters: 329c4762a1bSJed Brown ts - the timestep context 330c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 331c4762a1bSJed Brown initial condition) 332c4762a1bSJed Brown time - the current time 333c4762a1bSJed Brown u - the solution at this timestep 334c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 335c4762a1bSJed Brown In this case we use the application context which contains 336c4762a1bSJed Brown information about the problem size, workspace and the exact 337c4762a1bSJed Brown solution. 338c4762a1bSJed Brown */ 339c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 340c4762a1bSJed Brown { 341c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 342c4762a1bSJed Brown PetscErrorCode ierr; 343c4762a1bSJed Brown PetscReal en2,en2s,enmax; 344c4762a1bSJed Brown PetscDraw draw; 345c4762a1bSJed Brown 346c4762a1bSJed Brown /* 347*e1dfdf8eSBarry Smith We use the default X Windows viewer 348c4762a1bSJed Brown PETSC_VIEWER_DRAW_(appctx->comm) 349c4762a1bSJed Brown that is associated with the current communicator. This saves 350c4762a1bSJed Brown the effort of calling PetscViewerDrawOpen() to create the window. 351c4762a1bSJed Brown Note that if we wished to plot several items in separate windows we 352c4762a1bSJed Brown would create each viewer with PetscViewerDrawOpen() and store them in 353c4762a1bSJed Brown the application context, appctx. 354c4762a1bSJed Brown 355c4762a1bSJed Brown PetscReal buffering makes graphics look better. 356c4762a1bSJed Brown */ 357c4762a1bSJed Brown ierr = PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_DRAW_(appctx->comm));CHKERRQ(ierr); 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* 362c4762a1bSJed Brown Compute the exact solution at this timestep 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 ierr = PetscPrintf(appctx->comm,"Computed solution vector\n");CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = VecView(u,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = PetscPrintf(appctx->comm,"Exact solution vector\n");CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);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,&en2);CHKERRQ(ierr); 381c4762a1bSJed Brown en2s = PetscSqrtReal(appctx->h)*en2; /* scale the 2-norm by the grid spacing */ 382c4762a1bSJed Brown ierr = VecNorm(appctx->solution,NORM_MAX,&enmax);CHKERRQ(ierr); 383c4762a1bSJed Brown 384c4762a1bSJed Brown /* 385c4762a1bSJed Brown PetscPrintf() causes only the first processor in this 386c4762a1bSJed Brown communicator to print the timestep information. 387c4762a1bSJed Brown */ 388c4762a1bSJed Brown ierr = PetscPrintf(appctx->comm,"Timestep %D: time = %g 2-norm error = %g max norm error = %g\n",step,(double)time,(double)en2s,(double)enmax);CHKERRQ(ierr); 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* 391c4762a1bSJed Brown Print debugging information if desired 392c4762a1bSJed Brown */ 393c4762a1bSJed Brown if (appctx->debug) { 394c4762a1bSJed Brown ierr = PetscPrintf(appctx->comm,"Error vector\n");CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 396c4762a1bSJed Brown } 397c4762a1bSJed Brown return 0; 398c4762a1bSJed Brown } 399c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 400c4762a1bSJed Brown /* 401c4762a1bSJed Brown RHSFunction - User-provided routine that evalues the right-hand-side 402c4762a1bSJed Brown function of the ODE. This routine is set in the main program by 403c4762a1bSJed Brown calling TSSetRHSFunction(). We compute: 404c4762a1bSJed Brown global_out = F(global_in) 405c4762a1bSJed Brown 406c4762a1bSJed Brown Input Parameters: 407c4762a1bSJed Brown ts - timesteping context 408c4762a1bSJed Brown t - current time 409c4762a1bSJed Brown global_in - vector containing the current iterate 410c4762a1bSJed Brown ctx - (optional) user-provided context for function evaluation. 411c4762a1bSJed Brown In this case we use the appctx defined above. 412c4762a1bSJed Brown 413c4762a1bSJed Brown Output Parameter: 414c4762a1bSJed Brown global_out - vector containing the newly evaluated function 415c4762a1bSJed Brown */ 416c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx) 417c4762a1bSJed Brown { 418c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 419c4762a1bSJed Brown DM da = appctx->da; /* distributed array */ 420c4762a1bSJed Brown Vec local_in = appctx->u_local; /* local ghosted input vector */ 421c4762a1bSJed Brown Vec localwork = appctx->localwork; /* local ghosted work vector */ 422c4762a1bSJed Brown PetscErrorCode ierr; 423c4762a1bSJed Brown PetscInt i,localsize; 424c4762a1bSJed Brown PetscMPIInt rank,size; 425c4762a1bSJed Brown PetscScalar *copyptr,sc; 426c4762a1bSJed Brown const PetscScalar *localptr; 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 429c4762a1bSJed Brown Get ready for local function computations 430c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 431c4762a1bSJed Brown /* 432c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 433c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 434c4762a1bSJed Brown By placing code between these two statements, computations can be 435c4762a1bSJed Brown done while messages are in transition. 436c4762a1bSJed Brown */ 437c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr); 438c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr); 439c4762a1bSJed Brown 440c4762a1bSJed Brown /* 441c4762a1bSJed Brown Access directly the values in our local INPUT work array 442c4762a1bSJed Brown */ 443c4762a1bSJed Brown ierr = VecGetArrayRead(local_in,&localptr);CHKERRQ(ierr); 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown Access directly the values in our local OUTPUT work array 447c4762a1bSJed Brown */ 448c4762a1bSJed Brown ierr = VecGetArray(localwork,©ptr);CHKERRQ(ierr); 449c4762a1bSJed Brown 450c4762a1bSJed Brown sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t)); 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* 453c4762a1bSJed Brown Evaluate our function on the nodes owned by this processor 454c4762a1bSJed Brown */ 455c4762a1bSJed Brown ierr = VecGetLocalSize(local_in,&localsize);CHKERRQ(ierr); 456c4762a1bSJed Brown 457c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 458c4762a1bSJed Brown Compute entries for the locally owned part 459c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 460c4762a1bSJed Brown 461c4762a1bSJed Brown /* 462c4762a1bSJed Brown Handle boundary conditions: This is done by using the boundary condition 463c4762a1bSJed Brown u(t,boundary) = g(t,boundary) 464c4762a1bSJed Brown for some function g. Now take the derivative with respect to t to obtain 465c4762a1bSJed Brown u_{t}(t,boundary) = g_{t}(t,boundary) 466c4762a1bSJed Brown 467c4762a1bSJed Brown In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1 468c4762a1bSJed Brown and u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2 469c4762a1bSJed Brown */ 470ffc4695bSBarry Smith ierr = MPI_Comm_rank(appctx->comm,&rank);CHKERRMPI(ierr); 471ffc4695bSBarry Smith ierr = MPI_Comm_size(appctx->comm,&size);CHKERRMPI(ierr); 472dd400576SPatrick Sanan if (rank == 0) copyptr[0] = 1.0; 473c4762a1bSJed Brown if (rank == size-1) copyptr[localsize-1] = 2.0; 474c4762a1bSJed Brown 475c4762a1bSJed Brown /* 476c4762a1bSJed Brown Handle the interior nodes where the PDE is replace by finite 477c4762a1bSJed Brown difference operators. 478c4762a1bSJed Brown */ 479c4762a1bSJed Brown for (i=1; i<localsize-1; i++) copyptr[i] = localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]); 480c4762a1bSJed Brown 481c4762a1bSJed Brown /* 482c4762a1bSJed Brown Restore vectors 483c4762a1bSJed Brown */ 484c4762a1bSJed Brown ierr = VecRestoreArrayRead(local_in,&localptr);CHKERRQ(ierr); 485c4762a1bSJed Brown ierr = VecRestoreArray(localwork,©ptr);CHKERRQ(ierr); 486c4762a1bSJed Brown 487c4762a1bSJed Brown /* 488c4762a1bSJed Brown Insert values from the local OUTPUT vector into the global 489c4762a1bSJed Brown output vector 490c4762a1bSJed Brown */ 491c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out);CHKERRQ(ierr); 492c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out);CHKERRQ(ierr); 493c4762a1bSJed Brown 494c4762a1bSJed Brown /* Print debugging information if desired */ 495c4762a1bSJed Brown if (appctx->debug) { 496c4762a1bSJed Brown ierr = PetscPrintf(appctx->comm,"RHS function vector\n");CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = VecView(global_out,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 498c4762a1bSJed Brown } 499c4762a1bSJed Brown 500c4762a1bSJed Brown return 0; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 503c4762a1bSJed Brown /* 504c4762a1bSJed Brown RHSJacobian - User-provided routine to compute the Jacobian of 505c4762a1bSJed Brown the nonlinear right-hand-side function of the ODE. 506c4762a1bSJed Brown 507c4762a1bSJed Brown Input Parameters: 508c4762a1bSJed Brown ts - the TS context 509c4762a1bSJed Brown t - current time 510c4762a1bSJed Brown global_in - global input vector 511c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 512c4762a1bSJed Brown 513c4762a1bSJed Brown Output Parameters: 514c4762a1bSJed Brown AA - Jacobian matrix 515c4762a1bSJed Brown BB - optionally different preconditioning matrix 516c4762a1bSJed Brown str - flag indicating matrix structure 517c4762a1bSJed Brown 518c4762a1bSJed Brown Notes: 519c4762a1bSJed Brown RHSJacobian computes entries for the locally owned part of the Jacobian. 520c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 521c4762a1bSJed Brown contiguous chunks of rows across the processors. 522c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 523c4762a1bSJed Brown locally (but any non-local elements will be sent to the 524c4762a1bSJed Brown appropriate processor during matrix assembly). 525c4762a1bSJed Brown - Always specify global row and columns of matrix entries when 526c4762a1bSJed Brown using MatSetValues(). 527c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 528c4762a1bSJed Brown - Note that MatSetValues() uses 0-based row and column numbers 529c4762a1bSJed Brown in Fortran as well as in C. 530c4762a1bSJed Brown */ 531c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat BB,void *ctx) 532c4762a1bSJed Brown { 533c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 534c4762a1bSJed Brown Vec local_in = appctx->u_local; /* local ghosted input vector */ 535c4762a1bSJed Brown DM da = appctx->da; /* distributed array */ 536c4762a1bSJed Brown PetscScalar v[3],sc; 537c4762a1bSJed Brown const PetscScalar *localptr; 538c4762a1bSJed Brown PetscErrorCode ierr; 539c4762a1bSJed Brown PetscInt i,mstart,mend,mstarts,mends,idx[3],is; 540c4762a1bSJed Brown 541c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 542c4762a1bSJed Brown Get ready for local Jacobian computations 543c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 544c4762a1bSJed Brown /* 545c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 546c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 547c4762a1bSJed Brown By placing code between these two statements, computations can be 548c4762a1bSJed Brown done while messages are in transition. 549c4762a1bSJed Brown */ 550c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr); 551c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr); 552c4762a1bSJed Brown 553c4762a1bSJed Brown /* 554c4762a1bSJed Brown Get pointer to vector data 555c4762a1bSJed Brown */ 556c4762a1bSJed Brown ierr = VecGetArrayRead(local_in,&localptr);CHKERRQ(ierr); 557c4762a1bSJed Brown 558c4762a1bSJed Brown /* 559c4762a1bSJed Brown Get starting and ending locally owned rows of the matrix 560c4762a1bSJed Brown */ 561c4762a1bSJed Brown ierr = MatGetOwnershipRange(BB,&mstarts,&mends);CHKERRQ(ierr); 562c4762a1bSJed Brown mstart = mstarts; mend = mends; 563c4762a1bSJed Brown 564c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 565c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 566c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 567c4762a1bSJed Brown contiguous chunks of rows across the processors. 568c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 569c4762a1bSJed Brown locally (but any non-local elements will be sent to the 570c4762a1bSJed Brown appropriate processor during matrix assembly). 571c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 572c4762a1bSJed Brown - We can set matrix entries either using either 573c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(). 574c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 575c4762a1bSJed Brown 576c4762a1bSJed Brown /* 577c4762a1bSJed Brown Set matrix rows corresponding to boundary data 578c4762a1bSJed Brown */ 579c4762a1bSJed Brown if (mstart == 0) { 580c4762a1bSJed Brown v[0] = 0.0; 581c4762a1bSJed Brown ierr = MatSetValues(BB,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr); 582c4762a1bSJed Brown mstart++; 583c4762a1bSJed Brown } 584c4762a1bSJed Brown if (mend == appctx->m) { 585c4762a1bSJed Brown mend--; 586c4762a1bSJed Brown v[0] = 0.0; 587c4762a1bSJed Brown ierr = MatSetValues(BB,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr); 588c4762a1bSJed Brown } 589c4762a1bSJed Brown 590c4762a1bSJed Brown /* 591c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 592c4762a1bSJed Brown matrix one row at a time. 593c4762a1bSJed Brown */ 594c4762a1bSJed Brown sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t)); 595c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 596c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 597c4762a1bSJed Brown is = i - mstart + 1; 598c4762a1bSJed Brown v[0] = sc*localptr[is]; 599c4762a1bSJed Brown v[1] = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]); 600c4762a1bSJed Brown v[2] = sc*localptr[is]; 601c4762a1bSJed Brown ierr = MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 602c4762a1bSJed Brown } 603c4762a1bSJed Brown 604c4762a1bSJed Brown /* 605c4762a1bSJed Brown Restore vector 606c4762a1bSJed Brown */ 607c4762a1bSJed Brown ierr = VecRestoreArrayRead(local_in,&localptr);CHKERRQ(ierr); 608c4762a1bSJed Brown 609c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 610c4762a1bSJed Brown Complete the matrix assembly process and set some options 611c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 612c4762a1bSJed Brown /* 613c4762a1bSJed Brown Assemble matrix, using the 2-step process: 614c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 615c4762a1bSJed Brown Computations can be done while messages are in transition 616c4762a1bSJed Brown by placing code between these two statements. 617c4762a1bSJed Brown */ 618c4762a1bSJed Brown ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 619c4762a1bSJed Brown ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 620c4762a1bSJed Brown if (BB != AA) { 621c4762a1bSJed Brown ierr = MatAssemblyBegin(AA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 622c4762a1bSJed Brown ierr = MatAssemblyEnd(AA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 623c4762a1bSJed Brown } 624c4762a1bSJed Brown 625c4762a1bSJed Brown /* 626c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 627c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 628c4762a1bSJed Brown */ 629c4762a1bSJed Brown ierr = MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 630c4762a1bSJed Brown 631c4762a1bSJed Brown return 0; 632c4762a1bSJed Brown } 633c4762a1bSJed Brown 634c4762a1bSJed Brown /*TEST 635c4762a1bSJed Brown 636c4762a1bSJed Brown test: 637c4762a1bSJed Brown args: -nox -ts_dt 10 -mymonitor 638c4762a1bSJed Brown nsize: 2 639c4762a1bSJed Brown requires: !single 640c4762a1bSJed Brown 641c4762a1bSJed Brown test: 642c4762a1bSJed Brown suffix: tut_1 643c4762a1bSJed Brown nsize: 1 644c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor 645c4762a1bSJed Brown 646c4762a1bSJed Brown test: 647c4762a1bSJed Brown suffix: tut_2 648c4762a1bSJed Brown nsize: 4 649c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -snes_monitor -ksp_monitor 650c4762a1bSJed Brown 651c4762a1bSJed Brown test: 652c4762a1bSJed Brown suffix: tut_3 653c4762a1bSJed Brown nsize: 4 654c4762a1bSJed Brown args: ./ex2 -ts_max_steps 10 -ts_monitor -M 128 655c4762a1bSJed Brown 656c4762a1bSJed Brown TEST*/ 657c4762a1bSJed Brown 658