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