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