xref: /petsc/src/ts/tutorials/ex21.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves a time-dependent nonlinear PDE with lower and upper bounds on the interior grid points. Uses implicit\n\
3c4762a1bSJed Brown timestepping.  Runtime options include:\n\
4c4762a1bSJed Brown   -M <xg>, where <xg> = number of grid points\n\
5c4762a1bSJed Brown   -debug : Activate debugging printouts\n\
6c4762a1bSJed Brown   -nox   : Deactivate x-window graphics\n\
7c4762a1bSJed Brown   -ul   : lower bound\n\
8c4762a1bSJed Brown   -uh  : upper bound\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /* ------------------------------------------------------------------------
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    This is a variation of ex2.c to solve the PDE
13c4762a1bSJed Brown 
14c4762a1bSJed Brown                u * u_xx
15c4762a1bSJed Brown          u_t = ---------
16c4762a1bSJed Brown                2*(t+1)^2
17c4762a1bSJed Brown 
18c4762a1bSJed Brown     with box constraints on the interior grid points
19c4762a1bSJed Brown     ul <= u(t,x) <= uh with x != 0,1
20c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
21c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
22c4762a1bSJed Brown     and initial condition
23c4762a1bSJed Brown          u(0,x) = 1 + x*x.
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     The exact solution is:
26c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
27c4762a1bSJed Brown 
28c4762a1bSJed Brown     We use by default the backward Euler method.
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   ------------------------------------------------------------------------- */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    Include "petscts.h" to use the PETSc timestepping routines. Note that
34c4762a1bSJed Brown    this file automatically includes "petscsys.h" and other lower-level
35c4762a1bSJed Brown    PETSc include files.
36c4762a1bSJed Brown 
37c4762a1bSJed Brown    Include the "petscdmda.h" to allow us to use the distributed array data
38c4762a1bSJed Brown    structures to manage the parallel grid.
39c4762a1bSJed Brown */
40c4762a1bSJed Brown #include <petscts.h>
41c4762a1bSJed Brown #include <petscdm.h>
42c4762a1bSJed Brown #include <petscdmda.h>
43c4762a1bSJed Brown #include <petscdraw.h>
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*
46c4762a1bSJed Brown    User-defined application context - contains data needed by the
47c4762a1bSJed Brown    application-provided callback routines.
48c4762a1bSJed Brown */
49c4762a1bSJed Brown typedef struct {
50c4762a1bSJed Brown   MPI_Comm  comm;      /* communicator */
51c4762a1bSJed Brown   DM        da;        /* distributed array data structure */
52c4762a1bSJed Brown   Vec       localwork; /* local ghosted work vector */
53c4762a1bSJed Brown   Vec       u_local;   /* local ghosted approximate solution vector */
54c4762a1bSJed Brown   Vec       solution;  /* global exact solution vector */
55c4762a1bSJed Brown   PetscInt  m;         /* total number of grid points */
56c4762a1bSJed Brown   PetscReal h;         /* mesh width: h = 1/(m-1) */
57c4762a1bSJed Brown   PetscBool debug;     /* flag (1 indicates activation of debugging printouts) */
58c4762a1bSJed Brown } AppCtx;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown /*
61c4762a1bSJed Brown    User-defined routines, provided below.
62c4762a1bSJed Brown */
63c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
64c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
65c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
66c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
67c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal, Vec, AppCtx *);
68c4762a1bSJed Brown extern PetscErrorCode SetBounds(Vec, Vec, PetscScalar, PetscScalar, AppCtx *);
69c4762a1bSJed Brown 
70d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
71d71ae5a4SJacob Faibussowitsch {
72c4762a1bSJed Brown   AppCtx      appctx;                /* user-defined application context */
73c4762a1bSJed Brown   TS          ts;                    /* timestepping context */
74c4762a1bSJed Brown   Mat         A;                     /* Jacobian matrix data structure */
75c4762a1bSJed Brown   Vec         u;                     /* approximate solution vector */
76c4762a1bSJed Brown   Vec         r;                     /* residual vector */
77c4762a1bSJed Brown   PetscInt    time_steps_max = 1000; /* default max timesteps */
78c4762a1bSJed Brown   PetscReal   dt;
79c4762a1bSJed Brown   PetscReal   time_total_max = 100.0; /* default max total time */
80c4762a1bSJed Brown   Vec         xl, xu;                 /* Lower and upper bounds on variables */
81c4762a1bSJed Brown   PetscScalar ul = 0.0, uh = 3.0;
82c4762a1bSJed Brown   PetscBool   mymonitor;
83c4762a1bSJed Brown   PetscReal   bounds[] = {1.0, 3.3};
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown      Initialize program and set problem parameters
87c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88c4762a1bSJed Brown 
89327415f7SBarry Smith   PetscFunctionBeginUser;
909566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
919566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, bounds));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
94c4762a1bSJed Brown   appctx.m    = 60;
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &appctx.m, NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ul", &ul, NULL));
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-uh", &uh, NULL));
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-debug", &appctx.debug));
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-mymonitor", &mymonitor));
100c4762a1bSJed Brown   appctx.h = 1.0 / (appctx.m - 1.0);
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown      Create vector data structures
104c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
108c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
109c4762a1bSJed Brown      total grid values spread equally among all the processors.
110c4762a1bSJed Brown   */
1119566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, appctx.m, 1, 1, NULL, &appctx.da));
1129566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1139566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*
116c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
117c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
118c4762a1bSJed Brown      have the same types.
119c4762a1bSJed Brown   */
1209566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da, &u));
1219566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(appctx.da, &appctx.u_local));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /*
124c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
125c4762a1bSJed Brown      create global work vector for storing exact solution.
126c4762a1bSJed Brown   */
1279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx.u_local, &appctx.localwork));
1289566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.solution));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Create residual vector */
1319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
132c4762a1bSJed Brown   /* Create lower and upper bound vectors */
1339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &xl));
1349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &xu));
1359566063dSJacob Faibussowitsch   PetscCall(SetBounds(xl, xu, ul, uh, &appctx));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
139c4762a1bSJed Brown      right-hand-side function evaluation.
140c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1439566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1449566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, r, RHSFunction, &appctx));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown      Set optional user-defined monitoring routine
148c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149c4762a1bSJed Brown 
15048a46eb9SPierre Jolivet   if (mymonitor) PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
154c4762a1bSJed Brown      routine (or use a finite differencing approximation).
155c4762a1bSJed Brown 
156c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
157c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158c4762a1bSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1609566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, appctx.m, appctx.m));
1619566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1629566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1639566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &appctx));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166c4762a1bSJed Brown      Set solution vector and initial timestep
167c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   dt = appctx.h / 2.0;
1709566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown      Customize timestepping solver:
174c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
175c4762a1bSJed Brown        - Set timestepping duration info
176c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
177c4762a1bSJed Brown      For example,
178c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
179c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
180c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1839566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps_max));
1849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, time_total_max));
1859566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
186c4762a1bSJed Brown   /* Set lower and upper bound on the solution vector for each time step */
1879566063dSJacob Faibussowitsch   PetscCall(TSVISetVariableBounds(ts, xl, xu));
1889566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191c4762a1bSJed Brown      Solve the problem
192c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /*
195c4762a1bSJed Brown      Evaluate initial conditions
196c4762a1bSJed Brown   */
1979566063dSJacob Faibussowitsch   PetscCall(InitialConditions(u, &appctx));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /*
200c4762a1bSJed Brown      Run the timestepping solver
201c4762a1bSJed Brown   */
2029566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
206c4762a1bSJed Brown      are no longer needed.
207c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208c4762a1bSJed Brown 
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
2129566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2159566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
2169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.localwork));
2179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
2189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.u_local));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /*
221c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
222c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
223c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
224c4762a1bSJed Brown          options are chosen (e.g., -log_view).
225c4762a1bSJed Brown   */
2269566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
227b122ec5aSJacob Faibussowitsch   return 0;
228c4762a1bSJed Brown }
229c4762a1bSJed Brown /* --------------------------------------------------------------------- */
230c4762a1bSJed Brown /*
231c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
232c4762a1bSJed Brown 
233c4762a1bSJed Brown    Input Parameters:
234c4762a1bSJed Brown    u - uninitialized solution vector (global)
235c4762a1bSJed Brown    appctx - user-defined application context
236c4762a1bSJed Brown 
237c4762a1bSJed Brown    Output Parameter:
238c4762a1bSJed Brown    u - vector with solution at initial time (global)
239c4762a1bSJed Brown */
240d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
241d71ae5a4SJacob Faibussowitsch {
242c4762a1bSJed Brown   PetscScalar *u_localptr, h = appctx->h, x;
243c4762a1bSJed Brown   PetscInt     i, mybase, myend;
244c4762a1bSJed Brown 
245*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
246c4762a1bSJed Brown   /*
247c4762a1bSJed Brown      Determine starting point of each processor's range of
248c4762a1bSJed Brown      grid values.
249c4762a1bSJed Brown   */
2509566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(u, &mybase, &myend));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /*
253c4762a1bSJed Brown     Get a pointer to vector data.
254c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
255c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
256c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
257c4762a1bSJed Brown       the array.
258c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
259c4762a1bSJed Brown       C version.  See the users manual for details.
260c4762a1bSJed Brown   */
2619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &u_localptr));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /*
264c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
265c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
266c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
267c4762a1bSJed Brown   */
268c4762a1bSJed Brown   for (i = mybase; i < myend; i++) {
269c4762a1bSJed Brown     x                      = h * (PetscReal)i; /* current location in global grid */
270c4762a1bSJed Brown     u_localptr[i - mybase] = 1.0 + x * x;
271c4762a1bSJed Brown   }
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   /*
274c4762a1bSJed Brown      Restore vector
275c4762a1bSJed Brown   */
2769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &u_localptr));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /*
279c4762a1bSJed Brown      Print debugging information if desired
280c4762a1bSJed Brown   */
281c4762a1bSJed Brown   if (appctx->debug) {
2829566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "initial guess vector\n"));
2839566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
284c4762a1bSJed Brown   }
285c4762a1bSJed Brown 
286*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287c4762a1bSJed Brown }
288c4762a1bSJed Brown 
289c4762a1bSJed Brown /* --------------------------------------------------------------------- */
290c4762a1bSJed Brown /*
291c4762a1bSJed Brown   SetBounds - Sets the lower and uper bounds on the interior points
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   Input parameters:
294c4762a1bSJed Brown   xl - vector of lower bounds
295c4762a1bSJed Brown   xu - vector of upper bounds
296c4762a1bSJed Brown   ul - constant lower bound for all variables
297c4762a1bSJed Brown   uh - constant upper bound for all variables
298c4762a1bSJed Brown   appctx - Application context
299c4762a1bSJed Brown  */
300d71ae5a4SJacob Faibussowitsch PetscErrorCode SetBounds(Vec xl, Vec xu, PetscScalar ul, PetscScalar uh, AppCtx *appctx)
301d71ae5a4SJacob Faibussowitsch {
302c4762a1bSJed Brown   PetscScalar *l, *u;
303c4762a1bSJed Brown   PetscMPIInt  rank, size;
304c4762a1bSJed Brown   PetscInt     localsize;
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   PetscFunctionBeginUser;
3079566063dSJacob Faibussowitsch   PetscCall(VecSet(xl, ul));
3089566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, uh));
3099566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xl, &localsize));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xl, &l));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xu, &u));
312c4762a1bSJed Brown 
3139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm, &rank));
3149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm, &size));
315dd400576SPatrick Sanan   if (rank == 0) {
316c4762a1bSJed Brown     l[0] = -PETSC_INFINITY;
317c4762a1bSJed Brown     u[0] = PETSC_INFINITY;
318c4762a1bSJed Brown   }
319c4762a1bSJed Brown   if (rank == size - 1) {
320c4762a1bSJed Brown     l[localsize - 1] = -PETSC_INFINITY;
321c4762a1bSJed Brown     u[localsize - 1] = PETSC_INFINITY;
322c4762a1bSJed Brown   }
3239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xl, &l));
3249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xu, &u));
325*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326c4762a1bSJed Brown }
327c4762a1bSJed Brown 
328c4762a1bSJed Brown /* --------------------------------------------------------------------- */
329c4762a1bSJed Brown /*
330c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
331c4762a1bSJed Brown 
332c4762a1bSJed Brown    Input Parameters:
333c4762a1bSJed Brown    t - current time
334c4762a1bSJed Brown    solution - vector in which exact solution will be computed
335c4762a1bSJed Brown    appctx - user-defined application context
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    Output Parameter:
338c4762a1bSJed Brown    solution - vector with the newly computed exact solution
339c4762a1bSJed Brown */
340d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(PetscReal t, Vec solution, AppCtx *appctx)
341d71ae5a4SJacob Faibussowitsch {
342c4762a1bSJed Brown   PetscScalar *s_localptr, h = appctx->h, x;
343c4762a1bSJed Brown   PetscInt     i, mybase, myend;
344c4762a1bSJed Brown 
345*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
346c4762a1bSJed Brown   /*
347c4762a1bSJed Brown      Determine starting and ending points of each processor's
348c4762a1bSJed Brown      range of grid values
349c4762a1bSJed Brown   */
3509566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(solution, &mybase, &myend));
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   /*
353c4762a1bSJed Brown      Get a pointer to vector data.
354c4762a1bSJed Brown   */
3559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(solution, &s_localptr));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /*
358c4762a1bSJed Brown      Simply write the solution directly into the array locations.
359c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
360c4762a1bSJed Brown   */
361c4762a1bSJed Brown   for (i = mybase; i < myend; i++) {
362c4762a1bSJed Brown     x                      = h * (PetscReal)i;
363c4762a1bSJed Brown     s_localptr[i - mybase] = (t + 1.0) * (1.0 + x * x);
364c4762a1bSJed Brown   }
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   /*
367c4762a1bSJed Brown      Restore vector
368c4762a1bSJed Brown   */
3699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(solution, &s_localptr));
370*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
371c4762a1bSJed Brown }
372c4762a1bSJed Brown /* --------------------------------------------------------------------- */
373c4762a1bSJed Brown /*
374c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
375c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
376c4762a1bSJed Brown    error in two different norms.
377c4762a1bSJed Brown 
378c4762a1bSJed Brown    Input Parameters:
379c4762a1bSJed Brown    ts     - the timestep context
380c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
381c4762a1bSJed Brown             initial condition)
382c4762a1bSJed Brown    time   - the current time
383c4762a1bSJed Brown    u      - the solution at this timestep
384c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
385c4762a1bSJed Brown             In this case we use the application context which contains
386c4762a1bSJed Brown             information about the problem size, workspace and the exact
387c4762a1bSJed Brown             solution.
388c4762a1bSJed Brown */
389d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx)
390d71ae5a4SJacob Faibussowitsch {
391c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx; /* user-defined application context */
392c4762a1bSJed Brown   PetscReal en2, en2s, enmax;
393c4762a1bSJed Brown   PetscDraw draw;
394c4762a1bSJed Brown 
395*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
396c4762a1bSJed Brown   /*
397c4762a1bSJed Brown      We use the default X windows viewer
398c4762a1bSJed Brown              PETSC_VIEWER_DRAW_(appctx->comm)
399c4762a1bSJed Brown      that is associated with the current communicator. This saves
400c4762a1bSJed Brown      the effort of calling PetscViewerDrawOpen() to create the window.
401c4762a1bSJed Brown      Note that if we wished to plot several items in separate windows we
402c4762a1bSJed Brown      would create each viewer with PetscViewerDrawOpen() and store them in
403c4762a1bSJed Brown      the application context, appctx.
404c4762a1bSJed Brown 
405c4762a1bSJed Brown      PetscReal buffering makes graphics look better.
406c4762a1bSJed Brown   */
4079566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm), 0, &draw));
4089566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
4099566063dSJacob Faibussowitsch   PetscCall(VecView(u, PETSC_VIEWER_DRAW_(appctx->comm)));
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   /*
412c4762a1bSJed Brown      Compute the exact solution at this timestep
413c4762a1bSJed Brown   */
4149566063dSJacob Faibussowitsch   PetscCall(ExactSolution(time, appctx->solution, appctx));
415c4762a1bSJed Brown 
416c4762a1bSJed Brown   /*
417c4762a1bSJed Brown      Print debugging information if desired
418c4762a1bSJed Brown   */
419c4762a1bSJed Brown   if (appctx->debug) {
4209566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Computed solution vector\n"));
4219566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
4229566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(appctx->comm, "Exact solution vector\n"));
4239566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_WORLD));
424c4762a1bSJed Brown   }
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   /*
427c4762a1bSJed Brown      Compute the 2-norm and max-norm of the error
428c4762a1bSJed Brown   */
4299566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
4309566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &en2));
431c4762a1bSJed Brown   en2s = PetscSqrtReal(appctx->h) * en2; /* scale the 2-norm by the grid spacing */
4329566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &enmax));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /*
435c4762a1bSJed Brown      PetscPrintf() causes only the first processor in this
436c4762a1bSJed Brown      communicator to print the timestep information.
437c4762a1bSJed Brown   */
43863a3b9bcSJacob 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));
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /*
441c4762a1bSJed Brown      Print debugging information if desired
442c4762a1bSJed Brown    */
443c4762a1bSJed Brown   /*  if (appctx->debug) {
4449566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"Error vector\n"));
4459566063dSJacob Faibussowitsch      PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD));
446c4762a1bSJed Brown    } */
447*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
448c4762a1bSJed Brown }
449c4762a1bSJed Brown /* --------------------------------------------------------------------- */
450c4762a1bSJed Brown /*
451c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
452c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
453c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
454c4762a1bSJed Brown           global_out = F(global_in)
455c4762a1bSJed Brown 
456c4762a1bSJed Brown    Input Parameters:
457c4762a1bSJed Brown    ts         - timesteping context
458c4762a1bSJed Brown    t          - current time
459c4762a1bSJed Brown    global_in  - vector containing the current iterate
460c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
461c4762a1bSJed Brown                 In this case we use the appctx defined above.
462c4762a1bSJed Brown 
463c4762a1bSJed Brown    Output Parameter:
464c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
465c4762a1bSJed Brown */
466d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec global_in, Vec global_out, void *ctx)
467d71ae5a4SJacob Faibussowitsch {
468c4762a1bSJed Brown   AppCtx            *appctx    = (AppCtx *)ctx;     /* user-defined application context */
469c4762a1bSJed Brown   DM                 da        = appctx->da;        /* distributed array */
470c4762a1bSJed Brown   Vec                local_in  = appctx->u_local;   /* local ghosted input vector */
471c4762a1bSJed Brown   Vec                localwork = appctx->localwork; /* local ghosted work vector */
472c4762a1bSJed Brown   PetscInt           i, localsize;
473c4762a1bSJed Brown   PetscMPIInt        rank, size;
474c4762a1bSJed Brown   PetscScalar       *copyptr, sc;
475c4762a1bSJed Brown   const PetscScalar *localptr;
476c4762a1bSJed Brown 
477*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
478c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
479c4762a1bSJed Brown      Get ready for local function computations
480c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
481c4762a1bSJed Brown   /*
482c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
483c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
484c4762a1bSJed Brown      By placing code between these two statements, computations can be
485c4762a1bSJed Brown      done while messages are in transition.
486c4762a1bSJed Brown   */
4879566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, global_in, INSERT_VALUES, local_in));
4889566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, global_in, INSERT_VALUES, local_in));
489c4762a1bSJed Brown 
490c4762a1bSJed Brown   /*
491c4762a1bSJed Brown       Access directly the values in our local INPUT work array
492c4762a1bSJed Brown   */
4939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in, &localptr));
494c4762a1bSJed Brown 
495c4762a1bSJed Brown   /*
496c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
497c4762a1bSJed Brown   */
4989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localwork, &copyptr));
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   sc = 1.0 / (appctx->h * appctx->h * 2.0 * (1.0 + t) * (1.0 + t));
501c4762a1bSJed Brown 
502c4762a1bSJed Brown   /*
503c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
504c4762a1bSJed Brown   */
5059566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local_in, &localsize));
506c4762a1bSJed Brown 
507c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
508c4762a1bSJed Brown      Compute entries for the locally owned part
509c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
510c4762a1bSJed Brown 
511c4762a1bSJed Brown   /*
512c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
513c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
514c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
515c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
516c4762a1bSJed Brown 
517c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
518c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
519c4762a1bSJed Brown   */
5209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(appctx->comm, &rank));
5219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(appctx->comm, &size));
522dd400576SPatrick Sanan   if (rank == 0) copyptr[0] = 1.0;
523c4762a1bSJed Brown   if (rank == size - 1) copyptr[localsize - 1] = (t < .5) ? 2.0 : 0.0;
524c4762a1bSJed Brown 
525c4762a1bSJed Brown   /*
526c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
527c4762a1bSJed Brown      difference operators.
528c4762a1bSJed Brown   */
529c4762a1bSJed Brown   for (i = 1; i < localsize - 1; i++) copyptr[i] = localptr[i] * sc * (localptr[i + 1] + localptr[i - 1] - 2.0 * localptr[i]);
530c4762a1bSJed Brown 
531c4762a1bSJed Brown   /*
532c4762a1bSJed Brown      Restore vectors
533c4762a1bSJed Brown   */
5349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in, &localptr));
5359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localwork, &copyptr));
536c4762a1bSJed Brown 
537c4762a1bSJed Brown   /*
538c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
539c4762a1bSJed Brown      output vector
540c4762a1bSJed Brown   */
5419566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, localwork, INSERT_VALUES, global_out));
5429566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, localwork, INSERT_VALUES, global_out));
543c4762a1bSJed Brown 
544c4762a1bSJed Brown   /* Print debugging information if desired */
545c4762a1bSJed Brown   /*  if (appctx->debug) {
5469566063dSJacob Faibussowitsch      PetscCall(PetscPrintf(appctx->comm,"RHS function vector\n"));
5479566063dSJacob Faibussowitsch      PetscCall(VecView(global_out,PETSC_VIEWER_STDOUT_WORLD));
548c4762a1bSJed Brown    } */
549c4762a1bSJed Brown 
550*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
551c4762a1bSJed Brown }
552c4762a1bSJed Brown /* --------------------------------------------------------------------- */
553c4762a1bSJed Brown /*
554c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
555c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
556c4762a1bSJed Brown 
557c4762a1bSJed Brown    Input Parameters:
558c4762a1bSJed Brown    ts - the TS context
559c4762a1bSJed Brown    t - current time
560c4762a1bSJed Brown    global_in - global input vector
561c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
562c4762a1bSJed Brown 
563c4762a1bSJed Brown    Output Parameters:
564c4762a1bSJed Brown    AA - Jacobian matrix
565c4762a1bSJed Brown    BB - optionally different preconditioning matrix
566c4762a1bSJed Brown    str - flag indicating matrix structure
567c4762a1bSJed Brown 
568c4762a1bSJed Brown   Notes:
569c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
570c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
571c4762a1bSJed Brown      contiguous chunks of rows across the processors.
572c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
573c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
574c4762a1bSJed Brown      appropriate processor during matrix assembly).
575c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
576c4762a1bSJed Brown      using MatSetValues().
577c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
578c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
579c4762a1bSJed Brown      in Fortran as well as in C.
580c4762a1bSJed Brown */
581d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec global_in, Mat AA, Mat B, void *ctx)
582d71ae5a4SJacob Faibussowitsch {
583c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx *)ctx;   /* user-defined application context */
584c4762a1bSJed Brown   Vec                local_in = appctx->u_local; /* local ghosted input vector */
585c4762a1bSJed Brown   DM                 da       = appctx->da;      /* distributed array */
586c4762a1bSJed Brown   PetscScalar        v[3], sc;
587c4762a1bSJed Brown   const PetscScalar *localptr;
588c4762a1bSJed Brown   PetscInt           i, mstart, mend, mstarts, mends, idx[3], is;
589c4762a1bSJed Brown 
590*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
591c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
592c4762a1bSJed Brown      Get ready for local Jacobian computations
593c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
594c4762a1bSJed Brown   /*
595c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
596c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
597c4762a1bSJed Brown      By placing code between these two statements, computations can be
598c4762a1bSJed Brown      done while messages are in transition.
599c4762a1bSJed Brown   */
6009566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, global_in, INSERT_VALUES, local_in));
6019566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, global_in, INSERT_VALUES, local_in));
602c4762a1bSJed Brown 
603c4762a1bSJed Brown   /*
604c4762a1bSJed Brown      Get pointer to vector data
605c4762a1bSJed Brown   */
6069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(local_in, &localptr));
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   /*
609c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
610c4762a1bSJed Brown   */
6119566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &mstarts, &mends));
6129371c9d4SSatish Balay   mstart = mstarts;
6139371c9d4SSatish Balay   mend   = mends;
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
616c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
617c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
618c4762a1bSJed Brown         contiguous chunks of rows across the processors.
619c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
620c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
621c4762a1bSJed Brown         appropriate processor during matrix assembly).
622c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
623c4762a1bSJed Brown       - We can set matrix entries either using either
624c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
625c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
626c4762a1bSJed Brown 
627c4762a1bSJed Brown   /*
628c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
629c4762a1bSJed Brown   */
630c4762a1bSJed Brown   if (mstart == 0) {
631c4762a1bSJed Brown     v[0] = 0.0;
6329566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &mstart, 1, &mstart, v, INSERT_VALUES));
633c4762a1bSJed Brown     mstart++;
634c4762a1bSJed Brown   }
635c4762a1bSJed Brown   if (mend == appctx->m) {
636c4762a1bSJed Brown     mend--;
637c4762a1bSJed Brown     v[0] = 0.0;
6389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &mend, 1, &mend, v, INSERT_VALUES));
639c4762a1bSJed Brown   }
640c4762a1bSJed Brown 
641c4762a1bSJed Brown   /*
642c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
643c4762a1bSJed Brown      matrix one row at a time.
644c4762a1bSJed Brown   */
645c4762a1bSJed Brown   sc = 1.0 / (appctx->h * appctx->h * 2.0 * (1.0 + t) * (1.0 + t));
646c4762a1bSJed Brown   for (i = mstart; i < mend; i++) {
6479371c9d4SSatish Balay     idx[0] = i - 1;
6489371c9d4SSatish Balay     idx[1] = i;
6499371c9d4SSatish Balay     idx[2] = i + 1;
650c4762a1bSJed Brown     is     = i - mstart + 1;
651c4762a1bSJed Brown     v[0]   = sc * localptr[is];
652c4762a1bSJed Brown     v[1]   = sc * (localptr[is + 1] + localptr[is - 1] - 4.0 * localptr[is]);
653c4762a1bSJed Brown     v[2]   = sc * localptr[is];
6549566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, idx, v, INSERT_VALUES));
655c4762a1bSJed Brown   }
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   /*
658c4762a1bSJed Brown      Restore vector
659c4762a1bSJed Brown   */
6609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(local_in, &localptr));
661c4762a1bSJed Brown 
662c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
663c4762a1bSJed Brown      Complete the matrix assembly process and set some options
664c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
665c4762a1bSJed Brown   /*
666c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
667c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
668c4762a1bSJed Brown      Computations can be done while messages are in transition
669c4762a1bSJed Brown      by placing code between these two statements.
670c4762a1bSJed Brown   */
6719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
6729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
673c4762a1bSJed Brown 
674c4762a1bSJed Brown   /*
675c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
676c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
677c4762a1bSJed Brown   */
6789566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
679c4762a1bSJed Brown 
680*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
681c4762a1bSJed Brown }
682c4762a1bSJed Brown 
683c4762a1bSJed Brown /*TEST
684c4762a1bSJed Brown 
685c4762a1bSJed Brown     test:
686c4762a1bSJed Brown       args: -snes_type vinewtonrsls -ts_type glee -mymonitor -ts_max_steps 10 -nox
687c4762a1bSJed Brown       requires: !single
688c4762a1bSJed Brown 
689c4762a1bSJed Brown TEST*/
690