1c4762a1bSJed Brown 2*c2eed0edSBarry Smith 3c4762a1bSJed Brown static char help[] ="Model Equations for Advection-Diffusion\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* 6c4762a1bSJed Brown Page 9, Section 1.2 Model Equations for Advection-Diffusion 7c4762a1bSJed Brown 8c4762a1bSJed Brown u_t = a u_x + d u_xx 9c4762a1bSJed Brown 10c4762a1bSJed Brown The initial conditions used here different then in the book. 11c4762a1bSJed Brown 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* 15c4762a1bSJed Brown Helpful runtime linear solver options: 16c4762a1bSJed Brown -pc_type mg -da_refine 2 -snes_monitor -ksp_monitor -ts_view (geometric multigrid with three levels) 17c4762a1bSJed Brown 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* 21c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 22c4762a1bSJed Brown automatically includes: 23c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 24c4762a1bSJed Brown petscmat.h - matrices 25c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 26c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 27c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 28c4762a1bSJed Brown */ 29c4762a1bSJed Brown 30c4762a1bSJed Brown #include <petscts.h> 31c4762a1bSJed Brown #include <petscdm.h> 32c4762a1bSJed Brown #include <petscdmda.h> 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown User-defined application context - contains data needed by the 36c4762a1bSJed Brown application-provided call-back routines. 37c4762a1bSJed Brown */ 38c4762a1bSJed Brown typedef struct { 39c4762a1bSJed Brown PetscScalar a,d; /* advection and diffusion strength */ 40c4762a1bSJed Brown PetscBool upwind; 41c4762a1bSJed Brown } AppCtx; 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* 44c4762a1bSJed Brown User-defined routines 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown extern PetscErrorCode InitialConditions(TS,Vec,AppCtx*); 47c4762a1bSJed Brown extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*); 48c4762a1bSJed Brown extern PetscErrorCode Solution(TS,PetscReal,Vec,AppCtx*); 49c4762a1bSJed Brown 50c4762a1bSJed Brown int main(int argc,char **argv) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 53c4762a1bSJed Brown TS ts; /* timestepping context */ 54c4762a1bSJed Brown Vec U; /* approximate solution vector */ 55c4762a1bSJed Brown PetscErrorCode ierr; 56c4762a1bSJed Brown PetscReal dt; 57c4762a1bSJed Brown DM da; 58c4762a1bSJed Brown PetscInt M; 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown Initialize program and set problem parameters 62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 65c4762a1bSJed Brown appctx.a = 1.0; 66c4762a1bSJed Brown appctx.d = 0.0; 67c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-a",&appctx.a,NULL);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-d",&appctx.d,NULL);CHKERRQ(ierr); 69c4762a1bSJed Brown appctx.upwind = PETSC_TRUE; 70c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-upwind",&appctx.upwind,NULL);CHKERRQ(ierr); 71c4762a1bSJed Brown 72c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, 60, 1, 1,NULL,&da);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 75c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 76c4762a1bSJed Brown Create vector data structures 77c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* 80c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 81c4762a1bSJed Brown */ 82c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&U);CHKERRQ(ierr); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 85c4762a1bSJed Brown Create timestepping solver context 86c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 87c4762a1bSJed Brown 88c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* 92c4762a1bSJed Brown For linear problems with a time-dependent f(U,t) in the equation 93c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 94c4762a1bSJed Brown as a time-dependent matrix. 95c4762a1bSJed Brown */ 96c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSMatrixHeat,&appctx);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = TSSetSolutionFunction(ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void*))Solution,&appctx);CHKERRQ(ierr); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Customize timestepping solver: 102c4762a1bSJed Brown - Set timestepping duration info 103c4762a1bSJed Brown Then set runtime options, which can override these defaults. 104c4762a1bSJed Brown For example, 105c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 106c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 107c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108c4762a1bSJed Brown 109c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 110c4762a1bSJed Brown dt = .48/(M*M); 111c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,1000);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = TSSetMaxTime(ts,100.0);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* 119c4762a1bSJed Brown Evaluate initial conditions 120c4762a1bSJed Brown */ 121c4762a1bSJed Brown ierr = InitialConditions(ts,U,&appctx);CHKERRQ(ierr); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* 124c4762a1bSJed Brown Run the timestepping solver 125c4762a1bSJed Brown */ 126c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 127c4762a1bSJed Brown 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 131c4762a1bSJed Brown are no longer needed. 132c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 133c4762a1bSJed Brown 134c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 140c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 141c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 142c4762a1bSJed Brown options are chosen (e.g., -log_view). 143c4762a1bSJed Brown */ 144c4762a1bSJed Brown ierr = PetscFinalize(); 145c4762a1bSJed Brown return ierr; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 148c4762a1bSJed Brown /* 149c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 150c4762a1bSJed Brown 151c4762a1bSJed Brown Input Parameter: 152c4762a1bSJed Brown u - uninitialized solution vector (global) 153c4762a1bSJed Brown appctx - user-defined application context 154c4762a1bSJed Brown 155c4762a1bSJed Brown Output Parameter: 156c4762a1bSJed Brown u - vector with solution at initial time (global) 157c4762a1bSJed Brown */ 158c4762a1bSJed Brown PetscErrorCode InitialConditions(TS ts,Vec U,AppCtx *appctx) 159c4762a1bSJed Brown { 160c4762a1bSJed Brown PetscScalar *u,h; 161c4762a1bSJed Brown PetscErrorCode ierr; 162c4762a1bSJed Brown PetscInt i,mstart,mend,xm,M; 163c4762a1bSJed Brown DM da; 164c4762a1bSJed Brown 165c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 168c4762a1bSJed Brown h = 1.0/M; 169c4762a1bSJed Brown mend = mstart + xm; 170c4762a1bSJed Brown /* 171c4762a1bSJed Brown Get a pointer to vector data. 172c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 173c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 174c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 175c4762a1bSJed Brown the array. 176c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 177c4762a1bSJed Brown C version. See the users manual for details. 178c4762a1bSJed Brown */ 179c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* 182c4762a1bSJed Brown We initialize the solution array by simply writing the solution 183c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 184c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 185c4762a1bSJed Brown */ 186c4762a1bSJed Brown for (i=mstart; i<mend; i++) u[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* 189c4762a1bSJed Brown Restore vector 190c4762a1bSJed Brown */ 191c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 192c4762a1bSJed Brown return 0; 193c4762a1bSJed Brown } 194c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 195c4762a1bSJed Brown /* 196c4762a1bSJed Brown Solution - Computes the exact solution at a given time. 197c4762a1bSJed Brown 198c4762a1bSJed Brown Input Parameters: 199c4762a1bSJed Brown t - current time 200c4762a1bSJed Brown solution - vector in which exact solution will be computed 201c4762a1bSJed Brown appctx - user-defined application context 202c4762a1bSJed Brown 203c4762a1bSJed Brown Output Parameter: 204c4762a1bSJed Brown solution - vector with the newly computed exact solution 205c4762a1bSJed Brown */ 206c4762a1bSJed Brown PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *appctx) 207c4762a1bSJed Brown { 208c4762a1bSJed Brown PetscScalar *u,ex1,ex2,sc1,sc2,h; 209c4762a1bSJed Brown PetscErrorCode ierr; 210c4762a1bSJed Brown PetscInt i,mstart,mend,xm,M; 211c4762a1bSJed Brown DM da; 212c4762a1bSJed Brown 213c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 216c4762a1bSJed Brown h = 1.0/M; 217c4762a1bSJed Brown mend = mstart + xm; 218c4762a1bSJed Brown /* 219c4762a1bSJed Brown Get a pointer to vector data. 220c4762a1bSJed Brown */ 221c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* 224c4762a1bSJed Brown Simply write the solution directly into the array locations. 225c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 226c4762a1bSJed Brown */ 227c4762a1bSJed Brown ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*appctx->d*t); 228c4762a1bSJed Brown ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*appctx->d*t); 229c4762a1bSJed Brown sc1 = PETSC_PI*6.*h; sc2 = PETSC_PI*2.*h; 230c4762a1bSJed Brown for (i=mstart; i<mend; i++) u[i] = PetscSinScalar(sc1*(PetscReal)i + appctx->a*PETSC_PI*6.*t)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i + appctx->a*PETSC_PI*2.*t)*ex2; 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* 233c4762a1bSJed Brown Restore vector 234c4762a1bSJed Brown */ 235c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 236c4762a1bSJed Brown return 0; 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 240c4762a1bSJed Brown /* 241c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 242c4762a1bSJed Brown matrix for the heat equation. 243c4762a1bSJed Brown 244c4762a1bSJed Brown Input Parameters: 245c4762a1bSJed Brown ts - the TS context 246c4762a1bSJed Brown t - current time 247c4762a1bSJed Brown global_in - global input vector 248c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 249c4762a1bSJed Brown 250c4762a1bSJed Brown Output Parameters: 251c4762a1bSJed Brown AA - Jacobian matrix 252c4762a1bSJed Brown BB - optionally different preconditioning matrix 253c4762a1bSJed Brown str - flag indicating matrix structure 254c4762a1bSJed Brown 255c4762a1bSJed Brown Notes: 256c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 257c4762a1bSJed Brown in Fortran as well as in C. 258c4762a1bSJed Brown */ 259c4762a1bSJed Brown PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec U,Mat AA,Mat BB,void *ctx) 260c4762a1bSJed Brown { 261c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 262c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 263c4762a1bSJed Brown PetscInt mstart, mend; 264c4762a1bSJed Brown PetscErrorCode ierr; 265c4762a1bSJed Brown PetscInt i,idx[3],M,xm; 266c4762a1bSJed Brown PetscScalar v[3],h; 267c4762a1bSJed Brown DM da; 268c4762a1bSJed Brown 269c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr); 272c4762a1bSJed Brown h = 1.0/M; 273c4762a1bSJed Brown mend = mstart + xm; 274c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 275c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 276c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 277c4762a1bSJed Brown /* 278c4762a1bSJed Brown Set matrix rows corresponding to boundary data 279c4762a1bSJed Brown */ 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* diffusion */ 282c4762a1bSJed Brown v[0] = appctx->d/(h*h); 283c4762a1bSJed Brown v[1] = -2.0*appctx->d/(h*h); 284c4762a1bSJed Brown v[2] = appctx->d/(h*h); 285c4762a1bSJed Brown if (!mstart) { 286c4762a1bSJed Brown idx[0] = M-1; idx[1] = 0; idx[2] = 1; 287c4762a1bSJed Brown ierr = MatSetValues(A,1,&mstart,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 288c4762a1bSJed Brown mstart++; 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown if (mend == M) { 292c4762a1bSJed Brown mend--; 293c4762a1bSJed Brown idx[0] = M-2; idx[1] = M-1; idx[2] = 0; 294c4762a1bSJed Brown ierr = MatSetValues(A,1,&mend,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 299c4762a1bSJed Brown matrix one row at a time. 300c4762a1bSJed Brown */ 301c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 302c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 303c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 307c4762a1bSJed Brown 308c4762a1bSJed Brown ierr = DMDAGetCorners(da,&mstart,0,0,&xm,0,0);CHKERRQ(ierr); 309c4762a1bSJed Brown mend = mstart + xm; 310c4762a1bSJed Brown if (!appctx->upwind) { 311c4762a1bSJed Brown /* advection -- centered differencing */ 312c4762a1bSJed Brown v[0] = -.5*appctx->a/(h); 313c4762a1bSJed Brown v[1] = .5*appctx->a/(h); 314c4762a1bSJed Brown if (!mstart) { 315c4762a1bSJed Brown idx[0] = M-1; idx[1] = 1; 316c4762a1bSJed Brown ierr = MatSetValues(A,1,&mstart,2,idx,v,ADD_VALUES);CHKERRQ(ierr); 317c4762a1bSJed Brown mstart++; 318c4762a1bSJed Brown } 319c4762a1bSJed Brown 320c4762a1bSJed Brown if (mend == M) { 321c4762a1bSJed Brown mend--; 322c4762a1bSJed Brown idx[0] = M-2; idx[1] = 0; 323c4762a1bSJed Brown ierr = MatSetValues(A,1,&mend,2,idx,v,ADD_VALUES);CHKERRQ(ierr); 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 326c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 327c4762a1bSJed Brown idx[0] = i-1; idx[1] = i+1; 328c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,idx,v,ADD_VALUES);CHKERRQ(ierr); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown } else { 331c4762a1bSJed Brown /* advection -- upwinding */ 332c4762a1bSJed Brown v[0] = -appctx->a/(h); 333c4762a1bSJed Brown v[1] = appctx->a/(h); 334c4762a1bSJed Brown if (!mstart) { 335c4762a1bSJed Brown idx[0] = 0; idx[1] = 1; 336c4762a1bSJed Brown ierr = MatSetValues(A,1,&mstart,2,idx,v,ADD_VALUES);CHKERRQ(ierr); 337c4762a1bSJed Brown mstart++; 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown if (mend == M) { 341c4762a1bSJed Brown mend--; 342c4762a1bSJed Brown idx[0] = M-1; idx[1] = 0; 343c4762a1bSJed Brown ierr = MatSetValues(A,1,&mend,2,idx,v,ADD_VALUES);CHKERRQ(ierr); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 346c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 347c4762a1bSJed Brown idx[0] = i; idx[1] = i+1; 348c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,idx,v,ADD_VALUES);CHKERRQ(ierr); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 351c4762a1bSJed Brown 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 354c4762a1bSJed Brown Complete the matrix assembly process and set some options 355c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown Assemble matrix, using the 2-step process: 358c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 359c4762a1bSJed Brown Computations can be done while messages are in transition 360c4762a1bSJed Brown by placing code between these two statements. 361c4762a1bSJed Brown */ 362c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 364c4762a1bSJed Brown 365c4762a1bSJed Brown /* 366c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 367c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 368c4762a1bSJed Brown */ 369c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 370c4762a1bSJed Brown return 0; 371c4762a1bSJed Brown } 372c4762a1bSJed Brown 373c4762a1bSJed Brown 374c4762a1bSJed Brown /*TEST 375c4762a1bSJed Brown 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown args: -pc_type mg -da_refine 2 -ts_view -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3 378c4762a1bSJed Brown requires: double 379*c2eed0edSBarry Smith filter: grep -v "total number of" 380c4762a1bSJed Brown 381c4762a1bSJed Brown test: 382c4762a1bSJed Brown suffix: 2 383c4762a1bSJed Brown args: -pc_type mg -da_refine 2 -ts_view -ts_monitor_draw_solution -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3 384c4762a1bSJed Brown requires: x 385c4762a1bSJed Brown output_file: output/ex3_1.out 386c4762a1bSJed Brown requires: double 387*c2eed0edSBarry Smith filter: grep -v "total number of" 388c4762a1bSJed Brown 389c4762a1bSJed Brown TEST*/ 390