xref: /petsc/src/snes/tests/ex17.c (revision 226f8a8a5081bc6ad7227cd631662400f0d6e2a0)
1c4762a1bSJed Brown static const char help[] = "Newton's method to solve a two-variable system, sequentially.\n"
2c4762a1bSJed Brown                            "The same problem is solved twice - i) fully assembled system + ii) block system\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers.  Note that this
6c4762a1bSJed Brown file automatically includes:
7c4762a1bSJed Brown petscsys.h       - base PETSc routines   petscvec.h - vectors
8c4762a1bSJed Brown petscsys.h    - system routines       petscmat.h - matrices
9c4762a1bSJed Brown petscis.h     - index sets            petscksp.h - Krylov subspace methods
10c4762a1bSJed Brown petscviewer.h - viewers               petscpc.h  - preconditioners
11c4762a1bSJed Brown petscksp.h   - linear solvers
12c4762a1bSJed Brown */
13c4762a1bSJed Brown #include <petscsnes.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown /*
16c4762a1bSJed Brown This example is block version of the test found at
17c4762a1bSJed Brown   ${PETSC_DIR}/src/snes/tutorials/ex1.c
18c4762a1bSJed Brown In this test we replace the Jacobian systems
19c4762a1bSJed Brown   [J]{x} = {F}
20c4762a1bSJed Brown where
21c4762a1bSJed Brown 
22c4762a1bSJed Brown [J] = (j_00, j_01),  {x} = (x_0, x_1)^T,   {F} = (f_0, f_1)^T
23c4762a1bSJed Brown       (j_10, j_11)
24c4762a1bSJed Brown where [J] \in \mathbb^{2 \times 2}, {x},{F} \in \mathbb^{2 \times 1},
25c4762a1bSJed Brown 
26c4762a1bSJed Brown with a block system in which each block is of length 1.
27c4762a1bSJed Brown i.e. The block system is thus
28c4762a1bSJed Brown 
29c4762a1bSJed Brown [J] = ([j00], [j01]),  {x} = ({x0}, {x1})^T, {F} = ({f0}, {f1})^T
30c4762a1bSJed Brown       ([j10], [j11])
31c4762a1bSJed Brown where
32c4762a1bSJed Brown [j00], [j01], [j10], [j11] \in \mathbb^{1 \times 1}
33c4762a1bSJed Brown {x0}, {x1}, {f0}, {f1} \in \mathbb^{1 \times 1}
34c4762a1bSJed Brown 
35c4762a1bSJed Brown In practice we would not bother defing blocks of size one, and would instead assemble the
36c4762a1bSJed Brown full system. This is just a simple test to illustrate how to manipulate the blocks and
37c4762a1bSJed Brown to confirm the implementation is correct.
38c4762a1bSJed Brown */
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /*
41c4762a1bSJed Brown User-defined routines
42c4762a1bSJed Brown */
43c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES, Vec, Mat, Mat, void *);
44c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES, Vec, Vec, void *);
45c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES, Vec, Mat, Mat, void *);
46c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES, Vec, Vec, void *);
47c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES, Vec, Mat, Mat, void *);
48c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES, Vec, Vec, void *);
49c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES, Vec, Mat, Mat, void *);
50c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES, Vec, Vec, void *);
51c4762a1bSJed Brown 
assembled_system(void)52d71ae5a4SJacob Faibussowitsch static PetscErrorCode assembled_system(void)
53d71ae5a4SJacob Faibussowitsch {
54c4762a1bSJed Brown   SNES        snes; /* nonlinear solver context */
55c4762a1bSJed Brown   KSP         ksp;  /* linear solver context */
56c4762a1bSJed Brown   PC          pc;   /* preconditioner context */
57c4762a1bSJed Brown   Vec         x, r; /* solution, residual vectors */
58c4762a1bSJed Brown   Mat         J;    /* Jacobian matrix */
59c4762a1bSJed Brown   PetscInt    its;
60c4762a1bSJed Brown   PetscScalar pfive = .5, *xx;
61c4762a1bSJed Brown   PetscBool   flg;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   PetscFunctionBeginUser;
649566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n"));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4762a1bSJed Brown   Create nonlinear solver context
68c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73c4762a1bSJed Brown   Create matrix and vector data structures; set corresponding routines
74c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /*
77c4762a1bSJed Brown   Create vectors for solution and nonlinear function
78c4762a1bSJed Brown   */
799566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &x));
809566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /*
83c4762a1bSJed Brown   Create Jacobian matrix data structure
84c4762a1bSJed Brown   */
859566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
879566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
889566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg));
91c4762a1bSJed Brown   if (!flg) {
92c4762a1bSJed Brown     /*
93c4762a1bSJed Brown     Set function evaluation routine and vector.
94c4762a1bSJed Brown     */
959566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown     /*
98c4762a1bSJed Brown     Set Jacobian matrix data structure and Jacobian evaluation routine
99c4762a1bSJed Brown     */
1009566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL));
101c4762a1bSJed Brown   } else {
1029566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction2, NULL));
1039566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2, NULL));
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown   Customize nonlinear solver; set runtime options
108c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /*
111c4762a1bSJed Brown   Set linear solver defaults for this problem. By extracting the
112c4762a1bSJed Brown   KSP, KSP, and PC contexts from the SNES context, we can then
113c4762a1bSJed Brown   directly call any KSP, KSP, and PC routines to set various options.
114c4762a1bSJed Brown   */
1159566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
1169566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
1179566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
118fb842aefSJose E. Roman   PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_CURRENT, PETSC_CURRENT, 20));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*
121c4762a1bSJed Brown   Set SNES/KSP/KSP/PC runtime options, e.g.,
122c4762a1bSJed Brown   -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
123c4762a1bSJed Brown   These options will override those specified above as long as
124c4762a1bSJed Brown   SNESSetFromOptions() is called _after_ any other customization
125c4762a1bSJed Brown   routines.
126c4762a1bSJed Brown   */
1279566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown   Evaluate initial guess; then solve nonlinear system
131c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4762a1bSJed Brown   if (!flg) {
1339566063dSJacob Faibussowitsch     PetscCall(VecSet(x, pfive));
134c4762a1bSJed Brown   } else {
1359566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &xx));
1369371c9d4SSatish Balay     xx[0] = 2.0;
1379371c9d4SSatish Balay     xx[1] = 3.0;
1389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &xx));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown   /*
141c4762a1bSJed Brown   Note: The user should initialize the vector, x, with the initial guess
142c4762a1bSJed Brown   for the nonlinear solver prior to calling SNESSolve().  In particular,
143c4762a1bSJed Brown   to employ an initial guess of zero, the user should explicitly set
144c4762a1bSJed Brown   this vector to zero by calling VecSet().
145c4762a1bSJed Brown   */
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1489566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
149c4762a1bSJed Brown   if (flg) {
150c4762a1bSJed Brown     Vec f;
1519566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1529566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &f, 0, 0));
1539566063dSJacob Faibussowitsch     PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
154c4762a1bSJed Brown   }
15563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown   Free work space.  All PETSc objects should be destroyed when they
159c4762a1bSJed Brown   are no longer needed.
160c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161c4762a1bSJed Brown 
1629371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1639371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1649371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1659371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown /*
170c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x).
171c4762a1bSJed Brown 
172c4762a1bSJed Brown Input Parameters:
173c4762a1bSJed Brown .  snes - the SNES context
174c4762a1bSJed Brown .  x - input vector
175c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
176c4762a1bSJed Brown 
177c4762a1bSJed Brown Output Parameter:
178c4762a1bSJed Brown .  f - function vector
179c4762a1bSJed Brown */
FormFunction1(SNES snes,Vec x,Vec f,void * dummy)180d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *dummy)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown   const PetscScalar *xx;
183c4762a1bSJed Brown   PetscScalar       *ff;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   PetscFunctionBeginUser;
186c4762a1bSJed Brown   /*
187c4762a1bSJed Brown   Get pointers to vector data.
188c4762a1bSJed Brown   - For default PETSc vectors, VecGetArray() returns a pointer to
189c4762a1bSJed Brown   the data array.  Otherwise, the routine is implementation dependent.
190c4762a1bSJed Brown   - You MUST call VecRestoreArray() when you no longer need access to
191c4762a1bSJed Brown   the array.
192c4762a1bSJed Brown   */
1939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /*
197c4762a1bSJed Brown   Compute function
198c4762a1bSJed Brown   */
199c4762a1bSJed Brown   ff[0] = xx[0] * xx[0] + xx[0] * xx[1] - 3.0;
200c4762a1bSJed Brown   ff[1] = xx[0] * xx[1] + xx[1] * xx[1] - 6.0;
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /*
203c4762a1bSJed Brown   Restore vectors
204c4762a1bSJed Brown   */
2059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208c4762a1bSJed Brown }
209f6dfbefdSBarry Smith 
210c4762a1bSJed Brown /*
211c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix.
212c4762a1bSJed Brown 
213c4762a1bSJed Brown Input Parameters:
214c4762a1bSJed Brown .  snes - the SNES context
215c4762a1bSJed Brown .  x - input vector
216c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
217c4762a1bSJed Brown 
218c4762a1bSJed Brown Output Parameters:
219c4762a1bSJed Brown .  jac - Jacobian matrix
220*7addb90fSBarry Smith .  B - optionally different matrix used to construct the preconditioner
2210b4b7b1cSBarry Smith 
222c4762a1bSJed Brown */
FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void * dummy)223d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
224d71ae5a4SJacob Faibussowitsch {
225c4762a1bSJed Brown   const PetscScalar *xx;
226c4762a1bSJed Brown   PetscScalar        A[4];
227c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   PetscFunctionBeginUser;
230c4762a1bSJed Brown   /*
231c4762a1bSJed Brown   Get pointer to vector data
232c4762a1bSJed Brown   */
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /*
236c4762a1bSJed Brown   Compute Jacobian entries and insert into matrix.
237c4762a1bSJed Brown   - Since this is such a small problem, we set all entries for
238c4762a1bSJed Brown   the matrix at once.
239c4762a1bSJed Brown   */
2409371c9d4SSatish Balay   A[0] = 2.0 * xx[0] + xx[1];
2419371c9d4SSatish Balay   A[1] = xx[0];
2429371c9d4SSatish Balay   A[2] = xx[1];
2439371c9d4SSatish Balay   A[3] = xx[0] + 2.0 * xx[1];
2449566063dSJacob Faibussowitsch   PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES));
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /*
247c4762a1bSJed Brown   Restore vector
248c4762a1bSJed Brown   */
2499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /*
252c4762a1bSJed Brown   Assemble matrix
253c4762a1bSJed Brown   */
2549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
257c4762a1bSJed Brown }
258c4762a1bSJed Brown 
FormFunction2(SNES snes,Vec x,Vec f,void * dummy)259d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunction2(SNES snes, Vec x, Vec f, void *dummy)
260d71ae5a4SJacob Faibussowitsch {
261c4762a1bSJed Brown   const PetscScalar *xx;
262c4762a1bSJed Brown   PetscScalar       *ff;
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   PetscFunctionBeginUser;
265c4762a1bSJed Brown   /*
266c4762a1bSJed Brown   Get pointers to vector data.
267c4762a1bSJed Brown   - For default PETSc vectors, VecGetArray() returns a pointer to
268c4762a1bSJed Brown   the data array.  Otherwise, the routine is implementation dependent.
269c4762a1bSJed Brown   - You MUST call VecRestoreArray() when you no longer need access to
270c4762a1bSJed Brown   the array.
271c4762a1bSJed Brown   */
2729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /*
276c4762a1bSJed Brown   Compute function
277c4762a1bSJed Brown   */
278c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0];
279c4762a1bSJed Brown   ff[1] = xx[1];
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   /*
282c4762a1bSJed Brown   Restore vectors
283c4762a1bSJed Brown   */
2849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287c4762a1bSJed Brown }
288c4762a1bSJed Brown 
FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void * dummy)289d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobian2(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
290d71ae5a4SJacob Faibussowitsch {
291c4762a1bSJed Brown   const PetscScalar *xx;
292c4762a1bSJed Brown   PetscScalar        A[4];
293c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   PetscFunctionBeginUser;
296c4762a1bSJed Brown   /*
297c4762a1bSJed Brown   Get pointer to vector data
298c4762a1bSJed Brown   */
2999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   /*
302c4762a1bSJed Brown   Compute Jacobian entries and insert into matrix.
303c4762a1bSJed Brown   - Since this is such a small problem, we set all entries for
304c4762a1bSJed Brown   the matrix at once.
305c4762a1bSJed Brown   */
3069371c9d4SSatish Balay   A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0;
3079371c9d4SSatish Balay   A[1] = 0.0;
3089371c9d4SSatish Balay   A[2] = 0.0;
3099371c9d4SSatish Balay   A[3] = 1.0;
3109566063dSJacob Faibussowitsch   PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /*
313c4762a1bSJed Brown   Restore vector
314c4762a1bSJed Brown   */
3159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /*
318c4762a1bSJed Brown   Assemble matrix
319c4762a1bSJed Brown   */
3209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
3219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
block_system(void)325d71ae5a4SJacob Faibussowitsch static PetscErrorCode block_system(void)
326d71ae5a4SJacob Faibussowitsch {
327c4762a1bSJed Brown   SNES        snes; /* nonlinear solver context */
328c4762a1bSJed Brown   KSP         ksp;  /* linear solver context */
329c4762a1bSJed Brown   PC          pc;   /* preconditioner context */
330c4762a1bSJed Brown   Vec         x, r; /* solution, residual vectors */
331c4762a1bSJed Brown   Mat         J;    /* Jacobian matrix */
332c4762a1bSJed Brown   PetscInt    its;
333c4762a1bSJed Brown   PetscScalar pfive = .5;
334c4762a1bSJed Brown   PetscBool   flg;
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   Mat j11, j12, j21, j22;
337c4762a1bSJed Brown   Vec x1, x2, r1, r2;
338c4762a1bSJed Brown   Vec bv;
339c4762a1bSJed Brown   Vec bx[2];
340c4762a1bSJed Brown   Mat bA[2][2];
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   PetscFunctionBeginUser;
3439566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n"));
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
348c4762a1bSJed Brown   Create matrix and vector data structures; set corresponding routines
349c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   /*
352c4762a1bSJed Brown   Create sub vectors for solution and nonlinear function
353c4762a1bSJed Brown   */
3549566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &x1));
3559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x1, &r1));
356c4762a1bSJed Brown 
3579566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &x2));
3589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x2, &r2));
359c4762a1bSJed Brown 
360c4762a1bSJed Brown   /*
361c4762a1bSJed Brown   Create the block vectors
362c4762a1bSJed Brown   */
363c4762a1bSJed Brown   bx[0] = x1;
364c4762a1bSJed Brown   bx[1] = x2;
3659566063dSJacob Faibussowitsch   PetscCall(VecCreateNest(PETSC_COMM_WORLD, 2, NULL, bx, &x));
3669566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
3679566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
3689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x1));
3699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x2));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   bx[0] = r1;
372c4762a1bSJed Brown   bx[1] = r2;
3739566063dSJacob Faibussowitsch   PetscCall(VecCreateNest(PETSC_COMM_WORLD, 2, NULL, bx, &r));
3749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r1));
3759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r2));
3769566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(r));
3779566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(r));
378c4762a1bSJed Brown 
379c4762a1bSJed Brown   /*
380c4762a1bSJed Brown   Create sub Jacobian matrix data structure
381c4762a1bSJed Brown   */
3829566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j11));
3839566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j11, 1, 1, 1, 1));
3849566063dSJacob Faibussowitsch   PetscCall(MatSetType(j11, MATSEQAIJ));
3859566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j11));
386c4762a1bSJed Brown 
3879566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j12));
3889566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j12, 1, 1, 1, 1));
3899566063dSJacob Faibussowitsch   PetscCall(MatSetType(j12, MATSEQAIJ));
3909566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j12));
391c4762a1bSJed Brown 
3929566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j21));
3939566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j21, 1, 1, 1, 1));
3949566063dSJacob Faibussowitsch   PetscCall(MatSetType(j21, MATSEQAIJ));
3959566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j21));
396c4762a1bSJed Brown 
3979566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j22));
3989566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
3999566063dSJacob Faibussowitsch   PetscCall(MatSetType(j22, MATSEQAIJ));
4009566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j22));
401c4762a1bSJed Brown   /*
402c4762a1bSJed Brown   Create block Jacobian matrix data structure
403c4762a1bSJed Brown   */
404c4762a1bSJed Brown   bA[0][0] = j11;
405c4762a1bSJed Brown   bA[0][1] = j12;
406c4762a1bSJed Brown   bA[1][0] = j21;
407c4762a1bSJed Brown   bA[1][1] = j22;
408c4762a1bSJed Brown 
4099566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, &bA[0][0], &J));
4109566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
4119566063dSJacob Faibussowitsch   PetscCall(MatNestSetVecType(J, VECNEST));
4129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j11));
4139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j12));
4149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j21));
4159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j22));
416c4762a1bSJed Brown 
4179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
419c4762a1bSJed Brown 
4209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg));
421c4762a1bSJed Brown   if (!flg) {
422c4762a1bSJed Brown     /*
423c4762a1bSJed Brown     Set function evaluation routine and vector.
424c4762a1bSJed Brown     */
4259566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction1_block, NULL));
426c4762a1bSJed Brown 
427c4762a1bSJed Brown     /*
428c4762a1bSJed Brown     Set Jacobian matrix data structure and Jacobian evaluation routine
429c4762a1bSJed Brown     */
4309566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1_block, NULL));
431c4762a1bSJed Brown   } else {
4329566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction2_block, NULL));
4339566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2_block, NULL));
434c4762a1bSJed Brown   }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
437c4762a1bSJed Brown   Customize nonlinear solver; set runtime options
438c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /*
441c4762a1bSJed Brown   Set linear solver defaults for this problem. By extracting the
442c4762a1bSJed Brown   KSP, KSP, and PC contexts from the SNES context, we can then
443c4762a1bSJed Brown   directly call any KSP, KSP, and PC routines to set various options.
444c4762a1bSJed Brown   */
4459566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
4469566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
4479566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
448fb842aefSJose E. Roman   PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_CURRENT, PETSC_CURRENT, 20));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   /*
451c4762a1bSJed Brown   Set SNES/KSP/KSP/PC runtime options, e.g.,
452c4762a1bSJed Brown   -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
453c4762a1bSJed Brown   These options will override those specified above as long as
454c4762a1bSJed Brown   SNESSetFromOptions() is called _after_ any other customization
455c4762a1bSJed Brown   routines.
456c4762a1bSJed Brown   */
4579566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
458c4762a1bSJed Brown 
459c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
460c4762a1bSJed Brown   Evaluate initial guess; then solve nonlinear system
461c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
462c4762a1bSJed Brown   if (!flg) {
4639566063dSJacob Faibussowitsch     PetscCall(VecSet(x, pfive));
464c4762a1bSJed Brown   } else {
465c4762a1bSJed Brown     Vec *vecs;
4669566063dSJacob Faibussowitsch     PetscCall(VecNestGetSubVecs(x, NULL, &vecs));
467c4762a1bSJed Brown     bv = vecs[0];
4689566063dSJacob Faibussowitsch     /*    PetscCall(VecBlockGetSubVec(x, 0, &bv)); */
4699566063dSJacob Faibussowitsch     PetscCall(VecSetValue(bv, 0, 2.0, INSERT_VALUES)); /* xx[0] = 2.0; */
4709566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(bv));
4719566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(bv));
472c4762a1bSJed Brown 
4739566063dSJacob Faibussowitsch     /*    PetscCall(VecBlockGetSubVec(x, 1, &bv)); */
474c4762a1bSJed Brown     bv = vecs[1];
4759566063dSJacob Faibussowitsch     PetscCall(VecSetValue(bv, 0, 3.0, INSERT_VALUES)); /* xx[1] = 3.0; */
4769566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(bv));
4779566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(bv));
478c4762a1bSJed Brown   }
479c4762a1bSJed Brown   /*
480c4762a1bSJed Brown   Note: The user should initialize the vector, x, with the initial guess
481c4762a1bSJed Brown   for the nonlinear solver prior to calling SNESSolve().  In particular,
482c4762a1bSJed Brown   to employ an initial guess of zero, the user should explicitly set
483c4762a1bSJed Brown   this vector to zero by calling VecSet().
484c4762a1bSJed Brown   */
4859566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
4869566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
487c4762a1bSJed Brown   if (flg) {
488c4762a1bSJed Brown     Vec f;
4899566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
4909566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &f, 0, 0));
4919566063dSJacob Faibussowitsch     PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
492c4762a1bSJed Brown   }
493c4762a1bSJed Brown 
49463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
497c4762a1bSJed Brown   Free work space.  All PETSc objects should be destroyed when they
498c4762a1bSJed Brown   are no longer needed.
499c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
5009371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
5019371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
5029371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
5039371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
5043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
505c4762a1bSJed Brown }
506c4762a1bSJed Brown 
FormFunction1_block(SNES snes,Vec x,Vec f,void * dummy)507d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunction1_block(SNES snes, Vec x, Vec f, void *dummy)
508d71ae5a4SJacob Faibussowitsch {
509c4762a1bSJed Brown   Vec        *xx, *ff, x1, x2, f1, f2;
510c4762a1bSJed Brown   PetscScalar ff_0, ff_1;
511c4762a1bSJed Brown   PetscScalar xx_0, xx_1;
512c4762a1bSJed Brown   PetscInt    index, nb;
513c4762a1bSJed Brown 
514c4762a1bSJed Brown   PetscFunctionBeginUser;
515c4762a1bSJed Brown   /* get blocks for function */
5169566063dSJacob Faibussowitsch   PetscCall(VecNestGetSubVecs(f, &nb, &ff));
5179371c9d4SSatish Balay   f1 = ff[0];
5189371c9d4SSatish Balay   f2 = ff[1];
519c4762a1bSJed Brown 
520c4762a1bSJed Brown   /* get blocks for solution */
5219566063dSJacob Faibussowitsch   PetscCall(VecNestGetSubVecs(x, &nb, &xx));
5229371c9d4SSatish Balay   x1 = xx[0];
5239371c9d4SSatish Balay   x2 = xx[1];
524c4762a1bSJed Brown 
525c4762a1bSJed Brown   /* get solution values */
526c4762a1bSJed Brown   index = 0;
5279566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x1, 1, &index, &xx_0));
5289566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x2, 1, &index, &xx_1));
529c4762a1bSJed Brown 
530c4762a1bSJed Brown   /* Compute function */
531c4762a1bSJed Brown   ff_0 = xx_0 * xx_0 + xx_0 * xx_1 - 3.0;
532c4762a1bSJed Brown   ff_1 = xx_0 * xx_1 + xx_1 * xx_1 - 6.0;
533c4762a1bSJed Brown 
534c4762a1bSJed Brown   /* set function values */
5359566063dSJacob Faibussowitsch   PetscCall(VecSetValue(f1, index, ff_0, INSERT_VALUES));
536c4762a1bSJed Brown 
5379566063dSJacob Faibussowitsch   PetscCall(VecSetValue(f2, index, ff_1, INSERT_VALUES));
538c4762a1bSJed Brown 
5399566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(f));
5409566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(f));
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
542c4762a1bSJed Brown }
543c4762a1bSJed Brown 
FormJacobian1_block(SNES snes,Vec x,Mat jac,Mat B,void * dummy)544d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobian1_block(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
545d71ae5a4SJacob Faibussowitsch {
546c4762a1bSJed Brown   Vec        *xx, x1, x2;
547c4762a1bSJed Brown   PetscScalar xx_0, xx_1;
548c4762a1bSJed Brown   PetscInt    index, nb;
549c4762a1bSJed Brown   PetscScalar A_00, A_01, A_10, A_11;
550c4762a1bSJed Brown   Mat         j11, j12, j21, j22;
551c4762a1bSJed Brown   Mat       **mats;
552c4762a1bSJed Brown 
553c4762a1bSJed Brown   PetscFunctionBeginUser;
554c4762a1bSJed Brown   /* get blocks for solution */
5559566063dSJacob Faibussowitsch   PetscCall(VecNestGetSubVecs(x, &nb, &xx));
5569371c9d4SSatish Balay   x1 = xx[0];
5579371c9d4SSatish Balay   x2 = xx[1];
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   /* get solution values */
560c4762a1bSJed Brown   index = 0;
5619566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x1, 1, &index, &xx_0));
5629566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x2, 1, &index, &xx_1));
563c4762a1bSJed Brown 
564c4762a1bSJed Brown   /* get block matrices */
5659566063dSJacob Faibussowitsch   PetscCall(MatNestGetSubMats(jac, NULL, NULL, &mats));
566c4762a1bSJed Brown   j11 = mats[0][0];
567c4762a1bSJed Brown   j12 = mats[0][1];
568c4762a1bSJed Brown   j21 = mats[1][0];
569c4762a1bSJed Brown   j22 = mats[1][1];
570c4762a1bSJed Brown 
571c4762a1bSJed Brown   /* compute jacobian entries */
572c4762a1bSJed Brown   A_00 = 2.0 * xx_0 + xx_1;
573c4762a1bSJed Brown   A_01 = xx_0;
574c4762a1bSJed Brown   A_10 = xx_1;
575c4762a1bSJed Brown   A_11 = xx_0 + 2.0 * xx_1;
576c4762a1bSJed Brown 
577c4762a1bSJed Brown   /* set jacobian values */
5789566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j11, 0, 0, A_00, INSERT_VALUES));
5799566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j12, 0, 0, A_01, INSERT_VALUES));
5809566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j21, 0, 0, A_10, INSERT_VALUES));
5819566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j22, 0, 0, A_11, INSERT_VALUES));
582c4762a1bSJed Brown 
583c4762a1bSJed Brown   /* Assemble sub matrix */
5849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
5859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
587c4762a1bSJed Brown }
588c4762a1bSJed Brown 
FormFunction2_block(SNES snes,Vec x,Vec f,void * dummy)589d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunction2_block(SNES snes, Vec x, Vec f, void *dummy)
590d71ae5a4SJacob Faibussowitsch {
591c4762a1bSJed Brown   PetscScalar       *ff;
592c4762a1bSJed Brown   const PetscScalar *xx;
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   PetscFunctionBeginUser;
595c4762a1bSJed Brown   /*
596c4762a1bSJed Brown   Get pointers to vector data.
597c4762a1bSJed Brown   - For default PETSc vectors, VecGetArray() returns a pointer to
598c4762a1bSJed Brown   the data array.  Otherwise, the routine is implementation dependent.
599c4762a1bSJed Brown   - You MUST call VecRestoreArray() when you no longer need access to
600c4762a1bSJed Brown   the array.
601c4762a1bSJed Brown   */
6029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
6039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   /*
606c4762a1bSJed Brown   Compute function
607c4762a1bSJed Brown   */
608c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0];
609c4762a1bSJed Brown   ff[1] = xx[1];
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   /*
612c4762a1bSJed Brown   Restore vectors
613c4762a1bSJed Brown   */
6149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
6159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617c4762a1bSJed Brown }
618c4762a1bSJed Brown 
FormJacobian2_block(SNES snes,Vec x,Mat jac,Mat B,void * dummy)619d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobian2_block(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
620d71ae5a4SJacob Faibussowitsch {
621c4762a1bSJed Brown   const PetscScalar *xx;
622c4762a1bSJed Brown   PetscScalar        A[4];
623c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
624c4762a1bSJed Brown 
625c4762a1bSJed Brown   PetscFunctionBeginUser;
626c4762a1bSJed Brown   /*
627c4762a1bSJed Brown   Get pointer to vector data
628c4762a1bSJed Brown   */
6299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
630c4762a1bSJed Brown 
631c4762a1bSJed Brown   /*
632c4762a1bSJed Brown   Compute Jacobian entries and insert into matrix.
633c4762a1bSJed Brown   - Since this is such a small problem, we set all entries for
634c4762a1bSJed Brown   the matrix at once.
635c4762a1bSJed Brown   */
6369371c9d4SSatish Balay   A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0;
6379371c9d4SSatish Balay   A[1] = 0.0;
6389371c9d4SSatish Balay   A[2] = 0.0;
6399371c9d4SSatish Balay   A[3] = 1.0;
6409566063dSJacob Faibussowitsch   PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES));
641c4762a1bSJed Brown 
642c4762a1bSJed Brown   /*
643c4762a1bSJed Brown   Restore vector
644c4762a1bSJed Brown   */
6459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
646c4762a1bSJed Brown 
647c4762a1bSJed Brown   /*
648c4762a1bSJed Brown   Assemble matrix
649c4762a1bSJed Brown   */
6509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
6519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653c4762a1bSJed Brown }
654c4762a1bSJed Brown 
main(int argc,char ** argv)655d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
656d71ae5a4SJacob Faibussowitsch {
657c4762a1bSJed Brown   PetscMPIInt size;
658c4762a1bSJed Brown 
659327415f7SBarry Smith   PetscFunctionBeginUser;
660c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
6619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
662be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
6639566063dSJacob Faibussowitsch   PetscCall(assembled_system());
6649566063dSJacob Faibussowitsch   PetscCall(block_system());
6659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
666b122ec5aSJacob Faibussowitsch   return 0;
667c4762a1bSJed Brown }
668c4762a1bSJed Brown 
669c4762a1bSJed Brown /*TEST
670c4762a1bSJed Brown 
671c4762a1bSJed Brown    test:
672c4762a1bSJed Brown       args: -snes_monitor_short
673c4762a1bSJed Brown       requires: !single
674c4762a1bSJed Brown 
675c4762a1bSJed Brown TEST*/
676