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 49d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 50d71ae5a4SJacob Faibussowitsch { 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 62327415f7SBarry Smith PetscFunctionBeginUser; 639566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 64c4762a1bSJed Brown appctx.a = 1.0; 65c4762a1bSJed Brown appctx.d = 0.0; 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-a", &appctx.a, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-d", &appctx.d, NULL)); 68c4762a1bSJed Brown appctx.upwind = PETSC_TRUE; 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-upwind", &appctx.upwind, NULL)); 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 60, 1, 1, NULL, &da)); 729566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 739566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 74c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75c4762a1bSJed Brown Create vector data structures 76c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* 79c4762a1bSJed Brown Create vector data structures for approximate and exact solutions 80c4762a1bSJed Brown */ 819566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &U)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Create timestepping solver context 85c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 889566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown For linear problems with a time-dependent f(U,t) in the equation 92c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 93c4762a1bSJed Brown as a time-dependent matrix. 94c4762a1bSJed Brown */ 959566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, &appctx)); 969566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSMatrixHeat, &appctx)); 979566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))Solution, &appctx)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Customize timestepping solver: 101c4762a1bSJed Brown - Set timestepping duration info 102c4762a1bSJed Brown Then set runtime options, which can override these defaults. 103c4762a1bSJed Brown For example, 104c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 105c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 106c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 109c4762a1bSJed Brown dt = .48 / (M * M); 1109566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 1119566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1000)); 1129566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 100.0)); 1139566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1149566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX)); 1159566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* 118c4762a1bSJed Brown Evaluate initial conditions 119c4762a1bSJed Brown */ 1209566063dSJacob Faibussowitsch PetscCall(InitialConditions(ts, U, &appctx)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* 123c4762a1bSJed Brown Run the timestepping solver 124c4762a1bSJed Brown */ 1259566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 129c4762a1bSJed Brown are no longer needed. 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131c4762a1bSJed Brown 1329566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* 137c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 138c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 139c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 140c4762a1bSJed Brown options are chosen (e.g., -log_view). 141c4762a1bSJed Brown */ 1429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 143b122ec5aSJacob Faibussowitsch return 0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 146c4762a1bSJed Brown /* 147c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 148c4762a1bSJed Brown 149c4762a1bSJed Brown Input Parameter: 150c4762a1bSJed Brown u - uninitialized solution vector (global) 151c4762a1bSJed Brown appctx - user-defined application context 152c4762a1bSJed Brown 153c4762a1bSJed Brown Output Parameter: 154c4762a1bSJed Brown u - vector with solution at initial time (global) 155c4762a1bSJed Brown */ 156d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(TS ts, Vec U, AppCtx *appctx) 157d71ae5a4SJacob Faibussowitsch { 158c4762a1bSJed Brown PetscScalar *u, h; 159c4762a1bSJed Brown PetscInt i, mstart, mend, xm, M; 160c4762a1bSJed Brown DM da; 161c4762a1bSJed Brown 162*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1639566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 1659566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 166c4762a1bSJed Brown h = 1.0 / M; 167c4762a1bSJed Brown mend = mstart + xm; 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown Get a pointer to vector data. 170c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 171c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 172c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 173c4762a1bSJed Brown the array. 174c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 175c4762a1bSJed Brown C version. See the users manual for details. 176c4762a1bSJed Brown */ 1779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* 180c4762a1bSJed Brown We initialize the solution array by simply writing the solution 181c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 182c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 183c4762a1bSJed Brown */ 184c4762a1bSJed Brown for (i = mstart; i < mend; i++) u[i] = PetscSinScalar(PETSC_PI * i * 6. * h) + 3. * PetscSinScalar(PETSC_PI * i * 2. * h); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown Restore vector 188c4762a1bSJed Brown */ 1899566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 190*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 193c4762a1bSJed Brown /* 194c4762a1bSJed Brown Solution - Computes the exact solution at a given time. 195c4762a1bSJed Brown 196c4762a1bSJed Brown Input Parameters: 197c4762a1bSJed Brown t - current time 198c4762a1bSJed Brown solution - vector in which exact solution will be computed 199c4762a1bSJed Brown appctx - user-defined application context 200c4762a1bSJed Brown 201c4762a1bSJed Brown Output Parameter: 202c4762a1bSJed Brown solution - vector with the newly computed exact solution 203c4762a1bSJed Brown */ 204d71ae5a4SJacob Faibussowitsch PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *appctx) 205d71ae5a4SJacob Faibussowitsch { 206c4762a1bSJed Brown PetscScalar *u, ex1, ex2, sc1, sc2, h; 207c4762a1bSJed Brown PetscInt i, mstart, mend, xm, M; 208c4762a1bSJed Brown DM da; 209c4762a1bSJed Brown 210*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2119566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 2139566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 214c4762a1bSJed Brown h = 1.0 / M; 215c4762a1bSJed Brown mend = mstart + xm; 216c4762a1bSJed Brown /* 217c4762a1bSJed Brown Get a pointer to vector data. 218c4762a1bSJed Brown */ 2199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* 222c4762a1bSJed Brown Simply write the solution directly into the array locations. 223c4762a1bSJed Brown Alternatively, we culd use VecSetValues() or VecSetValuesLocal(). 224c4762a1bSJed Brown */ 225c4762a1bSJed Brown ex1 = PetscExpScalar(-36. * PETSC_PI * PETSC_PI * appctx->d * t); 226c4762a1bSJed Brown ex2 = PetscExpScalar(-4. * PETSC_PI * PETSC_PI * appctx->d * t); 2279371c9d4SSatish Balay sc1 = PETSC_PI * 6. * h; 2289371c9d4SSatish Balay sc2 = PETSC_PI * 2. * h; 229c4762a1bSJed 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; 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* 232c4762a1bSJed Brown Restore vector 233c4762a1bSJed Brown */ 2349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 235*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 239c4762a1bSJed Brown /* 240c4762a1bSJed Brown RHSMatrixHeat - User-provided routine to compute the right-hand-side 241c4762a1bSJed Brown matrix for the heat equation. 242c4762a1bSJed Brown 243c4762a1bSJed Brown Input Parameters: 244c4762a1bSJed Brown ts - the TS context 245c4762a1bSJed Brown t - current time 246c4762a1bSJed Brown global_in - global input vector 247c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 248c4762a1bSJed Brown 249c4762a1bSJed Brown Output Parameters: 250c4762a1bSJed Brown AA - Jacobian matrix 251c4762a1bSJed Brown BB - optionally different preconditioning matrix 252c4762a1bSJed Brown str - flag indicating matrix structure 253c4762a1bSJed Brown 254c4762a1bSJed Brown Notes: 255c4762a1bSJed Brown Recall that MatSetValues() uses 0-based row and column numbers 256c4762a1bSJed Brown in Fortran as well as in C. 257c4762a1bSJed Brown */ 258d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixHeat(TS ts, PetscReal t, Vec U, Mat AA, Mat BB, void *ctx) 259d71ae5a4SJacob Faibussowitsch { 260c4762a1bSJed Brown Mat A = AA; /* Jacobian matrix */ 261c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 262c4762a1bSJed Brown PetscInt mstart, mend; 263c4762a1bSJed Brown PetscInt i, idx[3], M, xm; 264c4762a1bSJed Brown PetscScalar v[3], h; 265c4762a1bSJed Brown DM da; 266c4762a1bSJed Brown 267*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2689566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2699566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 2709566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 271c4762a1bSJed Brown h = 1.0 / M; 272c4762a1bSJed Brown mend = mstart + xm; 273c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 274c4762a1bSJed Brown Compute entries for the locally owned part of the matrix 275c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 276c4762a1bSJed Brown /* 277c4762a1bSJed Brown Set matrix rows corresponding to boundary data 278c4762a1bSJed Brown */ 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* diffusion */ 281c4762a1bSJed Brown v[0] = appctx->d / (h * h); 282c4762a1bSJed Brown v[1] = -2.0 * appctx->d / (h * h); 283c4762a1bSJed Brown v[2] = appctx->d / (h * h); 284c4762a1bSJed Brown if (!mstart) { 2859371c9d4SSatish Balay idx[0] = M - 1; 2869371c9d4SSatish Balay idx[1] = 0; 2879371c9d4SSatish Balay idx[2] = 1; 2889566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mstart, 3, idx, v, INSERT_VALUES)); 289c4762a1bSJed Brown mstart++; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown if (mend == M) { 293c4762a1bSJed Brown mend--; 2949371c9d4SSatish Balay idx[0] = M - 2; 2959371c9d4SSatish Balay idx[1] = M - 1; 2969371c9d4SSatish Balay idx[2] = 0; 2979566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mend, 3, idx, v, INSERT_VALUES)); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* 301c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 302c4762a1bSJed Brown matrix one row at a time. 303c4762a1bSJed Brown */ 304c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 3059371c9d4SSatish Balay idx[0] = i - 1; 3069371c9d4SSatish Balay idx[1] = i; 3079371c9d4SSatish Balay idx[2] = i + 1; 3089566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES)); 309c4762a1bSJed Brown } 3109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FLUSH_ASSEMBLY)); 3119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FLUSH_ASSEMBLY)); 312c4762a1bSJed Brown 3139566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &mstart, 0, 0, &xm, 0, 0)); 314c4762a1bSJed Brown mend = mstart + xm; 315c4762a1bSJed Brown if (!appctx->upwind) { 316c4762a1bSJed Brown /* advection -- centered differencing */ 317c4762a1bSJed Brown v[0] = -.5 * appctx->a / (h); 318c4762a1bSJed Brown v[1] = .5 * appctx->a / (h); 319c4762a1bSJed Brown if (!mstart) { 3209371c9d4SSatish Balay idx[0] = M - 1; 3219371c9d4SSatish Balay idx[1] = 1; 3229566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES)); 323c4762a1bSJed Brown mstart++; 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 326c4762a1bSJed Brown if (mend == M) { 327c4762a1bSJed Brown mend--; 3289371c9d4SSatish Balay idx[0] = M - 2; 3299371c9d4SSatish Balay idx[1] = 0; 3309566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES)); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 3349371c9d4SSatish Balay idx[0] = i - 1; 3359371c9d4SSatish Balay idx[1] = i + 1; 3369566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES)); 337c4762a1bSJed Brown } 338c4762a1bSJed Brown } else { 339c4762a1bSJed Brown /* advection -- upwinding */ 340c4762a1bSJed Brown v[0] = -appctx->a / (h); 341c4762a1bSJed Brown v[1] = appctx->a / (h); 342c4762a1bSJed Brown if (!mstart) { 3439371c9d4SSatish Balay idx[0] = 0; 3449371c9d4SSatish Balay idx[1] = 1; 3459566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mstart, 2, idx, v, ADD_VALUES)); 346c4762a1bSJed Brown mstart++; 347c4762a1bSJed Brown } 348c4762a1bSJed Brown 349c4762a1bSJed Brown if (mend == M) { 350c4762a1bSJed Brown mend--; 3519371c9d4SSatish Balay idx[0] = M - 1; 3529371c9d4SSatish Balay idx[1] = 0; 3539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &mend, 2, idx, v, ADD_VALUES)); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 356c4762a1bSJed Brown for (i = mstart; i < mend; i++) { 3579371c9d4SSatish Balay idx[0] = i; 3589371c9d4SSatish Balay idx[1] = i + 1; 3599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, idx, v, ADD_VALUES)); 360c4762a1bSJed Brown } 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 364c4762a1bSJed Brown Complete the matrix assembly process and set some options 365c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 366c4762a1bSJed Brown /* 367c4762a1bSJed Brown Assemble matrix, using the 2-step process: 368c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 369c4762a1bSJed Brown Computations can be done while messages are in transition 370c4762a1bSJed Brown by placing code between these two statements. 371c4762a1bSJed Brown */ 3729566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3739566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* 376c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 377c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 378c4762a1bSJed Brown */ 3799566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 380*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown 383c4762a1bSJed Brown /*TEST 384c4762a1bSJed Brown 385c4762a1bSJed Brown test: 386c4762a1bSJed Brown args: -pc_type mg -da_refine 2 -ts_view -ts_monitor -ts_max_time .3 -mg_levels_ksp_max_it 3 387c4762a1bSJed Brown requires: double 388c2eed0edSBarry Smith filter: grep -v "total number of" 389c4762a1bSJed Brown 390c4762a1bSJed Brown test: 391c4762a1bSJed Brown suffix: 2 392c4762a1bSJed 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 393c4762a1bSJed Brown requires: x 394c4762a1bSJed Brown output_file: output/ex3_1.out 395c4762a1bSJed Brown requires: double 396c2eed0edSBarry Smith filter: grep -v "total number of" 397c4762a1bSJed Brown 398c4762a1bSJed Brown TEST*/ 399