xref: /petsc/src/snes/tests/ex17.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
529371c9d4SSatish Balay static PetscErrorCode assembled_system(void) {
53c4762a1bSJed Brown   SNES        snes; /* nonlinear solver context */
54c4762a1bSJed Brown   KSP         ksp;  /* linear solver context */
55c4762a1bSJed Brown   PC          pc;   /* preconditioner context */
56c4762a1bSJed Brown   Vec         x, r; /* solution, residual vectors */
57c4762a1bSJed Brown   Mat         J;    /* Jacobian matrix */
58c4762a1bSJed Brown   PetscInt    its;
59c4762a1bSJed Brown   PetscScalar pfive = .5, *xx;
60c4762a1bSJed Brown   PetscBool   flg;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBeginUser;
639566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n"));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown   Create nonlinear solver context
67c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72c4762a1bSJed Brown   Create matrix and vector data structures; set corresponding routines
73c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /*
76c4762a1bSJed Brown   Create vectors for solution and nonlinear function
77c4762a1bSJed Brown   */
789566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &x));
799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /*
82c4762a1bSJed Brown   Create Jacobian matrix data structure
83c4762a1bSJed Brown   */
849566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
859566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
869566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
879566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg));
90c4762a1bSJed Brown   if (!flg) {
91c4762a1bSJed Brown     /*
92c4762a1bSJed Brown     Set function evaluation routine and vector.
93c4762a1bSJed Brown     */
949566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown     /*
97c4762a1bSJed Brown     Set Jacobian matrix data structure and Jacobian evaluation routine
98c4762a1bSJed Brown     */
999566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL));
100c4762a1bSJed Brown   } else {
1019566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction2, NULL));
1029566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2, NULL));
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown   Customize nonlinear solver; set runtime options
107c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /*
110c4762a1bSJed Brown   Set linear solver defaults for this problem. By extracting the
111c4762a1bSJed Brown   KSP, KSP, and PC contexts from the SNES context, we can then
112c4762a1bSJed Brown   directly call any KSP, KSP, and PC routines to set various options.
113c4762a1bSJed Brown   */
1149566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
1159566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
1169566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
1179566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /*
120c4762a1bSJed Brown   Set SNES/KSP/KSP/PC runtime options, e.g.,
121c4762a1bSJed Brown   -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
122c4762a1bSJed Brown   These options will override those specified above as long as
123c4762a1bSJed Brown   SNESSetFromOptions() is called _after_ any other customization
124c4762a1bSJed Brown   routines.
125c4762a1bSJed Brown   */
1269566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129c4762a1bSJed Brown   Evaluate initial guess; then solve nonlinear system
130c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131c4762a1bSJed Brown   if (!flg) {
1329566063dSJacob Faibussowitsch     PetscCall(VecSet(x, pfive));
133c4762a1bSJed Brown   } else {
1349566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &xx));
1359371c9d4SSatish Balay     xx[0] = 2.0;
1369371c9d4SSatish Balay     xx[1] = 3.0;
1379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &xx));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown   /*
140c4762a1bSJed Brown   Note: The user should initialize the vector, x, with the initial guess
141c4762a1bSJed Brown   for the nonlinear solver prior to calling SNESSolve().  In particular,
142c4762a1bSJed Brown   to employ an initial guess of zero, the user should explicitly set
143c4762a1bSJed Brown   this vector to zero by calling VecSet().
144c4762a1bSJed Brown   */
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
1479566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
148c4762a1bSJed Brown   if (flg) {
149c4762a1bSJed Brown     Vec f;
1509566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1519566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &f, 0, 0));
1529566063dSJacob Faibussowitsch     PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
153c4762a1bSJed Brown   }
15463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157c4762a1bSJed Brown   Free work space.  All PETSc objects should be destroyed when they
158c4762a1bSJed Brown   are no longer needed.
159c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160c4762a1bSJed Brown 
1619371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1629371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
1639371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
1649371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
165c4762a1bSJed Brown   PetscFunctionReturn(0);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown /*
169c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x).
170c4762a1bSJed Brown 
171c4762a1bSJed Brown Input Parameters:
172c4762a1bSJed Brown .  snes - the SNES context
173c4762a1bSJed Brown .  x - input vector
174c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
175c4762a1bSJed Brown 
176c4762a1bSJed Brown Output Parameter:
177c4762a1bSJed Brown .  f - function vector
178c4762a1bSJed Brown */
1799371c9d4SSatish Balay static PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *dummy) {
180c4762a1bSJed Brown   const PetscScalar *xx;
181c4762a1bSJed Brown   PetscScalar       *ff;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   PetscFunctionBeginUser;
184c4762a1bSJed Brown   /*
185c4762a1bSJed Brown   Get pointers to vector data.
186c4762a1bSJed Brown   - For default PETSc vectors, VecGetArray() returns a pointer to
187c4762a1bSJed Brown   the data array.  Otherwise, the routine is implementation dependent.
188c4762a1bSJed Brown   - You MUST call VecRestoreArray() when you no longer need access to
189c4762a1bSJed Brown   the array.
190c4762a1bSJed Brown   */
1919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /*
195c4762a1bSJed Brown   Compute function
196c4762a1bSJed Brown   */
197c4762a1bSJed Brown   ff[0] = xx[0] * xx[0] + xx[0] * xx[1] - 3.0;
198c4762a1bSJed Brown   ff[1] = xx[0] * xx[1] + xx[1] * xx[1] - 6.0;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /*
201c4762a1bSJed Brown   Restore vectors
202c4762a1bSJed Brown   */
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207*f6dfbefdSBarry Smith 
208c4762a1bSJed Brown /*
209c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix.
210c4762a1bSJed Brown 
211c4762a1bSJed Brown Input Parameters:
212c4762a1bSJed Brown .  snes - the SNES context
213c4762a1bSJed Brown .  x - input vector
214c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
215c4762a1bSJed Brown 
216c4762a1bSJed Brown Output Parameters:
217c4762a1bSJed Brown .  jac - Jacobian matrix
218c4762a1bSJed Brown .  B - optionally different preconditioning matrix
219c4762a1bSJed Brown .  flag - flag indicating matrix structure
220c4762a1bSJed Brown */
2219371c9d4SSatish Balay static PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy) {
222c4762a1bSJed Brown   const PetscScalar *xx;
223c4762a1bSJed Brown   PetscScalar        A[4];
224c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   PetscFunctionBeginUser;
227c4762a1bSJed Brown   /*
228c4762a1bSJed Brown   Get pointer to vector data
229c4762a1bSJed Brown   */
2309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown   Compute Jacobian entries and insert into matrix.
234c4762a1bSJed Brown   - Since this is such a small problem, we set all entries for
235c4762a1bSJed Brown   the matrix at once.
236c4762a1bSJed Brown   */
2379371c9d4SSatish Balay   A[0] = 2.0 * xx[0] + xx[1];
2389371c9d4SSatish Balay   A[1] = xx[0];
2399371c9d4SSatish Balay   A[2] = xx[1];
2409371c9d4SSatish Balay   A[3] = xx[0] + 2.0 * xx[1];
2419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /*
244c4762a1bSJed Brown   Restore vector
245c4762a1bSJed Brown   */
2469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /*
249c4762a1bSJed Brown   Assemble matrix
250c4762a1bSJed Brown   */
2519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
2569371c9d4SSatish Balay static PetscErrorCode FormFunction2(SNES snes, Vec x, Vec f, void *dummy) {
257c4762a1bSJed Brown   const PetscScalar *xx;
258c4762a1bSJed Brown   PetscScalar       *ff;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   PetscFunctionBeginUser;
261c4762a1bSJed Brown   /*
262c4762a1bSJed Brown   Get pointers to vector data.
263c4762a1bSJed Brown   - For default PETSc vectors, VecGetArray() returns a pointer to
264c4762a1bSJed Brown   the data array.  Otherwise, the routine is implementation dependent.
265c4762a1bSJed Brown   - You MUST call VecRestoreArray() when you no longer need access to
266c4762a1bSJed Brown   the array.
267c4762a1bSJed Brown   */
2689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /*
272c4762a1bSJed Brown   Compute function
273c4762a1bSJed Brown   */
274c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0];
275c4762a1bSJed Brown   ff[1] = xx[1];
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /*
278c4762a1bSJed Brown   Restore vectors
279c4762a1bSJed Brown   */
2809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
282c4762a1bSJed Brown   PetscFunctionReturn(0);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
2859371c9d4SSatish Balay static PetscErrorCode FormJacobian2(SNES snes, Vec x, Mat jac, Mat B, void *dummy) {
286c4762a1bSJed Brown   const PetscScalar *xx;
287c4762a1bSJed Brown   PetscScalar        A[4];
288c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   PetscFunctionBeginUser;
291c4762a1bSJed Brown   /*
292c4762a1bSJed Brown   Get pointer to vector data
293c4762a1bSJed Brown   */
2949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   /*
297c4762a1bSJed Brown   Compute Jacobian entries and insert into matrix.
298c4762a1bSJed Brown   - Since this is such a small problem, we set all entries for
299c4762a1bSJed Brown   the matrix at once.
300c4762a1bSJed Brown   */
3019371c9d4SSatish Balay   A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0;
3029371c9d4SSatish Balay   A[1] = 0.0;
3039371c9d4SSatish Balay   A[2] = 0.0;
3049371c9d4SSatish Balay   A[3] = 1.0;
3059566063dSJacob Faibussowitsch   PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /*
308c4762a1bSJed Brown   Restore vector
309c4762a1bSJed Brown   */
3109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /*
313c4762a1bSJed Brown   Assemble matrix
314c4762a1bSJed Brown   */
3159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
3169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
317c4762a1bSJed Brown   PetscFunctionReturn(0);
318c4762a1bSJed Brown }
319c4762a1bSJed Brown 
3209371c9d4SSatish Balay static PetscErrorCode block_system(void) {
321c4762a1bSJed Brown   SNES        snes; /* nonlinear solver context */
322c4762a1bSJed Brown   KSP         ksp;  /* linear solver context */
323c4762a1bSJed Brown   PC          pc;   /* preconditioner context */
324c4762a1bSJed Brown   Vec         x, r; /* solution, residual vectors */
325c4762a1bSJed Brown   Mat         J;    /* Jacobian matrix */
326c4762a1bSJed Brown   PetscInt    its;
327c4762a1bSJed Brown   PetscScalar pfive = .5;
328c4762a1bSJed Brown   PetscBool   flg;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   Mat j11, j12, j21, j22;
331c4762a1bSJed Brown   Vec x1, x2, r1, r2;
332c4762a1bSJed Brown   Vec bv;
333c4762a1bSJed Brown   Vec bx[2];
334c4762a1bSJed Brown   Mat bA[2][2];
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   PetscFunctionBeginUser;
3379566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n"));
338c4762a1bSJed Brown 
3399566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
342c4762a1bSJed Brown   Create matrix and vector data structures; set corresponding routines
343c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
344c4762a1bSJed Brown 
345c4762a1bSJed Brown   /*
346c4762a1bSJed Brown   Create sub vectors for solution and nonlinear function
347c4762a1bSJed Brown   */
3489566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &x1));
3499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x1, &r1));
350c4762a1bSJed Brown 
3519566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &x2));
3529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x2, &r2));
353c4762a1bSJed Brown 
354c4762a1bSJed Brown   /*
355c4762a1bSJed Brown   Create the block vectors
356c4762a1bSJed Brown   */
357c4762a1bSJed Brown   bx[0] = x1;
358c4762a1bSJed Brown   bx[1] = x2;
3599566063dSJacob Faibussowitsch   PetscCall(VecCreateNest(PETSC_COMM_WORLD, 2, NULL, bx, &x));
3609566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
3619566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
3629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x1));
3639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x2));
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   bx[0] = r1;
366c4762a1bSJed Brown   bx[1] = r2;
3679566063dSJacob Faibussowitsch   PetscCall(VecCreateNest(PETSC_COMM_WORLD, 2, NULL, bx, &r));
3689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r1));
3699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r2));
3709566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(r));
3719566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(r));
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   /*
374c4762a1bSJed Brown   Create sub Jacobian matrix data structure
375c4762a1bSJed Brown   */
3769566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j11));
3779566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j11, 1, 1, 1, 1));
3789566063dSJacob Faibussowitsch   PetscCall(MatSetType(j11, MATSEQAIJ));
3799566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j11));
380c4762a1bSJed Brown 
3819566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j12));
3829566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j12, 1, 1, 1, 1));
3839566063dSJacob Faibussowitsch   PetscCall(MatSetType(j12, MATSEQAIJ));
3849566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j12));
385c4762a1bSJed Brown 
3869566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j21));
3879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j21, 1, 1, 1, 1));
3889566063dSJacob Faibussowitsch   PetscCall(MatSetType(j21, MATSEQAIJ));
3899566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j21));
390c4762a1bSJed Brown 
3919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &j22));
3929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
3939566063dSJacob Faibussowitsch   PetscCall(MatSetType(j22, MATSEQAIJ));
3949566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j22));
395c4762a1bSJed Brown   /*
396c4762a1bSJed Brown   Create block Jacobian matrix data structure
397c4762a1bSJed Brown   */
398c4762a1bSJed Brown   bA[0][0] = j11;
399c4762a1bSJed Brown   bA[0][1] = j12;
400c4762a1bSJed Brown   bA[1][0] = j21;
401c4762a1bSJed Brown   bA[1][1] = j22;
402c4762a1bSJed Brown 
4039566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, &bA[0][0], &J));
4049566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
4059566063dSJacob Faibussowitsch   PetscCall(MatNestSetVecType(J, VECNEST));
4069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j11));
4079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j12));
4089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j21));
4099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j22));
410c4762a1bSJed Brown 
4119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
413c4762a1bSJed Brown 
4149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg));
415c4762a1bSJed Brown   if (!flg) {
416c4762a1bSJed Brown     /*
417c4762a1bSJed Brown     Set function evaluation routine and vector.
418c4762a1bSJed Brown     */
4199566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction1_block, NULL));
420c4762a1bSJed Brown 
421c4762a1bSJed Brown     /*
422c4762a1bSJed Brown     Set Jacobian matrix data structure and Jacobian evaluation routine
423c4762a1bSJed Brown     */
4249566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1_block, NULL));
425c4762a1bSJed Brown   } else {
4269566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, r, FormFunction2_block, NULL));
4279566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2_block, NULL));
428c4762a1bSJed Brown   }
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
431c4762a1bSJed Brown   Customize nonlinear solver; set runtime options
432c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /*
435c4762a1bSJed Brown   Set linear solver defaults for this problem. By extracting the
436c4762a1bSJed Brown   KSP, KSP, and PC contexts from the SNES context, we can then
437c4762a1bSJed Brown   directly call any KSP, KSP, and PC routines to set various options.
438c4762a1bSJed Brown   */
4399566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
4409566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
4419566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
4429566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20));
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   /*
445c4762a1bSJed Brown   Set SNES/KSP/KSP/PC runtime options, e.g.,
446c4762a1bSJed Brown   -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
447c4762a1bSJed Brown   These options will override those specified above as long as
448c4762a1bSJed Brown   SNESSetFromOptions() is called _after_ any other customization
449c4762a1bSJed Brown   routines.
450c4762a1bSJed Brown   */
4519566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
454c4762a1bSJed Brown   Evaluate initial guess; then solve nonlinear system
455c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
456c4762a1bSJed Brown   if (!flg) {
4579566063dSJacob Faibussowitsch     PetscCall(VecSet(x, pfive));
458c4762a1bSJed Brown   } else {
459c4762a1bSJed Brown     Vec *vecs;
4609566063dSJacob Faibussowitsch     PetscCall(VecNestGetSubVecs(x, NULL, &vecs));
461c4762a1bSJed Brown     bv = vecs[0];
4629566063dSJacob Faibussowitsch     /*    PetscCall(VecBlockGetSubVec(x, 0, &bv)); */
4639566063dSJacob Faibussowitsch     PetscCall(VecSetValue(bv, 0, 2.0, INSERT_VALUES)); /* xx[0] = 2.0; */
4649566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(bv));
4659566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(bv));
466c4762a1bSJed Brown 
4679566063dSJacob Faibussowitsch     /*    PetscCall(VecBlockGetSubVec(x, 1, &bv)); */
468c4762a1bSJed Brown     bv = vecs[1];
4699566063dSJacob Faibussowitsch     PetscCall(VecSetValue(bv, 0, 3.0, INSERT_VALUES)); /* xx[1] = 3.0; */
4709566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(bv));
4719566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(bv));
472c4762a1bSJed Brown   }
473c4762a1bSJed Brown   /*
474c4762a1bSJed Brown   Note: The user should initialize the vector, x, with the initial guess
475c4762a1bSJed Brown   for the nonlinear solver prior to calling SNESSolve().  In particular,
476c4762a1bSJed Brown   to employ an initial guess of zero, the user should explicitly set
477c4762a1bSJed Brown   this vector to zero by calling VecSet().
478c4762a1bSJed Brown   */
4799566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
4809566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
481c4762a1bSJed Brown   if (flg) {
482c4762a1bSJed Brown     Vec f;
4839566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
4849566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &f, 0, 0));
4859566063dSJacob Faibussowitsch     PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
486c4762a1bSJed Brown   }
487c4762a1bSJed Brown 
48863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
489c4762a1bSJed Brown 
490c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
491c4762a1bSJed Brown   Free work space.  All PETSc objects should be destroyed when they
492c4762a1bSJed Brown   are no longer needed.
493c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4949371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
4959371c9d4SSatish Balay   PetscCall(VecDestroy(&r));
4969371c9d4SSatish Balay   PetscCall(MatDestroy(&J));
4979371c9d4SSatish Balay   PetscCall(SNESDestroy(&snes));
498c4762a1bSJed Brown   PetscFunctionReturn(0);
499c4762a1bSJed Brown }
500c4762a1bSJed Brown 
5019371c9d4SSatish Balay static PetscErrorCode FormFunction1_block(SNES snes, Vec x, Vec f, void *dummy) {
502c4762a1bSJed Brown   Vec        *xx, *ff, x1, x2, f1, f2;
503c4762a1bSJed Brown   PetscScalar ff_0, ff_1;
504c4762a1bSJed Brown   PetscScalar xx_0, xx_1;
505c4762a1bSJed Brown   PetscInt    index, nb;
506c4762a1bSJed Brown 
507c4762a1bSJed Brown   PetscFunctionBeginUser;
508c4762a1bSJed Brown   /* get blocks for function */
5099566063dSJacob Faibussowitsch   PetscCall(VecNestGetSubVecs(f, &nb, &ff));
5109371c9d4SSatish Balay   f1 = ff[0];
5119371c9d4SSatish Balay   f2 = ff[1];
512c4762a1bSJed Brown 
513c4762a1bSJed Brown   /* get blocks for solution */
5149566063dSJacob Faibussowitsch   PetscCall(VecNestGetSubVecs(x, &nb, &xx));
5159371c9d4SSatish Balay   x1 = xx[0];
5169371c9d4SSatish Balay   x2 = xx[1];
517c4762a1bSJed Brown 
518c4762a1bSJed Brown   /* get solution values */
519c4762a1bSJed Brown   index = 0;
5209566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x1, 1, &index, &xx_0));
5219566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x2, 1, &index, &xx_1));
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   /* Compute function */
524c4762a1bSJed Brown   ff_0 = xx_0 * xx_0 + xx_0 * xx_1 - 3.0;
525c4762a1bSJed Brown   ff_1 = xx_0 * xx_1 + xx_1 * xx_1 - 6.0;
526c4762a1bSJed Brown 
527c4762a1bSJed Brown   /* set function values */
5289566063dSJacob Faibussowitsch   PetscCall(VecSetValue(f1, index, ff_0, INSERT_VALUES));
529c4762a1bSJed Brown 
5309566063dSJacob Faibussowitsch   PetscCall(VecSetValue(f2, index, ff_1, INSERT_VALUES));
531c4762a1bSJed Brown 
5329566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(f));
5339566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(f));
534c4762a1bSJed Brown   PetscFunctionReturn(0);
535c4762a1bSJed Brown }
536c4762a1bSJed Brown 
5379371c9d4SSatish Balay static PetscErrorCode FormJacobian1_block(SNES snes, Vec x, Mat jac, Mat B, void *dummy) {
538c4762a1bSJed Brown   Vec        *xx, x1, x2;
539c4762a1bSJed Brown   PetscScalar xx_0, xx_1;
540c4762a1bSJed Brown   PetscInt    index, nb;
541c4762a1bSJed Brown   PetscScalar A_00, A_01, A_10, A_11;
542c4762a1bSJed Brown   Mat         j11, j12, j21, j22;
543c4762a1bSJed Brown   Mat       **mats;
544c4762a1bSJed Brown 
545c4762a1bSJed Brown   PetscFunctionBeginUser;
546c4762a1bSJed Brown   /* get blocks for solution */
5479566063dSJacob Faibussowitsch   PetscCall(VecNestGetSubVecs(x, &nb, &xx));
5489371c9d4SSatish Balay   x1 = xx[0];
5499371c9d4SSatish Balay   x2 = xx[1];
550c4762a1bSJed Brown 
551c4762a1bSJed Brown   /* get solution values */
552c4762a1bSJed Brown   index = 0;
5539566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x1, 1, &index, &xx_0));
5549566063dSJacob Faibussowitsch   PetscCall(VecGetValues(x2, 1, &index, &xx_1));
555c4762a1bSJed Brown 
556c4762a1bSJed Brown   /* get block matrices */
5579566063dSJacob Faibussowitsch   PetscCall(MatNestGetSubMats(jac, NULL, NULL, &mats));
558c4762a1bSJed Brown   j11 = mats[0][0];
559c4762a1bSJed Brown   j12 = mats[0][1];
560c4762a1bSJed Brown   j21 = mats[1][0];
561c4762a1bSJed Brown   j22 = mats[1][1];
562c4762a1bSJed Brown 
563c4762a1bSJed Brown   /* compute jacobian entries */
564c4762a1bSJed Brown   A_00 = 2.0 * xx_0 + xx_1;
565c4762a1bSJed Brown   A_01 = xx_0;
566c4762a1bSJed Brown   A_10 = xx_1;
567c4762a1bSJed Brown   A_11 = xx_0 + 2.0 * xx_1;
568c4762a1bSJed Brown 
569c4762a1bSJed Brown   /* set jacobian values */
5709566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j11, 0, 0, A_00, INSERT_VALUES));
5719566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j12, 0, 0, A_01, INSERT_VALUES));
5729566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j21, 0, 0, A_10, INSERT_VALUES));
5739566063dSJacob Faibussowitsch   PetscCall(MatSetValue(j22, 0, 0, A_11, INSERT_VALUES));
574c4762a1bSJed Brown 
575c4762a1bSJed Brown   /* Assemble sub matrix */
5769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
5779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
578c4762a1bSJed Brown   PetscFunctionReturn(0);
579c4762a1bSJed Brown }
580c4762a1bSJed Brown 
5819371c9d4SSatish Balay static PetscErrorCode FormFunction2_block(SNES snes, Vec x, Vec f, void *dummy) {
582c4762a1bSJed Brown   PetscScalar       *ff;
583c4762a1bSJed Brown   const PetscScalar *xx;
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   PetscFunctionBeginUser;
586c4762a1bSJed Brown   /*
587c4762a1bSJed Brown   Get pointers to vector data.
588c4762a1bSJed Brown   - For default PETSc vectors, VecGetArray() returns a pointer to
589c4762a1bSJed Brown   the data array.  Otherwise, the routine is implementation dependent.
590c4762a1bSJed Brown   - You MUST call VecRestoreArray() when you no longer need access to
591c4762a1bSJed Brown   the array.
592c4762a1bSJed Brown   */
5939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
5949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
595c4762a1bSJed Brown 
596c4762a1bSJed Brown   /*
597c4762a1bSJed Brown   Compute function
598c4762a1bSJed Brown   */
599c4762a1bSJed Brown   ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0];
600c4762a1bSJed Brown   ff[1] = xx[1];
601c4762a1bSJed Brown 
602c4762a1bSJed Brown   /*
603c4762a1bSJed Brown   Restore vectors
604c4762a1bSJed Brown   */
6059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
6069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
607c4762a1bSJed Brown   PetscFunctionReturn(0);
608c4762a1bSJed Brown }
609c4762a1bSJed Brown 
6109371c9d4SSatish Balay static PetscErrorCode FormJacobian2_block(SNES snes, Vec x, Mat jac, Mat B, void *dummy) {
611c4762a1bSJed Brown   const PetscScalar *xx;
612c4762a1bSJed Brown   PetscScalar        A[4];
613c4762a1bSJed Brown   PetscInt           idx[2] = {0, 1};
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   PetscFunctionBeginUser;
616c4762a1bSJed Brown   /*
617c4762a1bSJed Brown   Get pointer to vector data
618c4762a1bSJed Brown   */
6199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
620c4762a1bSJed Brown 
621c4762a1bSJed Brown   /*
622c4762a1bSJed Brown   Compute Jacobian entries and insert into matrix.
623c4762a1bSJed Brown   - Since this is such a small problem, we set all entries for
624c4762a1bSJed Brown   the matrix at once.
625c4762a1bSJed Brown   */
6269371c9d4SSatish Balay   A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0;
6279371c9d4SSatish Balay   A[1] = 0.0;
6289371c9d4SSatish Balay   A[2] = 0.0;
6299371c9d4SSatish Balay   A[3] = 1.0;
6309566063dSJacob Faibussowitsch   PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES));
631c4762a1bSJed Brown 
632c4762a1bSJed Brown   /*
633c4762a1bSJed Brown   Restore vector
634c4762a1bSJed Brown   */
6359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   /*
638c4762a1bSJed Brown   Assemble matrix
639c4762a1bSJed Brown   */
6409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
6419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
642c4762a1bSJed Brown   PetscFunctionReturn(0);
643c4762a1bSJed Brown }
644c4762a1bSJed Brown 
6459371c9d4SSatish Balay int main(int argc, char **argv) {
646c4762a1bSJed Brown   PetscMPIInt size;
647c4762a1bSJed Brown 
648327415f7SBarry Smith   PetscFunctionBeginUser;
6499566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
6509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
651be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
6529566063dSJacob Faibussowitsch   PetscCall(assembled_system());
6539566063dSJacob Faibussowitsch   PetscCall(block_system());
6549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
655b122ec5aSJacob Faibussowitsch   return 0;
656c4762a1bSJed Brown }
657c4762a1bSJed Brown 
658c4762a1bSJed Brown /*TEST
659c4762a1bSJed Brown 
660c4762a1bSJed Brown    test:
661c4762a1bSJed Brown       args: -snes_monitor_short
662c4762a1bSJed Brown       requires: !single
663c4762a1bSJed Brown 
664c4762a1bSJed Brown TEST*/
665