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