xref: /petsc/src/ts/tutorials/ex21.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Solves a time-dependent nonlinear PDE with lower and upper bounds on the interior grid points. Uses implicit\n\
2c4762a1bSJed Brown timestepping.  Runtime options include:\n\
3c4762a1bSJed Brown   -M <xg>, where <xg> = number of grid points\n\
4c4762a1bSJed Brown   -debug : Activate debugging printouts\n\
5c4762a1bSJed Brown   -nox   : Deactivate x-window graphics\n\
6c4762a1bSJed Brown   -ul   : lower bound\n\
7c4762a1bSJed Brown   -uh  : upper bound\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    This is a variation of ex2.c to solve the PDE
12c4762a1bSJed Brown 
13c4762a1bSJed Brown                u * u_xx
14c4762a1bSJed Brown          u_t = ---------
15c4762a1bSJed Brown                2*(t+1)^2
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     with box constraints on the interior grid points
18c4762a1bSJed Brown     ul <= u(t,x) <= uh with x != 0,1
19c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
20c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
21c4762a1bSJed Brown     and initial condition
22c4762a1bSJed Brown          u(0,x) = 1 + x*x.
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     The exact solution is:
25c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     We use by default the backward Euler method.
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   ------------------------------------------------------------------------- */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown    Include "petscts.h" to use the PETSc timestepping routines. Note that
33c4762a1bSJed Brown    this file automatically includes "petscsys.h" and other lower-level
34c4762a1bSJed Brown    PETSc include files.
35c4762a1bSJed Brown 
36c4762a1bSJed Brown    Include the "petscdmda.h" to allow us to use the distributed array data
37c4762a1bSJed Brown    structures to manage the parallel grid.
38c4762a1bSJed Brown */
39c4762a1bSJed Brown #include <petscts.h>
40c4762a1bSJed Brown #include <petscdm.h>
41c4762a1bSJed Brown #include <petscdmda.h>
42c4762a1bSJed Brown #include <petscdraw.h>
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*
45c4762a1bSJed Brown    User-defined application context - contains data needed by the
46c4762a1bSJed Brown    application-provided callback routines.
47c4762a1bSJed Brown */
48c4762a1bSJed Brown typedef struct {
49c4762a1bSJed Brown   MPI_Comm  comm;      /* communicator */
50c4762a1bSJed Brown   DM        da;        /* distributed array data structure */
51c4762a1bSJed Brown   Vec       localwork; /* local ghosted work vector */
52c4762a1bSJed Brown   Vec       u_local;   /* local ghosted approximate solution vector */
53c4762a1bSJed Brown   Vec       solution;  /* global exact solution vector */
54c4762a1bSJed Brown   PetscInt  m;         /* total number of grid points */
55c4762a1bSJed Brown   PetscReal h;         /* mesh width: h = 1/(m-1) */
56c4762a1bSJed Brown   PetscBool debug;     /* flag (1 indicates activation of debugging printouts) */
57c4762a1bSJed Brown } AppCtx;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*
60c4762a1bSJed Brown    User-defined routines, provided below.
61c4762a1bSJed Brown */
62c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
63c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
64c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
65c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
66c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
67c4762a1bSJed Brown extern PetscErrorCode SetBounds(Vec, Vec, PetscScalar, PetscScalar, AppCtx *);
68c4762a1bSJed Brown 
main(int argc,char ** argv)69d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
70d71ae5a4SJacob Faibussowitsch {
71c4762a1bSJed Brown   AppCtx      appctx;                /* user-defined application context */
72c4762a1bSJed Brown   TS          ts;                    /* timestepping context */
73c4762a1bSJed Brown   Mat         A;                     /* Jacobian matrix data structure */
74c4762a1bSJed Brown   Vec         u;                     /* approximate solution vector */
75c4762a1bSJed Brown   Vec         r;                     /* residual vector */
76c4762a1bSJed Brown   PetscInt    time_steps_max = 1000; /* default max timesteps */
77c4762a1bSJed Brown   PetscReal   dt;
78c4762a1bSJed Brown   PetscReal   time_total_max = 100.0; /* default max total time */
79c4762a1bSJed Brown   Vec         xl, xu;                 /* Lower and upper bounds on variables */
80c4762a1bSJed Brown   PetscScalar ul = 0.0, uh = 3.0;
81c4762a1bSJed Brown   PetscBool   mymonitor;
82c4762a1bSJed Brown   PetscReal   bounds[] = {1.0, 3.3};
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown      Initialize program and set problem parameters
86c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87c4762a1bSJed Brown 
88327415f7SBarry Smith   PetscFunctionBeginUser;
89c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
909566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, bounds));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
93c4762a1bSJed Brown   appctx.m    = 60;
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &appctx.m, NULL));
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ul", &ul, NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-uh", &uh, NULL));
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-mymonitor", &mymonitor));
99c4762a1bSJed Brown   appctx.h = 1.0 / (appctx.m - 1.0);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102c4762a1bSJed Brown      Create vector data structures
103c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /*
106c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
107c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
108c4762a1bSJed Brown      total grid values spread equally among all the processors.
109c4762a1bSJed Brown   */
1109566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, appctx.m, 1, 1, NULL, &appctx.da));
1119566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1129566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /*
115c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
116c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
117c4762a1bSJed Brown      have the same types.
118c4762a1bSJed Brown   */
1199566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da, &u));
1209566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(appctx.da, &appctx.u_local));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /*
123c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
124c4762a1bSJed Brown      create global work vector for storing exact solution.
125c4762a1bSJed Brown   */
1269566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx.u_local, &appctx.localwork));
1279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Create residual vector */
1309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
131c4762a1bSJed Brown   /* Create lower and upper bound vectors */
1329566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &xl));
1339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &xu));
1349566063dSJacob Faibussowitsch   PetscCall(SetBounds(xl, xu, ul, uh, &appctx));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
138c4762a1bSJed Brown      right-hand-side function evaluation.
139c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1429566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1439566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, r, RHSFunction, &appctx));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146c4762a1bSJed Brown      Set optional user-defined monitoring routine
147c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148c4762a1bSJed Brown 
14948a46eb9SPierre Jolivet   if (mymonitor) PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
153c4762a1bSJed Brown      routine (or use a finite differencing approximation).
154c4762a1bSJed Brown 
155c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
156c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, appctx.m, appctx.m));
1609566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1619566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1629566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &appctx));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165c4762a1bSJed Brown      Set solution vector and initial timestep
166c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   dt = appctx.h / 2.0;
1699566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172c4762a1bSJed Brown      Customize timestepping solver:
173c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
174c4762a1bSJed Brown        - Set timestepping duration info
175c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
176c4762a1bSJed Brown      For example,
177c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
178c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
179c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180c4762a1bSJed Brown 
1819566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1829566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
1839566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
1849566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
185c4762a1bSJed Brown   /* Set lower and upper bound on the solution vector for each time step */
1869566063dSJacob Faibussowitsch   PetscCall(TSVISetVariableBounds(ts, xl, xu));
1879566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown      Solve the problem
191c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /*
194c4762a1bSJed Brown      Evaluate initial conditions
195c4762a1bSJed Brown   */
1969566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /*
199c4762a1bSJed Brown      Run the timestepping solver
200c4762a1bSJed Brown   */
2019566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
205c4762a1bSJed Brown      are no longer needed.
206c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
2109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
2119566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2149566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
2159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.localwork));
2169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.u_local));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /*
220c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
221c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
222c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
223c4762a1bSJed Brown          options are chosen (e.g., -log_view).
224c4762a1bSJed Brown   */
2259566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
226b122ec5aSJacob Faibussowitsch   return 0;
227c4762a1bSJed Brown }
228c4762a1bSJed Brown /* --------------------------------------------------------------------- */
229c4762a1bSJed Brown /*
230c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
231c4762a1bSJed Brown 
232c4762a1bSJed Brown    Input Parameters:
233c4762a1bSJed Brown    u - uninitialized solution vector (global)
234c4762a1bSJed Brown    appctx - user-defined application context
235c4762a1bSJed Brown 
236c4762a1bSJed Brown    Output Parameter:
237c4762a1bSJed Brown    u - vector with solution at initial time (global)
238c4762a1bSJed Brown */
InitialConditions(Vec u,AppCtx * appctx)239d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
240d71ae5a4SJacob Faibussowitsch {
241c4762a1bSJed Brown   PetscScalar *u_localptr, h = appctx->h, x;
242c4762a1bSJed Brown   PetscInt     i, mybase, myend;
243c4762a1bSJed Brown 
2443ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
245c4762a1bSJed Brown   /*
246c4762a1bSJed Brown      Determine starting point of each processor's range of
247c4762a1bSJed Brown      grid values.
248c4762a1bSJed Brown   */
2499566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(u, &mybase, &myend));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /*
252c4762a1bSJed Brown     Get a pointer to vector data.
253c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
254c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
255c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
256c4762a1bSJed Brown       the array.
257c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
258c4762a1bSJed Brown       C version.  See the users manual for details.
259c4762a1bSJed Brown   */
2609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &u_localptr));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /*
263c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
264c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
265c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
266c4762a1bSJed Brown   */
267c4762a1bSJed Brown   for (i = mybase; i < myend; i++) {
268c4762a1bSJed Brown     x                      = h * (PetscReal)i; /* current location in global grid */
269c4762a1bSJed Brown     u_localptr[i - mybase] = 1.0 + x * x;
270c4762a1bSJed Brown   }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   /*
273c4762a1bSJed Brown      Restore vector
274c4762a1bSJed Brown   */
2759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &u_localptr));
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /*
278c4762a1bSJed Brown      Print debugging information if desired
279c4762a1bSJed Brown   */
280c4762a1bSJed Brown   if (appctx->debug) {
2819566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "initial guess vector\n"));
2829566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
283c4762a1bSJed Brown   }
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown /* --------------------------------------------------------------------- */
288c4762a1bSJed Brown /*
28969d47153SPierre Jolivet   SetBounds - Sets the lower and upper bounds on the interior points
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   Input parameters:
292c4762a1bSJed Brown   xl - vector of lower bounds
293c4762a1bSJed Brown   xu - vector of upper bounds
294c4762a1bSJed Brown   ul - constant lower bound for all variables
295c4762a1bSJed Brown   uh - constant upper bound for all variables
296c4762a1bSJed Brown   appctx - Application context
297c4762a1bSJed Brown  */
SetBounds(Vec xl,Vec xu,PetscScalar ul,PetscScalar uh,AppCtx * appctx)298d71ae5a4SJacob Faibussowitsch PetscErrorCode SetBounds(Vec xl, Vec xu, PetscScalar ul, PetscScalar uh, AppCtx *appctx)
299d71ae5a4SJacob Faibussowitsch {
300c4762a1bSJed Brown   PetscScalar *l, *u;
301c4762a1bSJed Brown   PetscMPIInt  rank, size;
302c4762a1bSJed Brown   PetscInt     localsize;
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   PetscFunctionBeginUser;
3059566063dSJacob Faibussowitsch   PetscCall(VecSet(xl, ul));
3069566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, uh));
3079566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &localsize));
3089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xl, &l));
3099566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xu, &u));
310c4762a1bSJed Brown 
3119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm, &rank));
3129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm, &size));
313dd400576SPatrick Sanan   if (rank == 0) {
314c4762a1bSJed Brown     l[0] = -PETSC_INFINITY;
315c4762a1bSJed Brown     u[0] = PETSC_INFINITY;
316c4762a1bSJed Brown   }
317c4762a1bSJed Brown   if (rank == size - 1) {
318c4762a1bSJed Brown     l[localsize - 1] = -PETSC_INFINITY;
319c4762a1bSJed Brown     u[localsize - 1] = PETSC_INFINITY;
320c4762a1bSJed Brown   }
3219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xl, &l));
3229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xu, &u));
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
324c4762a1bSJed Brown }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown /* --------------------------------------------------------------------- */
327c4762a1bSJed Brown /*
328c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
329c4762a1bSJed Brown 
330c4762a1bSJed Brown    Input Parameters:
331c4762a1bSJed Brown    t - current time
332c4762a1bSJed Brown    solution - vector in which exact solution will be computed
333c4762a1bSJed Brown    appctx - user-defined application context
334c4762a1bSJed Brown 
335c4762a1bSJed Brown    Output Parameter:
336c4762a1bSJed Brown    solution - vector with the newly computed exact solution
337c4762a1bSJed Brown */
ExactSolution(PetscReal t,Vec solution,AppCtx * appctx)338d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx)
339d71ae5a4SJacob Faibussowitsch {
340c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, x;
341c4762a1bSJed Brown   PetscInt     i, mybase, myend;
342c4762a1bSJed Brown 
3433ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
344c4762a1bSJed Brown   /*
345c4762a1bSJed Brown      Determine starting and ending points of each processor's
346c4762a1bSJed Brown      range of grid values
347c4762a1bSJed Brown   */
3489566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(solution, &mybase, &myend));
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   /*
351c4762a1bSJed Brown      Get a pointer to vector data.
352c4762a1bSJed Brown   */
3539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution, &s_localptr));
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   /*
356c4762a1bSJed Brown      Simply write the solution directly into the array locations.
357c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
358c4762a1bSJed Brown   */
359c4762a1bSJed Brown   for (i = mybase; i < myend; i++) {
360c4762a1bSJed Brown     x                      = h * (PetscReal)i;
361c4762a1bSJed Brown     s_localptr[i - mybase] = (t + 1.0) * (1.0 + x * x);
362c4762a1bSJed Brown   }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /*
365c4762a1bSJed Brown      Restore vector
366c4762a1bSJed Brown   */
3679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution, &s_localptr));
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369c4762a1bSJed Brown }
370c4762a1bSJed Brown /* --------------------------------------------------------------------- */
371c4762a1bSJed Brown /*
372c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
373c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
374c4762a1bSJed Brown    error in two different norms.
375c4762a1bSJed Brown 
376c4762a1bSJed Brown    Input Parameters:
377c4762a1bSJed Brown    ts     - the timestep context
378c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
379c4762a1bSJed Brown             initial condition)
380c4762a1bSJed Brown    time   - the current time
381c4762a1bSJed Brown    u      - the solution at this timestep
382c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
383c4762a1bSJed Brown             In this case we use the application context which contains
384c4762a1bSJed Brown             information about the problem size, workspace and the exact
385c4762a1bSJed Brown             solution.
386c4762a1bSJed Brown */
Monitor(TS ts,PetscInt step,PetscReal time,Vec u,PetscCtx ctx)387*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, PetscCtx ctx)
388d71ae5a4SJacob Faibussowitsch {
389c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
390c4762a1bSJed Brown   PetscReal en2, en2s, enmax;
391c4762a1bSJed Brown   PetscDraw draw;
392c4762a1bSJed Brown 
3933ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
394c4762a1bSJed Brown   /*
395c4762a1bSJed Brown      We use the default X windows viewer
396c4762a1bSJed Brown              PETSC_VIEWER_DRAW_(appctx->comm)
397c4762a1bSJed Brown      that is associated with the current communicator. This saves
398c4762a1bSJed Brown      the effort of calling PetscViewerDrawOpen() to create the window.
399c4762a1bSJed Brown      Note that if we wished to plot several items in separate windows we
400c4762a1bSJed Brown      would create each viewer with PetscViewerDrawOpen() and store them in
401c4762a1bSJed Brown      the application context, appctx.
402c4762a1bSJed Brown 
403c4762a1bSJed Brown      PetscReal buffering makes graphics look better.
404c4762a1bSJed Brown   */
4059566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm), 0, &draw));
4069566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
4079566063dSJacob Faibussowitsch   PetscCall(VecView(u, PETSC_VIEWER_DRAW_(appctx->comm)));
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   /*
410c4762a1bSJed Brown      Compute the exact solution at this timestep
411c4762a1bSJed Brown   */
4129566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time, appctx->solution, appctx));
413c4762a1bSJed Brown 
414c4762a1bSJed Brown   /*
415c4762a1bSJed Brown      Print debugging information if desired
416c4762a1bSJed Brown   */
417c4762a1bSJed Brown   if (appctx->debug) {
4189566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Computed solution vector\n"));
4199566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
4209566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Exact solution vector\n"));
4219566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD));
422c4762a1bSJed Brown   }
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /*
425c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
426c4762a1bSJed Brown   */
4279566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
4289566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &en2));
429c4762a1bSJed Brown   en2s = PetscSqrtReal(appctx->h) * en2; /* scale the 2-norm by the grid spacing */
4309566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &enmax));
431c4762a1bSJed Brown 
432c4762a1bSJed Brown   /*
433c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
434c4762a1bSJed Brown      communicator to print the timestep information.
435c4762a1bSJed Brown   */
43663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(appctx->comm, "Timestep %" PetscInt_FMT ": time = %g,2-norm error = %g, max norm error = %g\n", step, (double)time, (double)en2s, (double)enmax));
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   /*
439c4762a1bSJed Brown      Print debugging information if desired
440c4762a1bSJed Brown    */
441c4762a1bSJed Brown   /*  if (appctx->debug) {
4429566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"Error vector\n"));
4439566063dSJacob Faibussowitsch      PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD));
444c4762a1bSJed Brown    } */
4453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
446c4762a1bSJed Brown }
447c4762a1bSJed Brown /* --------------------------------------------------------------------- */
448c4762a1bSJed Brown /*
449c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
450c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
451c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
452c4762a1bSJed Brown           global_out = F(global_in)
453c4762a1bSJed Brown 
454c4762a1bSJed Brown    Input Parameters:
455c4762a1bSJed Brown    ts         - timesteping context
456c4762a1bSJed Brown    t          - current time
457c4762a1bSJed Brown    global_in  - vector containing the current iterate
458c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
459c4762a1bSJed Brown                 In this case we use the appctx defined above.
460c4762a1bSJed Brown 
461c4762a1bSJed Brown    Output Parameter:
462c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
463c4762a1bSJed Brown */
RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,PetscCtx ctx)464*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec global_in, Vec global_out, PetscCtx ctx)
465d71ae5a4SJacob Faibussowitsch {
466c4762a1bSJed Brown   AppCtx            *appctx    = (AppCtx *)ctx;     /* user-defined application context */
467c4762a1bSJed Brown   DM                 da        = appctx->da;        /* distributed array */
468c4762a1bSJed Brown   Vec                local_in  = appctx->u_local;   /* local ghosted input vector */
469c4762a1bSJed Brown   Vec                localwork = appctx->localwork; /* local ghosted work vector */
470c4762a1bSJed Brown   PetscInt           i, localsize;
471c4762a1bSJed Brown   PetscMPIInt        rank, size;
472c4762a1bSJed Brown   PetscScalar       *copyptr, sc;
473c4762a1bSJed Brown   const PetscScalar *localptr;
474c4762a1bSJed Brown 
4753ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
476c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
477c4762a1bSJed Brown      Get ready for local function computations
478c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
479c4762a1bSJed Brown   /*
480c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
481c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
482c4762a1bSJed Brown      By placing code between these two statements, computations can be
483c4762a1bSJed Brown      done while messages are in transition.
484c4762a1bSJed Brown   */
4859566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, global_in, INSERT_VALUES, local_in));
4869566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, global_in, INSERT_VALUES, local_in));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   /*
489c4762a1bSJed Brown       Access directly the values in our local INPUT work array
490c4762a1bSJed Brown   */
4919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in, &localptr));
492c4762a1bSJed Brown 
493c4762a1bSJed Brown   /*
494c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
495c4762a1bSJed Brown   */
4969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localwork, &copyptr));
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   sc = 1.0 / (appctx->h * appctx->h * 2.0 * (1.0 + t) * (1.0 + t));
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   /*
501c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
502c4762a1bSJed Brown   */
5039566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local_in, &localsize));
504c4762a1bSJed Brown 
505c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
506c4762a1bSJed Brown      Compute entries for the locally owned part
507c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
508c4762a1bSJed Brown 
509c4762a1bSJed Brown   /*
510c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
511c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
512c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
513c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
514c4762a1bSJed Brown 
515c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
516c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
517c4762a1bSJed Brown   */
5189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm, &rank));
5199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm, &size));
520dd400576SPatrick Sanan   if (rank == 0) copyptr[0] = 1.0;
521c4762a1bSJed Brown   if (rank == size - 1) copyptr[localsize - 1] = (t < .5) ? 2.0 : 0.0;
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   /*
524c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
525c4762a1bSJed Brown      difference operators.
526c4762a1bSJed Brown   */
527c4762a1bSJed Brown   for (i = 1; i < localsize - 1; i++) copyptr[i] = localptr[i] * sc * (localptr[i + 1] + localptr[i - 1] - 2.0 * localptr[i]);
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   /*
530c4762a1bSJed Brown      Restore vectors
531c4762a1bSJed Brown   */
5329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in, &localptr));
5339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localwork, &copyptr));
534c4762a1bSJed Brown 
535c4762a1bSJed Brown   /*
536c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
537c4762a1bSJed Brown      output vector
538c4762a1bSJed Brown   */
5399566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, localwork, INSERT_VALUES, global_out));
5409566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, localwork, INSERT_VALUES, global_out));
541c4762a1bSJed Brown 
542c4762a1bSJed Brown   /* Print debugging information if desired */
543c4762a1bSJed Brown   /*  if (appctx->debug) {
5449566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"RHS function vector\n"));
5459566063dSJacob Faibussowitsch      PetscCall(VecView(global_out,PETSC_VIEWER_STDOUT_WORLD));
546c4762a1bSJed Brown    } */
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
548c4762a1bSJed Brown }
549c4762a1bSJed Brown /* --------------------------------------------------------------------- */
550c4762a1bSJed Brown /*
551c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
552c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
553c4762a1bSJed Brown 
554c4762a1bSJed Brown    Input Parameters:
555c4762a1bSJed Brown    ts - the TS context
556c4762a1bSJed Brown    t - current time
557c4762a1bSJed Brown    global_in - global input vector
558c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
559c4762a1bSJed Brown 
560c4762a1bSJed Brown    Output Parameters:
561c4762a1bSJed Brown    AA - Jacobian matrix
5627addb90fSBarry Smith    BB - optionally different matrix used to construct the preconditioner
563c4762a1bSJed Brown 
564c4762a1bSJed Brown   Notes:
565c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
566c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
567c4762a1bSJed Brown      contiguous chunks of rows across the processors.
568c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
569c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
570c4762a1bSJed Brown      appropriate processor during matrix assembly).
571c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
572c4762a1bSJed Brown      using MatSetValues().
573c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
574c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
575c4762a1bSJed Brown      in Fortran as well as in C.
576c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat B,PetscCtx ctx)577*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec global_in, Mat AA, Mat B, PetscCtx ctx)
578d71ae5a4SJacob Faibussowitsch {
579c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx *)ctx;   /* user-defined application context */
580c4762a1bSJed Brown   Vec                local_in = appctx->u_local; /* local ghosted input vector */
581c4762a1bSJed Brown   DM                 da       = appctx->da;      /* distributed array */
582c4762a1bSJed Brown   PetscScalar        v[3], sc;
583c4762a1bSJed Brown   const PetscScalar *localptr;
584c4762a1bSJed Brown   PetscInt           i, mstart, mend, mstarts, mends, idx[3], is;
585c4762a1bSJed Brown 
5863ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
587c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
588c4762a1bSJed Brown      Get ready for local Jacobian computations
589c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
590c4762a1bSJed Brown   /*
591c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
592c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
593c4762a1bSJed Brown      By placing code between these two statements, computations can be
594c4762a1bSJed Brown      done while messages are in transition.
595c4762a1bSJed Brown   */
5969566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, global_in, INSERT_VALUES, local_in));
5979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, global_in, INSERT_VALUES, local_in));
598c4762a1bSJed Brown 
599c4762a1bSJed Brown   /*
600c4762a1bSJed Brown      Get pointer to vector data
601c4762a1bSJed Brown   */
6029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in, &localptr));
603c4762a1bSJed Brown 
604c4762a1bSJed Brown   /*
605c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
606c4762a1bSJed Brown   */
6079566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &mstarts, &mends));
6089371c9d4SSatish Balay   mstart = mstarts;
6099371c9d4SSatish Balay   mend   = mends;
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
612c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
613c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
614c4762a1bSJed Brown         contiguous chunks of rows across the processors.
615c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
616c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
617c4762a1bSJed Brown         appropriate processor during matrix assembly).
618c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
619c4762a1bSJed Brown       - We can set matrix entries either using either
620c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
621c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   /*
624c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
625c4762a1bSJed Brown   */
626c4762a1bSJed Brown   if (mstart == 0) {
627c4762a1bSJed Brown     v[0] = 0.0;
6289566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
629c4762a1bSJed Brown     mstart++;
630c4762a1bSJed Brown   }
631c4762a1bSJed Brown   if (mend == appctx->m) {
632c4762a1bSJed Brown     mend--;
633c4762a1bSJed Brown     v[0] = 0.0;
6349566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &mend, 1, &mend, v, INSERT_VALUES));
635c4762a1bSJed Brown   }
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   /*
638c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
639c4762a1bSJed Brown      matrix one row at a time.
640c4762a1bSJed Brown   */
641c4762a1bSJed Brown   sc = 1.0 / (appctx->h * appctx->h * 2.0 * (1.0 + t) * (1.0 + t));
642c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
6439371c9d4SSatish Balay     idx[0] = i - 1;
6449371c9d4SSatish Balay     idx[1] = i;
6459371c9d4SSatish Balay     idx[2] = i + 1;
646c4762a1bSJed Brown     is     = i - mstart + 1;
647c4762a1bSJed Brown     v[0]   = sc * localptr[is];
648c4762a1bSJed Brown     v[1]   = sc * (localptr[is + 1] + localptr[is - 1] - 4.0 * localptr[is]);
649c4762a1bSJed Brown     v[2]   = sc * localptr[is];
6509566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, idx, v, INSERT_VALUES));
651c4762a1bSJed Brown   }
652c4762a1bSJed Brown 
653c4762a1bSJed Brown   /*
654c4762a1bSJed Brown      Restore vector
655c4762a1bSJed Brown   */
6569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in, &localptr));
657c4762a1bSJed Brown 
658c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
659c4762a1bSJed Brown      Complete the matrix assembly process and set some options
660c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
661c4762a1bSJed Brown   /*
662c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
663c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
664c4762a1bSJed Brown      Computations can be done while messages are in transition
665c4762a1bSJed Brown      by placing code between these two statements.
666c4762a1bSJed Brown   */
6679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
6689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
669c4762a1bSJed Brown 
670c4762a1bSJed Brown   /*
671c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
672c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
673c4762a1bSJed Brown   */
6749566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
676c4762a1bSJed Brown }
677c4762a1bSJed Brown 
678c4762a1bSJed Brown /*TEST
679c4762a1bSJed Brown 
680c4762a1bSJed Brown     test:
681c4762a1bSJed Brown       args: -snes_type vinewtonrsls -ts_type glee -mymonitor -ts_max_steps 10 -nox
682c4762a1bSJed Brown       requires: !single
683c4762a1bSJed Brown 
684c4762a1bSJed Brown TEST*/
685