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