xref: /petsc/src/snes/tests/ex17.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 static const char help[] = "Newton's method to solve a two-variable system, sequentially.\n"
2                            "The same problem is solved twice - i) fully assembled system + ii) block system\n\n";
3 
4 /*T
5 Concepts: SNES^basic uniprocessor example, block objects
6 Processors: 1
7 T*/
8 
9 
10 
11 /*
12 Include "petscsnes.h" so that we can use SNES solvers.  Note that this
13 file automatically includes:
14 petscsys.h       - base PETSc routines   petscvec.h - vectors
15 petscsys.h    - system routines       petscmat.h - matrices
16 petscis.h     - index sets            petscksp.h - Krylov subspace methods
17 petscviewer.h - viewers               petscpc.h  - preconditioners
18 petscksp.h   - linear solvers
19 */
20 #include <petscsnes.h>
21 
22 /*
23 This example is block version of the test found at
24   ${PETSC_DIR}/src/snes/tutorials/ex1.c
25 In this test we replace the Jacobian systems
26   [J]{x} = {F}
27 where
28 
29 [J] = (j_00, j_01),  {x} = (x_0, x_1)^T,   {F} = (f_0, f_1)^T
30       (j_10, j_11)
31 where [J] \in \mathbb^{2 \times 2}, {x},{F} \in \mathbb^{2 \times 1},
32 
33 with a block system in which each block is of length 1.
34 i.e. The block system is thus
35 
36 [J] = ([j00], [j01]),  {x} = ({x0}, {x1})^T, {F} = ({f0}, {f1})^T
37       ([j10], [j11])
38 where
39 [j00], [j01], [j10], [j11] \in \mathbb^{1 \times 1}
40 {x0}, {x1}, {f0}, {f1} \in \mathbb^{1 \times 1}
41 
42 In practice we would not bother defing blocks of size one, and would instead assemble the
43 full system. This is just a simple test to illustrate how to manipulate the blocks and
44 to confirm the implementation is correct.
45 */
46 
47 /*
48 User-defined routines
49 */
50 static PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*);
51 static PetscErrorCode FormFunction1(SNES,Vec,Vec,void*);
52 static PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*);
53 static PetscErrorCode FormFunction2(SNES,Vec,Vec,void*);
54 static PetscErrorCode FormJacobian1_block(SNES,Vec,Mat,Mat,void*);
55 static PetscErrorCode FormFunction1_block(SNES,Vec,Vec,void*);
56 static PetscErrorCode FormJacobian2_block(SNES,Vec,Mat,Mat,void*);
57 static PetscErrorCode FormFunction2_block(SNES,Vec,Vec,void*);
58 
59 
60 static PetscErrorCode assembled_system(void)
61 {
62   SNES           snes;         /* nonlinear solver context */
63   KSP            ksp;         /* linear solver context */
64   PC             pc;           /* preconditioner context */
65   Vec            x,r;         /* solution, residual vectors */
66   Mat            J;            /* Jacobian matrix */
67   PetscErrorCode ierr;
68   PetscInt       its;
69   PetscScalar    pfive = .5,*xx;
70   PetscBool      flg;
71 
72   PetscFunctionBeginUser;
73   PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n");
74 
75   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76   Create nonlinear solver context
77   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78 
79   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
80 
81   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82   Create matrix and vector data structures; set corresponding routines
83   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84 
85   /*
86   Create vectors for solution and nonlinear function
87   */
88   ierr = VecCreateSeq(PETSC_COMM_SELF,2,&x);CHKERRQ(ierr);
89   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
90 
91   /*
92   Create Jacobian matrix data structure
93   */
94   ierr = MatCreate(PETSC_COMM_SELF,&J);CHKERRQ(ierr);
95   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
96   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
97   ierr = MatSetUp(J);CHKERRQ(ierr);
98 
99   ierr = PetscOptionsHasName(NULL,NULL,"-hard",&flg);CHKERRQ(ierr);
100   if (!flg) {
101     /*
102     Set function evaluation routine and vector.
103     */
104     ierr = SNESSetFunction(snes,r,FormFunction1,NULL);CHKERRQ(ierr);
105 
106     /*
107     Set Jacobian matrix data structure and Jacobian evaluation routine
108     */
109     ierr = SNESSetJacobian(snes,J,J,FormJacobian1,NULL);CHKERRQ(ierr);
110   } else {
111     ierr = SNESSetFunction(snes,r,FormFunction2,NULL);CHKERRQ(ierr);
112     ierr = SNESSetJacobian(snes,J,J,FormJacobian2,NULL);CHKERRQ(ierr);
113   }
114 
115   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116   Customize nonlinear solver; set runtime options
117   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118 
119   /*
120   Set linear solver defaults for this problem. By extracting the
121   KSP, KSP, and PC contexts from the SNES context, we can then
122   directly call any KSP, KSP, and PC routines to set various options.
123   */
124   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
125   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
126   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
127   ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr);
128 
129   /*
130   Set SNES/KSP/KSP/PC runtime options, e.g.,
131   -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
132   These options will override those specified above as long as
133   SNESSetFromOptions() is called _after_ any other customization
134   routines.
135   */
136   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
137 
138   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139   Evaluate initial guess; then solve nonlinear system
140   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141   if (!flg) {
142     ierr = VecSet(x,pfive);CHKERRQ(ierr);
143   } else {
144     ierr  = VecGetArray(x,&xx);CHKERRQ(ierr);
145     xx[0] = 2.0; xx[1] = 3.0;
146     ierr  = VecRestoreArray(x,&xx);CHKERRQ(ierr);
147   }
148   /*
149   Note: The user should initialize the vector, x, with the initial guess
150   for the nonlinear solver prior to calling SNESSolve().  In particular,
151   to employ an initial guess of zero, the user should explicitly set
152   this vector to zero by calling VecSet().
153   */
154 
155   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
156   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
157   if (flg) {
158     Vec f;
159     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
160     ierr = SNESGetFunction(snes,&f,0,0);CHKERRQ(ierr);
161     ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
162   }
163 
164   ierr = PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr);
165 
166   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167   Free work space.  All PETSc objects should be destroyed when they
168   are no longer needed.
169   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170 
171   ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr);
172   ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 /* ------------------------------------------------------------------- */
177 /*
178 FormFunction1 - Evaluates nonlinear function, F(x).
179 
180 Input Parameters:
181 .  snes - the SNES context
182 .  x - input vector
183 .  dummy - optional user-defined context (not used here)
184 
185 Output Parameter:
186 .  f - function vector
187 */
188 static PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *dummy)
189 {
190   PetscErrorCode    ierr;
191   const PetscScalar *xx;
192   PetscScalar       *ff;
193 
194   PetscFunctionBeginUser;
195   /*
196   Get pointers to vector data.
197   - For default PETSc vectors, VecGetArray() returns a pointer to
198   the data array.  Otherwise, the routine is implementation dependent.
199   - You MUST call VecRestoreArray() when you no longer need access to
200   the array.
201   */
202   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
203   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
204 
205   /*
206   Compute function
207   */
208   ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0;
209   ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0;
210 
211 
212   /*
213   Restore vectors
214   */
215   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
216   ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
217   PetscFunctionReturn(0);
218 }
219 /* ------------------------------------------------------------------- */
220 /*
221 FormJacobian1 - Evaluates Jacobian matrix.
222 
223 Input Parameters:
224 .  snes - the SNES context
225 .  x - input vector
226 .  dummy - optional user-defined context (not used here)
227 
228 Output Parameters:
229 .  jac - Jacobian matrix
230 .  B - optionally different preconditioning matrix
231 .  flag - flag indicating matrix structure
232 */
233 static PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
234 {
235   const PetscScalar *xx;
236   PetscScalar       A[4];
237   PetscErrorCode    ierr;
238   PetscInt          idx[2] = {0,1};
239 
240   PetscFunctionBeginUser;
241   /*
242   Get pointer to vector data
243   */
244   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
245 
246   /*
247   Compute Jacobian entries and insert into matrix.
248   - Since this is such a small problem, we set all entries for
249   the matrix at once.
250   */
251   A[0]  = 2.0*xx[0] + xx[1]; A[1] = xx[0];
252   A[2]  = xx[1]; A[3] = xx[0] + 2.0*xx[1];
253   ierr  = MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr);
254 
255   /*
256   Restore vector
257   */
258   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
259 
260   /*
261   Assemble matrix
262   */
263   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 
268 
269 /* ------------------------------------------------------------------- */
270 static PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy)
271 {
272   PetscErrorCode    ierr;
273   const PetscScalar *xx;
274   PetscScalar       *ff;
275 
276   PetscFunctionBeginUser;
277   /*
278   Get pointers to vector data.
279   - For default PETSc vectors, VecGetArray() returns a pointer to
280   the data array.  Otherwise, the routine is implementation dependent.
281   - You MUST call VecRestoreArray() when you no longer need access to
282   the array.
283   */
284   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
285   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
286 
287   /*
288   Compute function
289   */
290   ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0];
291   ff[1] = xx[1];
292 
293   /*
294   Restore vectors
295   */
296   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
297   ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
298   PetscFunctionReturn(0);
299 }
300 
301 /* ------------------------------------------------------------------- */
302 static PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
303 {
304   const PetscScalar *xx;
305   PetscScalar       A[4];
306   PetscErrorCode    ierr;
307   PetscInt          idx[2] = {0,1};
308 
309   PetscFunctionBeginUser;
310   /*
311   Get pointer to vector data
312   */
313   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
314 
315   /*
316   Compute Jacobian entries and insert into matrix.
317   - Since this is such a small problem, we set all entries for
318   the matrix at once.
319   */
320   A[0]  = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0;
321   A[2]  = 0.0;                     A[3] = 1.0;
322   ierr  = MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr);
323 
324   /*
325   Restore vector
326   */
327   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
328 
329   /*
330   Assemble matrix
331   */
332   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
333   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
334   PetscFunctionReturn(0);
335 }
336 
337 static int block_system(void)
338 {
339   SNES           snes;         /* nonlinear solver context */
340   KSP            ksp;         /* linear solver context */
341   PC             pc;           /* preconditioner context */
342   Vec            x,r;         /* solution, residual vectors */
343   Mat            J;            /* Jacobian matrix */
344   PetscErrorCode ierr;
345   PetscInt       its;
346   PetscScalar    pfive = .5;
347   PetscBool      flg;
348 
349   Mat j11, j12, j21, j22;
350   Vec x1, x2, r1, r2;
351   Vec bv;
352   Vec bx[2];
353   Mat bA[2][2];
354 
355   PetscFunctionBeginUser;
356   PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n");
357 
358   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
359 
360   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
361   Create matrix and vector data structures; set corresponding routines
362   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
363 
364   /*
365   Create sub vectors for solution and nonlinear function
366   */
367   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&x1);CHKERRQ(ierr);
368   ierr = VecDuplicate(x1,&r1);CHKERRQ(ierr);
369 
370   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&x2);CHKERRQ(ierr);
371   ierr = VecDuplicate(x2,&r2);CHKERRQ(ierr);
372 
373   /*
374   Create the block vectors
375   */
376   bx[0] = x1;
377   bx[1] = x2;
378   ierr  = VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&x);CHKERRQ(ierr);
379   ierr  = VecAssemblyBegin(x);CHKERRQ(ierr);
380   ierr  = VecAssemblyEnd(x);CHKERRQ(ierr);
381   ierr  = VecDestroy(&x1);CHKERRQ(ierr);
382   ierr  = VecDestroy(&x2);CHKERRQ(ierr);
383 
384   bx[0] = r1;
385   bx[1] = r2;
386   ierr  = VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&r);CHKERRQ(ierr);
387   ierr  = VecDestroy(&r1);CHKERRQ(ierr);
388   ierr  = VecDestroy(&r2);CHKERRQ(ierr);
389   ierr  = VecAssemblyBegin(r);CHKERRQ(ierr);
390   ierr  = VecAssemblyEnd(r);CHKERRQ(ierr);
391 
392   /*
393   Create sub Jacobian matrix data structure
394   */
395   ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr);
396   ierr = MatSetSizes(j11, 1, 1, 1, 1);CHKERRQ(ierr);
397   ierr = MatSetType(j11, MATSEQAIJ);CHKERRQ(ierr);
398   ierr = MatSetUp(j11);CHKERRQ(ierr);
399 
400   ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr);
401   ierr = MatSetSizes(j12, 1, 1, 1, 1);CHKERRQ(ierr);
402   ierr = MatSetType(j12, MATSEQAIJ);CHKERRQ(ierr);
403   ierr = MatSetUp(j12);CHKERRQ(ierr);
404 
405   ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr);
406   ierr = MatSetSizes(j21, 1, 1, 1, 1);CHKERRQ(ierr);
407   ierr = MatSetType(j21, MATSEQAIJ);CHKERRQ(ierr);
408   ierr = MatSetUp(j21);CHKERRQ(ierr);
409 
410   ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr);
411   ierr = MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1);CHKERRQ(ierr);
412   ierr = MatSetType(j22, MATSEQAIJ);CHKERRQ(ierr);
413   ierr = MatSetUp(j22);CHKERRQ(ierr);
414   /*
415   Create block Jacobian matrix data structure
416   */
417   bA[0][0] = j11;
418   bA[0][1] = j12;
419   bA[1][0] = j21;
420   bA[1][1] = j22;
421 
422   ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],&J);CHKERRQ(ierr);
423   ierr = MatSetUp(J);CHKERRQ(ierr);
424   ierr = MatNestSetVecType(J,VECNEST);CHKERRQ(ierr);
425   ierr = MatDestroy(&j11);CHKERRQ(ierr);
426   ierr = MatDestroy(&j12);CHKERRQ(ierr);
427   ierr = MatDestroy(&j21);CHKERRQ(ierr);
428   ierr = MatDestroy(&j22);CHKERRQ(ierr);
429 
430   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
431   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
432 
433   ierr = PetscOptionsHasName(NULL,NULL,"-hard",&flg);CHKERRQ(ierr);
434   if (!flg) {
435     /*
436     Set function evaluation routine and vector.
437     */
438     ierr = SNESSetFunction(snes,r,FormFunction1_block,NULL);CHKERRQ(ierr);
439 
440     /*
441     Set Jacobian matrix data structure and Jacobian evaluation routine
442     */
443     ierr = SNESSetJacobian(snes,J,J,FormJacobian1_block,NULL);CHKERRQ(ierr);
444   } else {
445     ierr = SNESSetFunction(snes,r,FormFunction2_block,NULL);CHKERRQ(ierr);
446     ierr = SNESSetJacobian(snes,J,J,FormJacobian2_block,NULL);CHKERRQ(ierr);
447   }
448 
449   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
450   Customize nonlinear solver; set runtime options
451   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
452 
453   /*
454   Set linear solver defaults for this problem. By extracting the
455   KSP, KSP, and PC contexts from the SNES context, we can then
456   directly call any KSP, KSP, and PC routines to set various options.
457   */
458   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
459   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
460   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
461   ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr);
462 
463   /*
464   Set SNES/KSP/KSP/PC runtime options, e.g.,
465   -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
466   These options will override those specified above as long as
467   SNESSetFromOptions() is called _after_ any other customization
468   routines.
469   */
470   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
471 
472   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
473   Evaluate initial guess; then solve nonlinear system
474   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
475   if (!flg) {
476     ierr = VecSet(x,pfive);CHKERRQ(ierr);
477   } else {
478     Vec *vecs;
479     ierr = VecNestGetSubVecs(x, NULL, &vecs);CHKERRQ(ierr);
480     bv   = vecs[0];
481 /*    ierr = VecBlockGetSubVec(x, 0, &bv);CHKERRQ(ierr); */
482     ierr = VecSetValue(bv, 0, 2.0, INSERT_VALUES);CHKERRQ(ierr);  /* xx[0] = 2.0; */
483     ierr = VecAssemblyBegin(bv);CHKERRQ(ierr);
484     ierr = VecAssemblyEnd(bv);CHKERRQ(ierr);
485 
486 /*    ierr = VecBlockGetSubVec(x, 1, &bv);CHKERRQ(ierr); */
487     bv   = vecs[1];
488     ierr = VecSetValue(bv, 0, 3.0, INSERT_VALUES);CHKERRQ(ierr);  /* xx[1] = 3.0; */
489     ierr = VecAssemblyBegin(bv);CHKERRQ(ierr);
490     ierr = VecAssemblyEnd(bv);CHKERRQ(ierr);
491   }
492   /*
493   Note: The user should initialize the vector, x, with the initial guess
494   for the nonlinear solver prior to calling SNESSolve().  In particular,
495   to employ an initial guess of zero, the user should explicitly set
496   this vector to zero by calling VecSet().
497   */
498   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
499   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
500   if (flg) {
501     Vec f;
502     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
503     ierr = SNESGetFunction(snes,&f,0,0);CHKERRQ(ierr);
504     ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
505   }
506 
507   ierr = PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr);
508 
509   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
510   Free work space.  All PETSc objects should be destroyed when they
511   are no longer needed.
512   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
513   ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr);
514   ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr);
515   PetscFunctionReturn(0);
516 }
517 
518 /* ------------------------------------------------------------------- */
519 static PetscErrorCode FormFunction1_block(SNES snes,Vec x,Vec f,void *dummy)
520 {
521   Vec            *xx, *ff, x1,x2, f1,f2;
522   PetscScalar    ff_0, ff_1;
523   PetscScalar    xx_0, xx_1;
524   PetscInt       index,nb;
525   PetscErrorCode ierr;
526 
527   PetscFunctionBeginUser;
528   /* get blocks for function */
529   ierr = VecNestGetSubVecs(f, &nb, &ff);CHKERRQ(ierr);
530   f1   = ff[0];  f2 = ff[1];
531 
532   /* get blocks for solution */
533   ierr = VecNestGetSubVecs(x, &nb, &xx);CHKERRQ(ierr);
534   x1   = xx[0];  x2 = xx[1];
535 
536   /* get solution values */
537   index = 0;
538   ierr  = VecGetValues(x1,1, &index, &xx_0);CHKERRQ(ierr);
539   ierr  = VecGetValues(x2,1, &index, &xx_1);CHKERRQ(ierr);
540 
541   /* Compute function */
542   ff_0 = xx_0*xx_0 + xx_0*xx_1 - 3.0;
543   ff_1 = xx_0*xx_1 + xx_1*xx_1 - 6.0;
544 
545   /* set function values */
546   ierr = VecSetValue(f1, index, ff_0, INSERT_VALUES);CHKERRQ(ierr);
547 
548   ierr = VecSetValue(f2, index, ff_1, INSERT_VALUES);CHKERRQ(ierr);
549 
550   ierr = VecAssemblyBegin(f);CHKERRQ(ierr);
551   ierr = VecAssemblyEnd(f);CHKERRQ(ierr);
552   PetscFunctionReturn(0);
553 }
554 
555 /* ------------------------------------------------------------------- */
556 static PetscErrorCode FormJacobian1_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
557 {
558   Vec            *xx, x1,x2;
559   PetscScalar    xx_0, xx_1;
560   PetscInt       index,nb;
561   PetscScalar    A_00, A_01, A_10, A_11;
562   Mat            j11, j12, j21, j22;
563   Mat            **mats;
564   PetscErrorCode ierr;
565 
566   PetscFunctionBeginUser;
567   /* get blocks for solution */
568   ierr = VecNestGetSubVecs(x, &nb, &xx);CHKERRQ(ierr);
569   x1   = xx[0];  x2 = xx[1];
570 
571   /* get solution values */
572   index = 0;
573   ierr  = VecGetValues(x1,1, &index, &xx_0);CHKERRQ(ierr);
574   ierr  = VecGetValues(x2,1, &index, &xx_1);CHKERRQ(ierr);
575 
576   /* get block matrices */
577   ierr = MatNestGetSubMats(jac,NULL,NULL,&mats);CHKERRQ(ierr);
578   j11  = mats[0][0];
579   j12  = mats[0][1];
580   j21  = mats[1][0];
581   j22  = mats[1][1];
582 
583   /* compute jacobian entries */
584   A_00 = 2.0*xx_0 + xx_1;
585   A_01 = xx_0;
586   A_10 = xx_1;
587   A_11 = xx_0 + 2.0*xx_1;
588 
589   /* set jacobian values */
590   ierr = MatSetValue(j11, 0,0, A_00, INSERT_VALUES);CHKERRQ(ierr);
591   ierr = MatSetValue(j12, 0,0, A_01, INSERT_VALUES);CHKERRQ(ierr);
592   ierr = MatSetValue(j21, 0,0, A_10, INSERT_VALUES);CHKERRQ(ierr);
593   ierr = MatSetValue(j22, 0,0, A_11, INSERT_VALUES);CHKERRQ(ierr);
594 
595   /* Assemble sub matrix */
596   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
597   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
598   PetscFunctionReturn(0);
599 }
600 
601 /* ------------------------------------------------------------------- */
602 static PetscErrorCode FormFunction2_block(SNES snes,Vec x,Vec f,void *dummy)
603 {
604   PetscErrorCode    ierr;
605   PetscScalar       *ff;
606   const PetscScalar *xx;
607 
608   PetscFunctionBeginUser;
609   /*
610   Get pointers to vector data.
611   - For default PETSc vectors, VecGetArray() returns a pointer to
612   the data array.  Otherwise, the routine is implementation dependent.
613   - You MUST call VecRestoreArray() when you no longer need access to
614   the array.
615   */
616   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
617   ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
618 
619   /*
620   Compute function
621   */
622   ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0];
623   ff[1] = xx[1];
624 
625   /*
626   Restore vectors
627   */
628   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
629   ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
630   PetscFunctionReturn(0);
631 }
632 
633 /* ------------------------------------------------------------------- */
634 static PetscErrorCode FormJacobian2_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
635 {
636   const PetscScalar *xx;
637   PetscScalar       A[4];
638   PetscErrorCode    ierr;
639   PetscInt          idx[2] = {0,1};
640 
641   PetscFunctionBeginUser;
642   /*
643   Get pointer to vector data
644   */
645   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
646 
647   /*
648   Compute Jacobian entries and insert into matrix.
649   - Since this is such a small problem, we set all entries for
650   the matrix at once.
651   */
652   A[0]  = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0;
653   A[2]  = 0.0;                     A[3] = 1.0;
654   ierr  = MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr);
655 
656   /*
657   Restore vector
658   */
659   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
660 
661   /*
662   Assemble matrix
663   */
664   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
665   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
666   PetscFunctionReturn(0);
667 }
668 
669 int main(int argc,char **argv)
670 {
671   PetscMPIInt    size;
672   PetscErrorCode ierr;
673 
674   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
675 
676   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
677   if (size != 1) SETERRQ(PETSC_COMM_WORLD, 1,"This is a uniprocessor example only!");
678 
679   ierr = assembled_system();CHKERRQ(ierr);
680 
681   ierr = block_system();CHKERRQ(ierr);
682 
683   ierr = PetscFinalize();
684   return ierr;
685 }
686 
687 
688 /*TEST
689 
690    test:
691       args: -snes_monitor_short
692       requires: !single
693 
694 TEST*/
695