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