xref: /petsc/src/snes/impls/multiblock/multiblock.c (revision 2c9ec6dfe874b911fa49ef7e759d29a8430d6aff)
1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
2 #include <petscdmcomposite.h>
3 
4 typedef struct _BlockDesc *BlockDesc;
5 struct _BlockDesc {
6   char       *name;     /* Block name */
7   PetscInt   nfields;   /* If block is defined on a DA, the number of DA fields */
8   PetscInt   *fields;   /* If block is defined on a DA, the list of DA fields */
9   IS         is;        /* Index sets defining the block */
10   VecScatter sctx;      /* Scatter mapping global Vec to blockVec */
11   SNES       snes;      /* Solver for this block */
12   Vec        x;
13   BlockDesc  next, previous;
14 };
15 
16 typedef struct {
17   PetscBool       issetup;       /* Flag is true after the all ISs and operators have been defined */
18   PetscBool       defined;       /* Flag is true after the blocks have been defined, to prevent more blocks from being added */
19   PetscBool       defaultblocks; /* Flag is true for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
20   PetscInt        numBlocks;     /* Number of blocks (can be fields, domains, etc.) */
21   PetscInt        bs;            /* Block size for IS, Vec and Mat structures */
22   PCCompositeType type;          /* Solver combination method (additive, multiplicative, etc.) */
23   BlockDesc       blocks;        /* Linked list of block descriptors */
24 } SNES_Multiblock;
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "SNESReset_Multiblock"
28 PetscErrorCode SNESReset_Multiblock(SNES snes)
29 {
30   SNES_Multiblock *mb    = (SNES_Multiblock*) snes->data;
31   BlockDesc       blocks = mb->blocks, next;
32   PetscErrorCode  ierr;
33 
34   PetscFunctionBegin;
35   while (blocks) {
36     ierr = SNESReset(blocks->snes);CHKERRQ(ierr);
37 #if 0
38     ierr = VecDestroy(&blocks->x);CHKERRQ(ierr);
39 #endif
40     ierr   = VecScatterDestroy(&blocks->sctx);CHKERRQ(ierr);
41     ierr   = ISDestroy(&blocks->is);CHKERRQ(ierr);
42     next   = blocks->next;
43     blocks = next;
44   }
45   PetscFunctionReturn(0);
46 }
47 
48 /*
49   SNESDestroy_Multiblock - Destroys the private SNES_Multiblock context that was created with SNESCreate_Multiblock().
50 
51   Input Parameter:
52 . snes - the SNES context
53 
54   Application Interface Routine: SNESDestroy()
55 */
56 #undef __FUNCT__
57 #define __FUNCT__ "SNESDestroy_Multiblock"
58 PetscErrorCode SNESDestroy_Multiblock(SNES snes)
59 {
60   SNES_Multiblock *mb    = (SNES_Multiblock*) snes->data;
61   BlockDesc       blocks = mb->blocks, next;
62   PetscErrorCode  ierr;
63 
64   PetscFunctionBegin;
65   ierr = SNESReset_Multiblock(snes);CHKERRQ(ierr);
66   while (blocks) {
67     next   = blocks->next;
68     ierr   = SNESDestroy(&blocks->snes);CHKERRQ(ierr);
69     ierr   = PetscFree(blocks->name);CHKERRQ(ierr);
70     ierr   = PetscFree(blocks->fields);CHKERRQ(ierr);
71     ierr   = PetscFree(blocks);CHKERRQ(ierr);
72     blocks = next;
73   }
74   ierr = PetscFree(snes->data);CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 
78 #undef __FUNCT__
79 #define __FUNCT__ "SNESMultiblockSetFieldsRuntime_Private"
80 /* Precondition: blocksize is set to a meaningful value */
81 static PetscErrorCode SNESMultiblockSetFieldsRuntime_Private(SNES snes)
82 {
83   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
84   PetscInt        *ifields;
85   PetscInt        i, nfields;
86   PetscBool       flg = PETSC_TRUE;
87   char            optionname[128], name[8];
88   PetscErrorCode  ierr;
89 
90   PetscFunctionBegin;
91   ierr = PetscMalloc1(mb->bs, &ifields);CHKERRQ(ierr);
92   for (i = 0;; ++i) {
93     ierr    = PetscSNPrintf(name, sizeof(name), "%D", i);CHKERRQ(ierr);
94     ierr    = PetscSNPrintf(optionname, sizeof(optionname), "-snes_multiblock_%D_fields", i);CHKERRQ(ierr);
95     nfields = mb->bs;
96     ierr    = PetscOptionsGetIntArray(((PetscObject) snes)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
97     if (!flg) break;
98     if (!nfields) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Cannot list zero fields");
99     ierr = SNESMultiblockSetFields(snes, name, nfields, ifields);CHKERRQ(ierr);
100   }
101   if (i > 0) {
102     /* Makes command-line setting of blocks take precedence over setting them in code.
103        Otherwise subsequent calls to SNESMultiblockSetIS() or SNESMultiblockSetFields() would
104        create new blocks, which would probably not be what the user wanted. */
105     mb->defined = PETSC_TRUE;
106   }
107   ierr = PetscFree(ifields);CHKERRQ(ierr);
108   PetscFunctionReturn(0);
109 }
110 
111 #undef __FUNCT__
112 #define __FUNCT__ "SNESMultiblockSetDefaults"
113 static PetscErrorCode SNESMultiblockSetDefaults(SNES snes)
114 {
115   SNES_Multiblock *mb    = (SNES_Multiblock*) snes->data;
116   BlockDesc       blocks = mb->blocks;
117   PetscInt        i;
118   PetscErrorCode  ierr;
119 
120   PetscFunctionBegin;
121   if (!blocks) {
122     if (snes->dm) {
123       PetscBool dmcomposite;
124 
125       ierr = PetscObjectTypeCompare((PetscObject) snes->dm, DMCOMPOSITE, &dmcomposite);CHKERRQ(ierr);
126       if (dmcomposite) {
127         PetscInt nDM;
128         IS       *fields;
129 
130         ierr = PetscInfo(snes,"Setting up physics based multiblock solver using the embedded DM\n");CHKERRQ(ierr);
131         ierr = DMCompositeGetNumberDM(snes->dm, &nDM);CHKERRQ(ierr);
132         ierr = DMCompositeGetGlobalISs(snes->dm, &fields);CHKERRQ(ierr);
133         for (i = 0; i < nDM; ++i) {
134           char name[8];
135 
136           ierr = PetscSNPrintf(name, sizeof(name), "%D", i);CHKERRQ(ierr);
137           ierr = SNESMultiblockSetIS(snes, name, fields[i]);CHKERRQ(ierr);
138           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
139         }
140         ierr = PetscFree(fields);CHKERRQ(ierr);
141       }
142     } else {
143       PetscBool flg    = PETSC_FALSE;
144       PetscBool stokes = PETSC_FALSE;
145 
146       if (mb->bs <= 0) {
147         if (snes->jacobian_pre) {
148           ierr = MatGetBlockSize(snes->jacobian_pre, &mb->bs);CHKERRQ(ierr);
149         } else mb->bs = 1;
150       }
151 
152       ierr = PetscOptionsGetBool(((PetscObject) snes)->prefix, "-snes_multiblock_default", &flg, NULL);CHKERRQ(ierr);
153       ierr = PetscOptionsGetBool(((PetscObject) snes)->prefix, "-snes_multiblock_detect_saddle_point", &stokes, NULL);CHKERRQ(ierr);
154       if (stokes) {
155         IS       zerodiags, rest;
156         PetscInt nmin, nmax;
157 
158         ierr = MatGetOwnershipRange(snes->jacobian_pre, &nmin, &nmax);CHKERRQ(ierr);
159         ierr = MatFindZeroDiagonals(snes->jacobian_pre, &zerodiags);CHKERRQ(ierr);
160         ierr = ISComplement(zerodiags, nmin, nmax, &rest);CHKERRQ(ierr);
161         ierr = SNESMultiblockSetIS(snes, "0", rest);CHKERRQ(ierr);
162         ierr = SNESMultiblockSetIS(snes, "1", zerodiags);CHKERRQ(ierr);
163         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
164         ierr = ISDestroy(&rest);CHKERRQ(ierr);
165       } else {
166         if (!flg) {
167           /* Allow user to set fields from command line, if bs was known at the time of SNESSetFromOptions_Multiblock()
168            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
169           ierr = SNESMultiblockSetFieldsRuntime_Private(snes);CHKERRQ(ierr);
170           if (mb->defined) {ierr = PetscInfo(snes, "Blocks defined using the options database\n");CHKERRQ(ierr);}
171         }
172         if (flg || !mb->defined) {
173           ierr = PetscInfo(snes, "Using default splitting of fields\n");CHKERRQ(ierr);
174           for (i = 0; i < mb->bs; ++i) {
175             char name[8];
176 
177             ierr = PetscSNPrintf(name, sizeof(name), "%D", i);CHKERRQ(ierr);
178             ierr = SNESMultiblockSetFields(snes, name, 1, &i);CHKERRQ(ierr);
179           }
180           mb->defaultblocks = PETSC_TRUE;
181         }
182       }
183     }
184   } else if (mb->numBlocks == 1) {
185     if (blocks->is) {
186       IS       is2;
187       PetscInt nmin, nmax;
188 
189       ierr = MatGetOwnershipRange(snes->jacobian_pre, &nmin, &nmax);CHKERRQ(ierr);
190       ierr = ISComplement(blocks->is, nmin, nmax, &is2);CHKERRQ(ierr);
191       ierr = SNESMultiblockSetIS(snes, "1", is2);CHKERRQ(ierr);
192       ierr = ISDestroy(&is2);CHKERRQ(ierr);
193     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least two sets of fields to SNES multiblock");
194   }
195   if (mb->numBlocks < 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unhandled case, must have at least two blocks");
196   PetscFunctionReturn(0);
197 }
198 
199 /*
200    SNESSetUp_Multiblock - Sets up the internal data structures for the later use
201    of the SNESMULTIBLOCK nonlinear solver.
202 
203    Input Parameters:
204 +  snes - the SNES context
205 -  x - the solution vector
206 
207    Application Interface Routine: SNESSetUp()
208 */
209 #undef __FUNCT__
210 #define __FUNCT__ "SNESSetUp_Multiblock"
211 PetscErrorCode SNESSetUp_Multiblock(SNES snes)
212 {
213   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
214   BlockDesc       blocks;
215   PetscInt        i, numBlocks;
216   PetscErrorCode  ierr;
217 
218   PetscFunctionBegin;
219   ierr      = SNESMultiblockSetDefaults(snes);CHKERRQ(ierr);
220   numBlocks = mb->numBlocks;
221   blocks    = mb->blocks;
222 
223   /* Create ISs */
224   if (!mb->issetup) {
225     PetscInt  ccsize, rstart, rend, nslots, bs;
226     PetscBool sorted;
227 
228     mb->issetup = PETSC_TRUE;
229     bs          = mb->bs;
230     ierr        = MatGetOwnershipRange(snes->jacobian_pre, &rstart, &rend);CHKERRQ(ierr);
231     ierr        = MatGetLocalSize(snes->jacobian_pre, NULL, &ccsize);CHKERRQ(ierr);
232     nslots      = (rend - rstart)/bs;
233     for (i = 0; i < numBlocks; ++i) {
234       if (mb->defaultblocks) {
235         ierr = ISCreateStride(PetscObjectComm((PetscObject)snes), nslots, rstart+i, numBlocks, &blocks->is);CHKERRQ(ierr);
236       } else if (!blocks->is) {
237         if (blocks->nfields > 1) {
238           PetscInt *ii, j, k, nfields = blocks->nfields, *fields = blocks->fields;
239 
240           ierr = PetscMalloc1(nfields*nslots, &ii);CHKERRQ(ierr);
241           for (j = 0; j < nslots; ++j) {
242             for (k = 0; k < nfields; ++k) {
243               ii[nfields*j + k] = rstart + bs*j + fields[k];
244             }
245           }
246           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)snes), nslots*nfields, ii, PETSC_OWN_POINTER, &blocks->is);CHKERRQ(ierr);
247         } else {
248           ierr = ISCreateStride(PetscObjectComm((PetscObject)snes), nslots, rstart+blocks->fields[0], bs, &blocks->is);CHKERRQ(ierr);
249         }
250       }
251       ierr = ISSorted(blocks->is, &sorted);CHKERRQ(ierr);
252       if (!sorted) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Fields must be sorted when creating split");
253       blocks = blocks->next;
254     }
255   }
256 
257 #if 0
258   /* Create matrices */
259   ilink = jac->head;
260   if (!jac->pmat) {
261     ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr);
262     for (i=0; i<nsplit; i++) {
263       ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
264       ilink = ilink->next;
265     }
266   } else {
267     for (i=0; i<nsplit; i++) {
268       ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
269       ilink = ilink->next;
270     }
271   }
272   if (jac->realdiagonal) {
273     ilink = jac->head;
274     if (!jac->mat) {
275       ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr);
276       for (i=0; i<nsplit; i++) {
277         ierr  = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
278         ilink = ilink->next;
279       }
280     } else {
281       for (i=0; i<nsplit; i++) {
282         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
283         ilink = ilink->next;
284       }
285     }
286   } else jac->mat = jac->pmat;
287 #endif
288 
289 #if 0
290   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
291     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
292     ilink = jac->head;
293     if (!jac->Afield) {
294       ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
295       for (i=0; i<nsplit; i++) {
296         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
297         ilink = ilink->next;
298       }
299     } else {
300       for (i=0; i<nsplit; i++) {
301         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
302         ilink = ilink->next;
303       }
304     }
305   }
306 #endif
307 
308   if (mb->type == PC_COMPOSITE_SCHUR) {
309 #if 0
310     IS       ccis;
311     PetscInt rstart,rend;
312     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
313 
314     /* When extracting off-diagonal submatrices, we take complements from this range */
315     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
316 
317     /* need to handle case when one is resetting up the preconditioner */
318     if (jac->schur) {
319       ilink = jac->head;
320       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
321       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
322       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
323       ilink = ilink->next;
324       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
325       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
326       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
327       ierr  = MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1]);CHKERRQ(ierr);
328       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
329 
330     } else {
331       KSP  ksp;
332       char schurprefix[256];
333 
334       /* extract the A01 and A10 matrices */
335       ilink = jac->head;
336       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
337       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
338       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
339       ilink = ilink->next;
340       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
341       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
342       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
343       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
344       ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
345       /* set tabbing and options prefix of KSP inside the MatSchur */
346       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
347       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
348       ierr = PetscSNPrintf(schurprefix,sizeof(schurprefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",jac->head->splitname);CHKERRQ(ierr);
349       ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
350       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
351 
352       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
353       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
354       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
355       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
356       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
357         PC pc;
358         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
359         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
360         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
361       }
362       ierr = PetscSNPrintf(schurprefix,sizeof(schurprefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
363       ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
364       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
365       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
366 
367       ierr     = PetscMalloc2(2,&jac->x,2,&jac->y);CHKERRQ(ierr);
368       ierr     = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
369       ierr     = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
370       ilink    = jac->head;
371       ilink->x = jac->x[0]; ilink->y = jac->y[0];
372       ilink    = ilink->next;
373       ilink->x = jac->x[1]; ilink->y = jac->y[1];
374     }
375 #endif
376   } else {
377     /* Set up the individual SNESs */
378     blocks = mb->blocks;
379     i      = 0;
380     while (blocks) {
381       /*TODO: Set these correctly */
382       /*ierr = SNESSetFunction(blocks->snes, blocks->x, func);CHKERRQ(ierr);*/
383       /*ierr = SNESSetJacobian(blocks->snes, blocks->x, jac);CHKERRQ(ierr);*/
384       ierr = VecDuplicate(blocks->snes->vec_sol, &blocks->x);CHKERRQ(ierr);
385       /* really want setfromoptions called in SNESSetFromOptions_Multiblock(), but it is not ready yet */
386       ierr   = SNESSetFromOptions(blocks->snes);CHKERRQ(ierr);
387       ierr   = SNESSetUp(blocks->snes);CHKERRQ(ierr);
388       blocks = blocks->next;
389       i++;
390     }
391   }
392 
393   /* Compute scatter contexts needed by multiplicative versions and non-default splits */
394   if (!mb->blocks->sctx) {
395     Vec xtmp;
396 
397     blocks = mb->blocks;
398     ierr   = MatGetVecs(snes->jacobian_pre, &xtmp, NULL);CHKERRQ(ierr);
399     while (blocks) {
400       ierr   = VecScatterCreate(xtmp, blocks->is, blocks->x, NULL, &blocks->sctx);CHKERRQ(ierr);
401       blocks = blocks->next;
402     }
403     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
404   }
405   PetscFunctionReturn(0);
406 }
407 
408 /*
409   SNESSetFromOptions_Multiblock - Sets various parameters for the SNESMULTIBLOCK method.
410 
411   Input Parameter:
412 . snes - the SNES context
413 
414   Application Interface Routine: SNESSetFromOptions()
415 */
416 #undef __FUNCT__
417 #define __FUNCT__ "SNESSetFromOptions_Multiblock"
418 static PetscErrorCode SNESSetFromOptions_Multiblock(SNES snes)
419 {
420   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
421   PCCompositeType ctype;
422   PetscInt        bs;
423   PetscBool       flg;
424   PetscErrorCode  ierr;
425 
426   PetscFunctionBegin;
427   ierr = PetscOptionsHead("SNES Multiblock options");CHKERRQ(ierr);
428   ierr = PetscOptionsInt("-snes_multiblock_block_size", "Blocksize that defines number of fields", "PCFieldSplitSetBlockSize", mb->bs, &bs, &flg);CHKERRQ(ierr);
429   if (flg) {ierr = SNESMultiblockSetBlockSize(snes, bs);CHKERRQ(ierr);}
430   ierr = PetscOptionsEnum("-snes_multiblock_type", "Type of composition", "PCFieldSplitSetType", PCCompositeTypes, (PetscEnum) mb->type, (PetscEnum*) &ctype, &flg);CHKERRQ(ierr);
431   if (flg) {
432     ierr = SNESMultiblockSetType(snes,ctype);CHKERRQ(ierr);
433   }
434   /* Only setup fields once */
435   if ((mb->bs > 0) && (mb->numBlocks == 0)) {
436     /* only allow user to set fields from command line if bs is already known, otherwise user can set them in SNESMultiblockSetDefaults() */
437     ierr = SNESMultiblockSetFieldsRuntime_Private(snes);CHKERRQ(ierr);
438     if (mb->defined) {ierr = PetscInfo(snes, "Blocks defined using the options database\n");CHKERRQ(ierr);}
439   }
440   ierr = PetscOptionsTail();CHKERRQ(ierr);
441   PetscFunctionReturn(0);
442 }
443 
444 /*
445   SNESView_Multiblock - Prints info from the SNESMULTIBLOCK data structure.
446 
447   Input Parameters:
448 + SNES - the SNES context
449 - viewer - visualization context
450 
451   Application Interface Routine: SNESView()
452 */
453 #undef __FUNCT__
454 #define __FUNCT__ "SNESView_Multiblock"
455 static PetscErrorCode SNESView_Multiblock(SNES snes, PetscViewer viewer)
456 {
457   SNES_Multiblock *mb    = (SNES_Multiblock*) snes->data;
458   BlockDesc       blocks = mb->blocks;
459   PetscBool       iascii;
460   PetscErrorCode  ierr;
461 
462   PetscFunctionBegin;
463   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
464   if (iascii) {
465     ierr = PetscViewerASCIIPrintf(viewer,"  Multiblock with %s composition: total blocks = %D, blocksize = %D\n", PCCompositeTypes[mb->type], mb->numBlocks, mb->bs);CHKERRQ(ierr);
466     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following SNES objects:\n");CHKERRQ(ierr);
467     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
468     while (blocks) {
469       if (blocks->fields) {
470         PetscInt j;
471 
472         ierr = PetscViewerASCIIPrintf(viewer, "Block %s Fields ", blocks->name);CHKERRQ(ierr);
473         ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr);
474         for (j = 0; j < blocks->nfields; ++j) {
475           if (j > 0) {
476             ierr = PetscViewerASCIIPrintf(viewer, ",");CHKERRQ(ierr);
477           }
478           ierr = PetscViewerASCIIPrintf(viewer, " %D", blocks->fields[j]);CHKERRQ(ierr);
479         }
480         ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
481         ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr);
482       } else {
483         ierr = PetscViewerASCIIPrintf(viewer, "Block %s Defined by IS\n", blocks->name);CHKERRQ(ierr);
484       }
485       ierr   = SNESView(blocks->snes, viewer);CHKERRQ(ierr);
486       blocks = blocks->next;
487     }
488     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
489   }
490   PetscFunctionReturn(0);
491 }
492 
493 /*
494   SNESSolve_Multiblock - Solves a nonlinear system with the Multiblock method.
495 
496   Input Parameters:
497 . snes - the SNES context
498 
499   Output Parameter:
500 . outits - number of iterations until termination
501 
502   Application Interface Routine: SNESSolve()
503 */
504 #undef __FUNCT__
505 #define __FUNCT__ "SNESSolve_Multiblock"
506 PetscErrorCode SNESSolve_Multiblock(SNES snes)
507 {
508   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
509   Vec             X, Y, F;
510   PetscReal       fnorm;
511   PetscInt        maxits, i;
512   PetscErrorCode  ierr;
513 
514   PetscFunctionBegin;
515   snes->reason = SNES_CONVERGED_ITERATING;
516 
517   maxits = snes->max_its;        /* maximum number of iterations */
518   X      = snes->vec_sol;        /* X^n */
519   Y      = snes->vec_sol_update; /* \tilde X */
520   F      = snes->vec_func;       /* residual vector */
521 
522   ierr       = VecSetBlockSize(X, mb->bs);CHKERRQ(ierr);
523   ierr       = VecSetBlockSize(Y, mb->bs);CHKERRQ(ierr);
524   ierr       = VecSetBlockSize(F, mb->bs);CHKERRQ(ierr);
525   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
526   snes->iter = 0;
527   snes->norm = 0.;
528   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
529 
530   if (!snes->vec_func_init_set) {
531     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
532     if (snes->domainerror) {
533       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
534       PetscFunctionReturn(0);
535     }
536   } else snes->vec_func_init_set = PETSC_FALSE;
537 
538   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
539   if (PetscIsInfOrNanReal(fnorm)) {
540     snes->reason = SNES_DIVERGED_FNORM_NAN;
541     PetscFunctionReturn(0);
542   }
543   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
544   snes->norm = fnorm;
545   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
546   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
547   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
548 
549   /* test convergence */
550   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
551   if (snes->reason) PetscFunctionReturn(0);
552 
553   for (i = 0; i < maxits; i++) {
554     /* Call general purpose update function */
555     if (snes->ops->update) {
556       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
557     }
558     /* Compute X^{new} from subsolves */
559     if (mb->type == PC_COMPOSITE_ADDITIVE) {
560       BlockDesc blocks = mb->blocks;
561 
562       if (mb->defaultblocks) {
563         /*TODO: Make an array of Vecs for this */
564         /*ierr = VecStrideGatherAll(X, mb->x, INSERT_VALUES);CHKERRQ(ierr);*/
565         while (blocks) {
566           ierr   = SNESSolve(blocks->snes, NULL, blocks->x);CHKERRQ(ierr);
567           blocks = blocks->next;
568         }
569         /*ierr = VecStrideScatterAll(mb->x, X, INSERT_VALUES);CHKERRQ(ierr);*/
570       } else {
571         while (blocks) {
572           ierr   = VecScatterBegin(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
573           ierr   = VecScatterEnd(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
574           ierr   = SNESSolve(blocks->snes, NULL, blocks->x);CHKERRQ(ierr);
575           ierr   = VecScatterBegin(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
576           ierr   = VecScatterEnd(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
577           blocks = blocks->next;
578         }
579       }
580     } else SETERRQ1(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unsupported or unknown composition", (int) mb->type);
581     /* Compute F(X^{new}) */
582     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
583     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
584     if (PetscIsInfOrNanReal(fnorm)) {
585       snes->reason = SNES_DIVERGED_FNORM_NAN;
586       PetscFunctionReturn(0);
587     }
588 
589     if (snes->nfuncs >= snes->max_funcs) {
590       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
591       break;
592     }
593     if (snes->domainerror) {
594       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
595       PetscFunctionReturn(0);
596     }
597     /* Monitor convergence */
598     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
599     snes->iter = i+1;
600     snes->norm = fnorm;
601     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
602     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
603     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
604     /* Test for convergence */
605     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
606     if (snes->reason) break;
607   }
608   if (i == maxits) {
609     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
610     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
611   }
612   PetscFunctionReturn(0);
613 }
614 
615 #undef __FUNCT__
616 #define __FUNCT__ "SNESMultiblockSetFields_Default"
617 PetscErrorCode SNESMultiblockSetFields_Default(SNES snes, const char name[], PetscInt n, const PetscInt fields[])
618 {
619   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
620   BlockDesc       newblock, next = mb->blocks;
621   char            prefix[128];
622   PetscInt        i;
623   PetscErrorCode  ierr;
624 
625   PetscFunctionBegin;
626   if (mb->defined) {
627     ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr);
628     PetscFunctionReturn(0);
629   }
630   for (i = 0; i < n; ++i) {
631     if (fields[i] >= mb->bs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %D requested but only %D exist", fields[i], mb->bs);
632     if (fields[i] < 0)       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %D requested", fields[i]);
633   }
634   ierr = PetscNew(&newblock);CHKERRQ(ierr);
635   if (name) {
636     ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr);
637   } else {
638     PetscInt len = floor(log10(mb->numBlocks))+1;
639 
640     ierr = PetscMalloc1((len+1), &newblock->name);CHKERRQ(ierr);
641     ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr);
642   }
643   newblock->nfields = n;
644 
645   ierr = PetscMalloc1(n, &newblock->fields);CHKERRQ(ierr);
646   ierr = PetscMemcpy(newblock->fields, fields, n*sizeof(PetscInt));CHKERRQ(ierr);
647 
648   newblock->next = NULL;
649 
650   ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &newblock->snes);CHKERRQ(ierr);
651   ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr);
652   ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr);
653   ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr);
654   ierr = PetscSNPrintf(prefix, sizeof(prefix), "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr);
655   ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr);
656 
657   if (!next) {
658     mb->blocks         = newblock;
659     newblock->previous = NULL;
660   } else {
661     while (next->next) {
662       next = next->next;
663     }
664     next->next         = newblock;
665     newblock->previous = next;
666   }
667   mb->numBlocks++;
668   PetscFunctionReturn(0);
669 }
670 
671 #undef __FUNCT__
672 #define __FUNCT__ "SNESMultiblockSetIS_Default"
673 PetscErrorCode SNESMultiblockSetIS_Default(SNES snes, const char name[], IS is)
674 {
675   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
676   BlockDesc       newblock, next = mb->blocks;
677   char            prefix[128];
678   PetscErrorCode  ierr;
679 
680   PetscFunctionBegin;
681   if (mb->defined) {
682     ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr);
683     PetscFunctionReturn(0);
684   }
685   ierr = PetscNew(&newblock);CHKERRQ(ierr);
686   if (name) {
687     ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr);
688   } else {
689     PetscInt len = floor(log10(mb->numBlocks))+1;
690 
691     ierr = PetscMalloc1((len+1), &newblock->name);CHKERRQ(ierr);
692     ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr);
693   }
694   newblock->is = is;
695 
696   ierr = PetscObjectReference((PetscObject) is);CHKERRQ(ierr);
697 
698   newblock->next = NULL;
699 
700   ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &newblock->snes);CHKERRQ(ierr);
701   ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr);
702   ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr);
703   ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr);
704   ierr = PetscSNPrintf(prefix, sizeof(prefix), "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr);
705   ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr);
706 
707   if (!next) {
708     mb->blocks         = newblock;
709     newblock->previous = NULL;
710   } else {
711     while (next->next) {
712       next = next->next;
713     }
714     next->next         = newblock;
715     newblock->previous = next;
716   }
717   mb->numBlocks++;
718   PetscFunctionReturn(0);
719 }
720 
721 #undef __FUNCT__
722 #define __FUNCT__ "SNESMultiblockSetBlockSize_Default"
723 PetscErrorCode  SNESMultiblockSetBlockSize_Default(SNES snes, PetscInt bs)
724 {
725   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
726 
727   PetscFunctionBegin;
728   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %D", bs);
729   if (mb->bs > 0 && mb->bs != bs) SETERRQ2(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot change blocksize from %D to %D after it has been set", mb->bs, bs);
730   mb->bs = bs;
731   PetscFunctionReturn(0);
732 }
733 
734 #undef __FUNCT__
735 #define __FUNCT__ "SNESMultiblockGetSubSNES_Default"
736 PetscErrorCode SNESMultiblockGetSubSNES_Default(SNES snes, PetscInt *n, SNES **subsnes)
737 {
738   SNES_Multiblock *mb    = (SNES_Multiblock*) snes->data;
739   BlockDesc       blocks = mb->blocks;
740   PetscInt        cnt    = 0;
741   PetscErrorCode  ierr;
742 
743   PetscFunctionBegin;
744   ierr = PetscMalloc1(mb->numBlocks, subsnes);CHKERRQ(ierr);
745   while (blocks) {
746     (*subsnes)[cnt++] = blocks->snes;
747     blocks            = blocks->next;
748   }
749   if (cnt != mb->numBlocks) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Corrupt SNESMULTIBLOCK object: number of blocks in linked list %D does not match number in object %D", cnt, mb->numBlocks);
750 
751   if (n) *n = mb->numBlocks;
752   PetscFunctionReturn(0);
753 }
754 
755 #undef __FUNCT__
756 #define __FUNCT__ "SNESMultiblockSetType_Default"
757 PetscErrorCode  SNESMultiblockSetType_Default(SNES snes, PCCompositeType type)
758 {
759   SNES_Multiblock *mb = (SNES_Multiblock*) snes->data;
760   PetscErrorCode  ierr;
761 
762   PetscFunctionBegin;
763   mb->type = type;
764   if (type == PC_COMPOSITE_SCHUR) {
765 #if 1
766     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "The Schur composite type is not yet supported");
767 #else
768     snes->ops->solve = SNESSolve_Multiblock_Schur;
769     snes->ops->view  = SNESView_Multiblock_Schur;
770 
771     ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockGetSubSNES_C", SNESMultiblockGetSubSNES_Schur);CHKERRQ(ierr);
772     ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", SNESMultiblockSchurPrecondition_Default);CHKERRQ(ierr);
773 #endif
774   } else {
775     snes->ops->solve = SNESSolve_Multiblock;
776     snes->ops->view  = SNESView_Multiblock;
777 
778     ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockGetSubSNES_C", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr);
779     ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", 0);CHKERRQ(ierr);
780   }
781   PetscFunctionReturn(0);
782 }
783 
784 #undef __FUNCT__
785 #define __FUNCT__ "SNESMultiblockSetFields"
786 /*@
787   SNESMultiblockSetFields - Sets the fields for one particular block in the solver
788 
789   Logically Collective on SNES
790 
791   Input Parameters:
792 + snes   - the solver
793 . name   - name of this block, if NULL the number of the block is used
794 . n      - the number of fields in this block
795 - fields - the fields in this block
796 
797   Level: intermediate
798 
799   Notes: Use SNESMultiblockSetIS() to set a completely general set of row indices as a block.
800 
801   The SNESMultiblockSetFields() is for defining blocks as a group of strided indices, or fields.
802   For example, if the vector block size is three then one can define a block as field 0, or
803   1 or 2, or field 0,1 or 0,2 or 1,2 which means
804     0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
805   where the numbered entries indicate what is in the block.
806 
807   This function is called once per block (it creates a new block each time). Solve options
808   for this block will be available under the prefix -multiblock_BLOCKNAME_.
809 
810 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetBlockSize(), SNESMultiblockSetIS()
811 @*/
812 PetscErrorCode SNESMultiblockSetFields(SNES snes, const char name[], PetscInt n, const PetscInt *fields)
813 {
814   PetscErrorCode ierr;
815 
816   PetscFunctionBegin;
817   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
818   PetscValidCharPointer(name, 2);
819   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Provided number of fields %D in split \"%s\" not positive", n, name);
820   PetscValidIntPointer(fields, 3);
821   ierr = PetscTryMethod(snes, "SNESMultiblockSetFields_C", (SNES, const char[], PetscInt, const PetscInt*), (snes, name, n, fields));CHKERRQ(ierr);
822   PetscFunctionReturn(0);
823 }
824 
825 #undef __FUNCT__
826 #define __FUNCT__ "SNESMultiblockSetIS"
827 /*@
828   SNESMultiblockSetIS - Sets the global row indices for the block
829 
830   Logically Collective on SNES
831 
832   Input Parameters:
833 + snes - the solver context
834 . name - name of this block, if NULL the number of the block is used
835 - is   - the index set that defines the global row indices in this block
836 
837   Notes:
838   Use SNESMultiblockSetFields(), for blocks defined by strides.
839 
840   This function is called once per block (it creates a new block each time). Solve options
841   for this block will be available under the prefix -multiblock_BLOCKNAME_.
842 
843   Level: intermediate
844 
845 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetBlockSize()
846 @*/
847 PetscErrorCode SNESMultiblockSetIS(SNES snes, const char name[], IS is)
848 {
849   PetscErrorCode ierr;
850 
851   PetscFunctionBegin;
852   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
853   PetscValidCharPointer(name, 2);
854   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
855   ierr = PetscTryMethod(snes, "SNESMultiblockSetIS_C", (SNES, const char[], IS), (snes, name, is));CHKERRQ(ierr);
856   PetscFunctionReturn(0);
857 }
858 
859 #undef __FUNCT__
860 #define __FUNCT__ "SNESMultiblockSetType"
861 /*@
862   SNESMultiblockSetType - Sets the type of block combination.
863 
864   Collective on SNES
865 
866   Input Parameters:
867 + snes - the solver context
868 - type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
869 
870   Options Database Key:
871 . -snes_multiblock_type <type: one of multiplicative, additive, symmetric_multiplicative> - Sets block combination type
872 
873   Level: Developer
874 
875 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
876 .seealso: PCCompositeSetType()
877 @*/
878 PetscErrorCode SNESMultiblockSetType(SNES snes, PCCompositeType type)
879 {
880   PetscErrorCode ierr;
881 
882   PetscFunctionBegin;
883   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
884   ierr = PetscTryMethod(snes, "SNESMultiblockSetType_C", (SNES, PCCompositeType), (snes, type));CHKERRQ(ierr);
885   PetscFunctionReturn(0);
886 }
887 
888 #undef __FUNCT__
889 #define __FUNCT__ "SNESMultiblockSetBlockSize"
890 /*@
891   SNESMultiblockSetBlockSize - Sets the block size for structured mesh block division. If not set the matrix block size is used.
892 
893   Logically Collective on SNES
894 
895   Input Parameters:
896 + snes - the solver context
897 - bs   - the block size
898 
899   Level: intermediate
900 
901 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetFields()
902 @*/
903 PetscErrorCode SNESMultiblockSetBlockSize(SNES snes, PetscInt bs)
904 {
905   PetscErrorCode ierr;
906 
907   PetscFunctionBegin;
908   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
909   PetscValidLogicalCollectiveInt(snes, bs, 2);
910   ierr = PetscTryMethod(snes, "SNESMultiblockSetBlockSize_C", (SNES, PetscInt), (snes,bs));CHKERRQ(ierr);
911   PetscFunctionReturn(0);
912 }
913 
914 #undef __FUNCT__
915 #define __FUNCT__ "SNESMultiblockGetSubSNES"
916 /*@C
917   SNESMultiblockGetSubSNES - Gets the SNES contexts for all blocks
918 
919   Collective on SNES
920 
921   Input Parameter:
922 . snes - the solver context
923 
924   Output Parameters:
925 + n       - the number of blocks
926 - subsnes - the array of SNES contexts
927 
928   Note:
929   After SNESMultiblockGetSubSNES() the array of SNESs MUST be freed by the user
930   (not each SNES, just the array that contains them).
931 
932   You must call SNESSetUp() before calling SNESMultiblockGetSubSNES().
933 
934   Level: advanced
935 
936 .seealso: SNESMULTIBLOCK
937 @*/
938 PetscErrorCode SNESMultiblockGetSubSNES(SNES snes, PetscInt *n, SNES *subsnes[])
939 {
940   PetscErrorCode ierr;
941 
942   PetscFunctionBegin;
943   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
944   if (n) PetscValidIntPointer(n, 2);
945   ierr = PetscUseMethod(snes, "SNESMultiblockGetSubSNES_C", (SNES, PetscInt*, SNES **), (snes, n, subsnes));CHKERRQ(ierr);
946   PetscFunctionReturn(0);
947 }
948 
949 /*MC
950   SNESMULTIBLOCK - Multiblock nonlinear solver that can use overlapping or nonoverlapping blocks, organized
951   additively (Jacobi) or multiplicatively (Gauss-Seidel).
952 
953   Level: beginner
954 
955 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNRICHARDSON
956 M*/
957 #undef __FUNCT__
958 #define __FUNCT__ "SNESCreate_Multiblock"
959 PETSC_EXTERN PetscErrorCode SNESCreate_Multiblock(SNES snes)
960 {
961   SNES_Multiblock *mb;
962   PetscErrorCode  ierr;
963 
964   PetscFunctionBegin;
965   snes->ops->destroy        = SNESDestroy_Multiblock;
966   snes->ops->setup          = SNESSetUp_Multiblock;
967   snes->ops->setfromoptions = SNESSetFromOptions_Multiblock;
968   snes->ops->view           = SNESView_Multiblock;
969   snes->ops->solve          = SNESSolve_Multiblock;
970   snes->ops->reset          = SNESReset_Multiblock;
971 
972   snes->usesksp = PETSC_FALSE;
973 
974   ierr          = PetscNewLog(snes,&mb);CHKERRQ(ierr);
975   snes->data    = (void*) mb;
976   mb->defined   = PETSC_FALSE;
977   mb->numBlocks = 0;
978   mb->bs        = -1;
979   mb->type      = PC_COMPOSITE_MULTIPLICATIVE;
980 
981   /* We attach functions so that they can be called on another PC without crashing the program */
982   ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetFields_C",    SNESMultiblockSetFields_Default);CHKERRQ(ierr);
983   ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetIS_C",        SNESMultiblockSetIS_Default);CHKERRQ(ierr);
984   ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetType_C",      SNESMultiblockSetType_Default);CHKERRQ(ierr);
985   ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetBlockSize_C", SNESMultiblockSetBlockSize_Default);CHKERRQ(ierr);
986   ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockGetSubSNES_C",   SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr);
987   PetscFunctionReturn(0);
988 }
989