Actual source code: multiblock.c

  1: #include <petsc/private/snesimpl.h>
  2: #include <petscdmcomposite.h>

  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: };

 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;

 26: PetscErrorCode SNESReset_Multiblock(SNES snes)
 27: {
 28:   SNES_Multiblock *mb     = (SNES_Multiblock *)snes->data;
 29:   BlockDesc        blocks = mb->blocks, next;

 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(PETSC_SUCCESS);
 43: }

 45: /*
 46:   SNESDestroy_Multiblock - Destroys the private SNES_Multiblock context that was created with SNESCreate_Multiblock().

 48:   Input Parameter:
 49: . snes - the SNES context

 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;

 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(PETSC_SUCCESS);
 70: }

 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];

 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(PETSC_SUCCESS);
100: }

102: static PetscErrorCode SNESMultiblockSetDefaults(SNES snes)
103: {
104:   SNES_Multiblock *mb     = (SNES_Multiblock *)snes->data;
105:   BlockDesc        blocks = mb->blocks;
106:   PetscInt         i;

108:   PetscFunctionBegin;
109:   if (!blocks) {
110:     if (snes->dm) {
111:       PetscBool dmcomposite;

113:       PetscCall(PetscObjectTypeCompare((PetscObject)snes->dm, DMCOMPOSITE, &dmcomposite));
114:       if (dmcomposite) {
115:         PetscInt nDM;
116:         IS      *fields;

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];

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;

134:       if (mb->bs <= 0) {
135:         if (snes->jacobian_pre) {
136:           PetscCall(MatGetBlockSize(snes->jacobian_pre, &mb->bs));
137:         } else mb->bs = 1;
138:       }

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;

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];

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;

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(PETSC_SUCCESS);
185: }

187: /*
188:    SNESSetUp_Multiblock - Sets up the internal data structures for the later use
189:    of the SNESMULTIBLOCK nonlinear solver.

191:    Input Parameters:
192: +  snes - the SNES context
193: -  x - the solution vector

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;

203:   PetscFunctionBegin;
204:   PetscCall(SNESMultiblockSetDefaults(snes));
205:   numBlocks = mb->numBlocks;
206:   blocks    = mb->blocks;

208:   /* Create ISs */
209:   if (!mb->issetup) {
210:     PetscInt  ccsize, rstart, rend, nslots, bs;
211:     PetscBool sorted;

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;

225:           PetscCall(PetscMalloc1(nfields * nslots, &ii));
226:           for (j = 0; j < nslots; ++j) {
227:             for (k = 0; k < nfields; ++k) ii[nfields * j + k] = rstart + bs * j + fields[k];
228:           }
229:           PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)snes), nslots * nfields, ii, PETSC_OWN_POINTER, &blocks->is));
230:         } else {
231:           PetscCall(ISCreateStride(PetscObjectComm((PetscObject)snes), nslots, rstart + blocks->fields[0], bs, &blocks->is));
232:         }
233:       }
234:       PetscCall(ISSorted(blocks->is, &sorted));
235:       PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_USER, "Fields must be sorted when creating split");
236:       blocks = blocks->next;
237:     }
238:   }

240: #if 0
241:   /* Create matrices */
242:   ilink = jac->head;
243:   if (!jac->pmat) {
244:     PetscCall(PetscMalloc1(nsplit,&jac->pmat));
245:     for (i=0; i<nsplit; i++) {
246:       PetscCall(MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]));
247:       ilink = ilink->next;
248:     }
249:   } else {
250:     for (i=0; i<nsplit; i++) {
251:       PetscCall(MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]));
252:       ilink = ilink->next;
253:     }
254:   }
255:   if (jac->realdiagonal) {
256:     ilink = jac->head;
257:     if (!jac->mat) {
258:       PetscCall(PetscMalloc1(nsplit,&jac->mat));
259:       for (i=0; i<nsplit; i++) {
260:         PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]));
261:         ilink = ilink->next;
262:       }
263:     } else {
264:       for (i=0; i<nsplit; i++) {
265:         if (jac->mat[i]) PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]));
266:         ilink = ilink->next;
267:       }
268:     }
269:   } else jac->mat = jac->pmat;
270: #endif

272: #if 0
273:   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
274:     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
275:     ilink = jac->head;
276:     if (!jac->Afield) {
277:       PetscCall(PetscMalloc1(nsplit,&jac->Afield));
278:       for (i=0; i<nsplit; i++) {
279:         PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]));
280:         ilink = ilink->next;
281:       }
282:     } else {
283:       for (i=0; i<nsplit; i++) {
284:         PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]));
285:         ilink = ilink->next;
286:       }
287:     }
288:   }
289: #endif

291:   if (mb->type == PC_COMPOSITE_SCHUR) {
292: #if 0
293:     IS       ccis;
294:     PetscInt rstart,rend;
295:     PetscCheck(nsplit == 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");

297:     /* When extracting off-diagonal submatrices, we take complements from this range */
298:     PetscCall(MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend));

300:     /* need to handle case when one is resetting up the preconditioner */
301:     if (jac->schur) {
302:       ilink = jac->head;
303:       PetscCall(ISComplement(ilink->is,rstart,rend,&ccis));
304:       PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B));
305:       PetscCall(ISDestroy(&ccis));
306:       ilink = ilink->next;
307:       PetscCall(ISComplement(ilink->is,rstart,rend,&ccis));
308:       PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C));
309:       PetscCall(ISDestroy(&ccis));
310:       PetscCall(MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1]));
311:       PetscCall(KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag));

313:     } else {
314:       KSP  ksp;
315:       char schurprefix[256];

317:       /* extract the A01 and A10 matrices */
318:       ilink = jac->head;
319:       PetscCall(ISComplement(ilink->is,rstart,rend,&ccis));
320:       PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B));
321:       PetscCall(ISDestroy(&ccis));
322:       ilink = ilink->next;
323:       PetscCall(ISComplement(ilink->is,rstart,rend,&ccis));
324:       PetscCall(MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C));
325:       PetscCall(ISDestroy(&ccis));
326:       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
327:       PetscCall(MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur));
328:       /* set tabbing and options prefix of KSP inside the MatSchur */
329:       PetscCall(MatSchurComplementGetKSP(jac->schur,&ksp));
330:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2));
331:       PetscCall(PetscSNPrintf(schurprefix,sizeof(schurprefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",jac->head->splitname));
332:       PetscCall(KSPSetOptionsPrefix(ksp,schurprefix));
333:       PetscCall(MatSetFromOptions(jac->schur));

335:       PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur));
336:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1));
337:       PetscCall(KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac)));
338:       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
339:         PC pc;
340:         PetscCall(KSPGetPC(jac->kspschur,&pc));
341:         PetscCall(PCSetType(pc,PCNONE));
342:         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
343:       }
344:       PetscCall(PetscSNPrintf(schurprefix,sizeof(schurprefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname));
345:       PetscCall(KSPSetOptionsPrefix(jac->kspschur,schurprefix));
346:       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
347:       PetscCall(KSPSetFromOptions(jac->kspschur));

349:       PetscCall(PetscMalloc2(2,&jac->x,2,&jac->y));
350:       PetscCall(MatCreateVecs(jac->pmat[0],&jac->x[0],&jac->y[0]));
351:       PetscCall(MatCreateVecs(jac->pmat[1],&jac->x[1],&jac->y[1]));
352:       ilink    = jac->head;
353:       ilink->x = jac->x[0]; ilink->y = jac->y[0];
354:       ilink    = ilink->next;
355:       ilink->x = jac->x[1]; ilink->y = jac->y[1];
356:     }
357: #endif
358:   } else {
359:     /* Set up the individual SNESs */
360:     blocks = mb->blocks;
361:     i      = 0;
362:     while (blocks) {
363:       /*TODO: Set these correctly */
364:       /* PetscCall(SNESSetFunction(blocks->snes, blocks->x, func)); */
365:       /* PetscCall(SNESSetJacobian(blocks->snes, blocks->x, jac)); */
366:       PetscCall(VecDuplicate(blocks->snes->vec_sol, &blocks->x));
367:       /* really want setfromoptions called in SNESSetFromOptions_Multiblock(), but it is not ready yet */
368:       PetscCall(SNESSetFromOptions(blocks->snes));
369:       PetscCall(SNESSetUp(blocks->snes));
370:       blocks = blocks->next;
371:       i++;
372:     }
373:   }

375:   /* Compute scatter contexts needed by multiplicative versions and non-default splits */
376:   if (!mb->blocks->sctx) {
377:     Vec xtmp;

379:     blocks = mb->blocks;
380:     PetscCall(MatCreateVecs(snes->jacobian_pre, &xtmp, NULL));
381:     while (blocks) {
382:       PetscCall(VecScatterCreate(xtmp, blocks->is, blocks->x, NULL, &blocks->sctx));
383:       blocks = blocks->next;
384:     }
385:     PetscCall(VecDestroy(&xtmp));
386:   }
387:   PetscFunctionReturn(PETSC_SUCCESS);
388: }

390: /*
391:   SNESSetFromOptions_Multiblock - Sets various parameters for the SNESMULTIBLOCK method.

393:   Input Parameter:
394: . snes - the SNES context

396:   Application Interface Routine: SNESSetFromOptions()
397: */
398: static PetscErrorCode SNESSetFromOptions_Multiblock(SNES snes, PetscOptionItems *PetscOptionsObject)
399: {
400:   SNES_Multiblock *mb = (SNES_Multiblock *)snes->data;
401:   PCCompositeType  ctype;
402:   PetscInt         bs;
403:   PetscBool        flg;

405:   PetscFunctionBegin;
406:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES Multiblock options");
407:   PetscCall(PetscOptionsInt("-snes_multiblock_block_size", "Blocksize that defines number of fields", "PCFieldSplitSetBlockSize", mb->bs, &bs, &flg));
408:   if (flg) PetscCall(SNESMultiblockSetBlockSize(snes, bs));
409:   PetscCall(PetscOptionsEnum("-snes_multiblock_type", "Type of composition", "PCFieldSplitSetType", PCCompositeTypes, (PetscEnum)mb->type, (PetscEnum *)&ctype, &flg));
410:   if (flg) PetscCall(SNESMultiblockSetType(snes, ctype));
411:   /* Only setup fields once */
412:   if ((mb->bs > 0) && (mb->numBlocks == 0)) {
413:     /* only allow user to set fields from command line if bs is already known, otherwise user can set them in SNESMultiblockSetDefaults() */
414:     PetscCall(SNESMultiblockSetFieldsRuntime_Private(snes));
415:     if (mb->defined) PetscCall(PetscInfo(snes, "Blocks defined using the options database\n"));
416:   }
417:   PetscOptionsHeadEnd();
418:   PetscFunctionReturn(PETSC_SUCCESS);
419: }

421: /*
422:   SNESView_Multiblock - Prints info from the SNESMULTIBLOCK data structure.

424:   Input Parameters:
425: + SNES - the SNES context
426: - viewer - visualization context

428:   Application Interface Routine: SNESView()
429: */
430: static PetscErrorCode SNESView_Multiblock(SNES snes, PetscViewer viewer)
431: {
432:   SNES_Multiblock *mb     = (SNES_Multiblock *)snes->data;
433:   BlockDesc        blocks = mb->blocks;
434:   PetscBool        iascii;

436:   PetscFunctionBegin;
437:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
438:   if (iascii) {
439:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Multiblock with %s composition: total blocks = %" PetscInt_FMT ", blocksize = %" PetscInt_FMT "\n", PCCompositeTypes[mb->type], mb->numBlocks, mb->bs));
440:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Solver info for each split is in the following SNES objects:\n"));
441:     PetscCall(PetscViewerASCIIPushTab(viewer));
442:     while (blocks) {
443:       if (blocks->fields) {
444:         PetscInt j;

446:         PetscCall(PetscViewerASCIIPrintf(viewer, "  Block %s Fields ", blocks->name));
447:         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
448:         for (j = 0; j < blocks->nfields; ++j) {
449:           if (j > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ","));
450:           PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, blocks->fields[j]));
451:         }
452:         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
453:         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
454:       } else {
455:         PetscCall(PetscViewerASCIIPrintf(viewer, "  Block %s Defined by IS\n", blocks->name));
456:       }
457:       PetscCall(SNESView(blocks->snes, viewer));
458:       blocks = blocks->next;
459:     }
460:     PetscCall(PetscViewerASCIIPopTab(viewer));
461:   }
462:   PetscFunctionReturn(PETSC_SUCCESS);
463: }

465: /*
466:   SNESSolve_Multiblock - Solves a nonlinear system with the Multiblock method.

468:   Input Parameters:
469: . snes - the SNES context

471:   Output Parameter:
472: . outits - number of iterations until termination

474:   Application Interface Routine: SNESSolve()
475: */
476: PetscErrorCode SNESSolve_Multiblock(SNES snes)
477: {
478:   SNES_Multiblock *mb = (SNES_Multiblock *)snes->data;
479:   Vec              X, Y, F;
480:   PetscReal        fnorm;
481:   PetscInt         maxits, i;

483:   PetscFunctionBegin;
484:   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);

486:   snes->reason = SNES_CONVERGED_ITERATING;

488:   maxits = snes->max_its;        /* maximum number of iterations */
489:   X      = snes->vec_sol;        /* X^n */
490:   Y      = snes->vec_sol_update; /* \tilde X */
491:   F      = snes->vec_func;       /* residual vector */

493:   PetscCall(VecSetBlockSize(X, mb->bs));
494:   PetscCall(VecSetBlockSize(Y, mb->bs));
495:   PetscCall(VecSetBlockSize(F, mb->bs));
496:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
497:   snes->iter = 0;
498:   snes->norm = 0.;
499:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

501:   if (!snes->vec_func_init_set) {
502:     PetscCall(SNESComputeFunction(snes, X, F));
503:   } else snes->vec_func_init_set = PETSC_FALSE;

505:   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
506:   SNESCheckFunctionNorm(snes, fnorm);
507:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
508:   snes->norm = fnorm;
509:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
510:   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
511:   PetscCall(SNESMonitor(snes, 0, fnorm));

513:   /* test convergence */
514:   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
515:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);

517:   for (i = 0; i < maxits; i++) {
518:     /* Call general purpose update function */
519:     PetscTryTypeMethod(snes, update, snes->iter);
520:     /* Compute X^{new} from subsolves */
521:     if (mb->type == PC_COMPOSITE_ADDITIVE) {
522:       BlockDesc blocks = mb->blocks;

524:       if (mb->defaultblocks) {
525:         /*TODO: Make an array of Vecs for this */
526:         /* PetscCall(VecStrideGatherAll(X, mb->x, INSERT_VALUES)); */
527:         while (blocks) {
528:           PetscCall(SNESSolve(blocks->snes, NULL, blocks->x));
529:           blocks = blocks->next;
530:         }
531:         /* PetscCall(VecStrideScatterAll(mb->x, X, INSERT_VALUES)); */
532:       } else {
533:         while (blocks) {
534:           PetscCall(VecScatterBegin(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD));
535:           PetscCall(VecScatterEnd(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD));
536:           PetscCall(SNESSolve(blocks->snes, NULL, blocks->x));
537:           PetscCall(VecScatterBegin(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE));
538:           PetscCall(VecScatterEnd(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE));
539:           blocks = blocks->next;
540:         }
541:       }
542:     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unsupported or unknown composition %d", (int)mb->type);
543:     /* Compute F(X^{new}) */
544:     PetscCall(SNESComputeFunction(snes, X, F));
545:     PetscCall(VecNorm(F, NORM_2, &fnorm));
546:     SNESCheckFunctionNorm(snes, fnorm);

548:     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
549:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
550:       break;
551:     }

553:     /* Monitor convergence */
554:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
555:     snes->iter = i + 1;
556:     snes->norm = fnorm;
557:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
558:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
559:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
560:     /* Test for convergence */
561:     PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
562:     if (snes->reason) break;
563:   }
564:   if (i == maxits) {
565:     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
566:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
567:   }
568:   PetscFunctionReturn(PETSC_SUCCESS);
569: }

571: PetscErrorCode SNESMultiblockSetFields_Default(SNES snes, const char name[], PetscInt n, const PetscInt fields[])
572: {
573:   SNES_Multiblock *mb = (SNES_Multiblock *)snes->data;
574:   BlockDesc        newblock, next = mb->blocks;
575:   char             prefix[128];
576:   PetscInt         i;

578:   PetscFunctionBegin;
579:   if (mb->defined) {
580:     PetscCall(PetscInfo(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name));
581:     PetscFunctionReturn(PETSC_SUCCESS);
582:   }
583:   for (i = 0; i < n; ++i) {
584:     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);
585:     PetscCheck(fields[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %" PetscInt_FMT " requested", fields[i]);
586:   }
587:   PetscCall(PetscNew(&newblock));
588:   if (name) {
589:     PetscCall(PetscStrallocpy(name, &newblock->name));
590:   } else {
591:     PetscInt len = floor(log10(mb->numBlocks)) + 1;

593:     PetscCall(PetscMalloc1(len + 1, &newblock->name));
594:     PetscCall(PetscSNPrintf(newblock->name, len, "%" PetscInt_FMT, mb->numBlocks));
595:   }
596:   newblock->nfields = n;

598:   PetscCall(PetscMalloc1(n, &newblock->fields));
599:   PetscCall(PetscArraycpy(newblock->fields, fields, n));

601:   newblock->next = NULL;

603:   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &newblock->snes));
604:   PetscCall(PetscObjectIncrementTabLevel((PetscObject)newblock->snes, (PetscObject)snes, 1));
605:   PetscCall(SNESSetType(newblock->snes, SNESNRICHARDSON));
606:   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%smultiblock_%s_", ((PetscObject)snes)->prefix ? ((PetscObject)snes)->prefix : "", newblock->name));
607:   PetscCall(SNESSetOptionsPrefix(newblock->snes, prefix));

609:   if (!next) {
610:     mb->blocks         = newblock;
611:     newblock->previous = NULL;
612:   } else {
613:     while (next->next) next = next->next;
614:     next->next         = newblock;
615:     newblock->previous = next;
616:   }
617:   mb->numBlocks++;
618:   PetscFunctionReturn(PETSC_SUCCESS);
619: }

621: PetscErrorCode SNESMultiblockSetIS_Default(SNES snes, const char name[], IS is)
622: {
623:   SNES_Multiblock *mb = (SNES_Multiblock *)snes->data;
624:   BlockDesc        newblock, next = mb->blocks;
625:   char             prefix[128];

627:   PetscFunctionBegin;
628:   if (mb->defined) {
629:     PetscCall(PetscInfo(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name));
630:     PetscFunctionReturn(PETSC_SUCCESS);
631:   }
632:   PetscCall(PetscNew(&newblock));
633:   if (name) {
634:     PetscCall(PetscStrallocpy(name, &newblock->name));
635:   } else {
636:     PetscInt len = floor(log10(mb->numBlocks)) + 1;

638:     PetscCall(PetscMalloc1(len + 1, &newblock->name));
639:     PetscCall(PetscSNPrintf(newblock->name, len, "%" PetscInt_FMT, mb->numBlocks));
640:   }
641:   newblock->is = is;

643:   PetscCall(PetscObjectReference((PetscObject)is));

645:   newblock->next = NULL;

647:   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &newblock->snes));
648:   PetscCall(PetscObjectIncrementTabLevel((PetscObject)newblock->snes, (PetscObject)snes, 1));
649:   PetscCall(SNESSetType(newblock->snes, SNESNRICHARDSON));
650:   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%smultiblock_%s_", ((PetscObject)snes)->prefix ? ((PetscObject)snes)->prefix : "", newblock->name));
651:   PetscCall(SNESSetOptionsPrefix(newblock->snes, prefix));

653:   if (!next) {
654:     mb->blocks         = newblock;
655:     newblock->previous = NULL;
656:   } else {
657:     while (next->next) next = next->next;
658:     next->next         = newblock;
659:     newblock->previous = next;
660:   }
661:   mb->numBlocks++;
662:   PetscFunctionReturn(PETSC_SUCCESS);
663: }

665: PetscErrorCode SNESMultiblockSetBlockSize_Default(SNES snes, PetscInt bs)
666: {
667:   SNES_Multiblock *mb = (SNES_Multiblock *)snes->data;

669:   PetscFunctionBegin;
670:   PetscCheck(bs >= 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %" PetscInt_FMT, bs);
671:   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);
672:   mb->bs = bs;
673:   PetscFunctionReturn(PETSC_SUCCESS);
674: }

676: PetscErrorCode SNESMultiblockGetSubSNES_Default(SNES snes, PetscInt *n, SNES **subsnes)
677: {
678:   SNES_Multiblock *mb     = (SNES_Multiblock *)snes->data;
679:   BlockDesc        blocks = mb->blocks;
680:   PetscInt         cnt    = 0;

682:   PetscFunctionBegin;
683:   PetscCall(PetscMalloc1(mb->numBlocks, subsnes));
684:   while (blocks) {
685:     (*subsnes)[cnt++] = blocks->snes;
686:     blocks            = blocks->next;
687:   }
688:   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);

690:   if (n) *n = mb->numBlocks;
691:   PetscFunctionReturn(PETSC_SUCCESS);
692: }

694: PetscErrorCode SNESMultiblockSetType_Default(SNES snes, PCCompositeType type)
695: {
696:   SNES_Multiblock *mb = (SNES_Multiblock *)snes->data;

698:   PetscFunctionBegin;
699:   mb->type = type;
700:   if (type == PC_COMPOSITE_SCHUR) {
701: #if 1
702:     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "The Schur composite type is not yet supported");
703: #else
704:     snes->ops->solve = SNESSolve_Multiblock_Schur;
705:     snes->ops->view  = SNESView_Multiblock_Schur;

707:     PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockGetSubSNES_C", SNESMultiblockGetSubSNES_Schur));
708:     PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockSchurPrecondition_C", SNESMultiblockSchurPrecondition_Default));
709: #endif
710:   } else {
711:     snes->ops->solve = SNESSolve_Multiblock;
712:     snes->ops->view  = SNESView_Multiblock;

714:     PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockGetSubSNES_C", SNESMultiblockGetSubSNES_Default));
715:     PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockSchurPrecondition_C", 0));
716:   }
717:   PetscFunctionReturn(PETSC_SUCCESS);
718: }

720: /*@
721:   SNESMultiblockSetFields - Sets the fields for one particular block in a `SNESMULTBLOCK` solver

723:   Logically Collective

725:   Input Parameters:
726: + snes   - the solver
727: . name   - name of this block, if NULL the number of the block is used
728: . n      - the number of fields in this block
729: - fields - the fields in this block

731:   Level: intermediate

733:   Notes:
734:     Use `SNESMultiblockSetIS()` to set a completely general set of row indices as a block.

736:   The `SNESMultiblockSetFields()` is for defining blocks as a group of strided indices, or fields.
737:   For example, if the vector block size is three then one can define a block as field 0, or
738:   1 or 2, or field 0,1 or 0,2 or 1,2 which means
739:     0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
740:   where the numbered entries indicate what is in the block.

742:   This function is called once per block (it creates a new block each time). Solve options
743:   for this block will be available under the prefix -multiblock_BLOCKNAME_.

745: .seealso: `SNESMULTBLOCK`, `SNESMultiblockGetSubSNES()`, `SNESMULTIBLOCK`, `SNESMultiblockSetBlockSize()`, `SNESMultiblockSetIS()`
746: @*/
747: PetscErrorCode SNESMultiblockSetFields(SNES snes, const char name[], PetscInt n, const PetscInt *fields)
748: {
749:   PetscFunctionBegin;
752:   PetscCheck(n >= 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Provided number of fields %" PetscInt_FMT " in split \"%s\" not positive", n, name);
754:   PetscTryMethod(snes, "SNESMultiblockSetFields_C", (SNES, const char[], PetscInt, const PetscInt *), (snes, name, n, fields));
755:   PetscFunctionReturn(PETSC_SUCCESS);
756: }

758: /*@
759:   SNESMultiblockSetIS - Sets the global row indices for one particular block in a `SNESMULTBLOCK` solver

761:   Logically Collective

763:   Input Parameters:
764: + snes - the solver context
765: . name - name of this block, if NULL the number of the block is used
766: - is   - the index set that defines the global row indices in this block

768:   Notes:
769:   Use `SNESMultiblockSetFields()`, for blocks defined by strides.

771:   This function is called once per block (it creates a new block each time). Solve options
772:   for this block will be available under the prefix -multiblock_BLOCKNAME_.

774:   Level: intermediate

776: .seealso: `SNESMULTBLOCK`, `SNESMultiblockGetSubSNES()`, `SNESMULTIBLOCK`, `SNESMultiblockSetBlockSize()`
777: @*/
778: PetscErrorCode SNESMultiblockSetIS(SNES snes, const char name[], IS is)
779: {
780:   PetscFunctionBegin;
784:   PetscTryMethod(snes, "SNESMultiblockSetIS_C", (SNES, const char[], IS), (snes, name, is));
785:   PetscFunctionReturn(PETSC_SUCCESS);
786: }

788: /*@
789:   SNESMultiblockSetType - Sets the type of block combination used for a `SNESMULTBLOCK` solver

791:   Logically Collective

793:   Input Parameters:
794: + snes - the solver context
795: - type - `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE` (default), `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`

797:   Options Database Key:
798: . -snes_multiblock_type <type: one of multiplicative, additive, symmetric_multiplicative> - Sets block combination type

800:   Level: advanced

802: .seealso: `SNESMULTBLOCK`, `PCCompositeSetType()`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`
803: @*/
804: PetscErrorCode SNESMultiblockSetType(SNES snes, PCCompositeType type)
805: {
806:   PetscFunctionBegin;
808:   PetscTryMethod(snes, "SNESMultiblockSetType_C", (SNES, PCCompositeType), (snes, type));
809:   PetscFunctionReturn(PETSC_SUCCESS);
810: }

812: /*@
813:   SNESMultiblockSetBlockSize - Sets the block size for structured block division in a `SNESMULTBLOCK` solver. If not set the matrix block size is used.

815:   Logically Collective

817:   Input Parameters:
818: + snes - the solver context
819: - bs   - the block size

821:   Level: intermediate

823: .seealso: `SNESMULTBLOCK`, `SNESMultiblockGetSubSNES()`, `SNESMULTIBLOCK`, `SNESMultiblockSetFields()`
824: @*/
825: PetscErrorCode SNESMultiblockSetBlockSize(SNES snes, PetscInt bs)
826: {
827:   PetscFunctionBegin;
830:   PetscTryMethod(snes, "SNESMultiblockSetBlockSize_C", (SNES, PetscInt), (snes, bs));
831:   PetscFunctionReturn(PETSC_SUCCESS);
832: }

834: /*@C
835:   SNESMultiblockGetSubSNES - Gets the `SNES` contexts for all blocks in a `SNESMULTBLOCK` solver.

837:   Not Collective but each `SNES` obtained is parallel

839:   Input Parameter:
840: . snes - the solver context

842:   Output Parameters:
843: + n       - the number of blocks
844: - subsnes - the array of `SNES` contexts

846:   Note:
847:   After `SNESMultiblockGetSubSNES()` the array of `SNES`s MUST be freed by the user
848:   (not each `SNES`, just the array that contains them).

850:   You must call `SNESSetUp()` before calling `SNESMultiblockGetSubSNES()`.

852:   Level: advanced

854: .seealso: `SNESMULTBLOCK`, `SNESMultiblockSetIS()`, `SNESMultiblockSetFields()`
855: @*/
856: PetscErrorCode SNESMultiblockGetSubSNES(SNES snes, PetscInt *n, SNES *subsnes[])
857: {
858:   PetscFunctionBegin;
861:   PetscUseMethod(snes, "SNESMultiblockGetSubSNES_C", (SNES, PetscInt *, SNES **), (snes, n, subsnes));
862:   PetscFunctionReturn(PETSC_SUCCESS);
863: }

865: /*MC
866:   SNESMULTIBLOCK - Multiblock nonlinear solver that can use overlapping or nonoverlapping blocks, organized
867:   additively (Jacobi) or multiplicatively (Gauss-Seidel).

869:   Level: beginner

871: .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNRICHARDSON`, `SNESMultiblockSetType()`,
872:           `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE`
873: M*/
874: PETSC_EXTERN PetscErrorCode SNESCreate_Multiblock(SNES snes)
875: {
876:   SNES_Multiblock *mb;

878:   PetscFunctionBegin;
879:   snes->ops->destroy        = SNESDestroy_Multiblock;
880:   snes->ops->setup          = SNESSetUp_Multiblock;
881:   snes->ops->setfromoptions = SNESSetFromOptions_Multiblock;
882:   snes->ops->view           = SNESView_Multiblock;
883:   snes->ops->solve          = SNESSolve_Multiblock;
884:   snes->ops->reset          = SNESReset_Multiblock;

886:   snes->usesksp = PETSC_FALSE;
887:   snes->usesnpc = PETSC_FALSE;

889:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

891:   PetscCall(PetscNew(&mb));
892:   snes->data    = (void *)mb;
893:   mb->defined   = PETSC_FALSE;
894:   mb->numBlocks = 0;
895:   mb->bs        = -1;
896:   mb->type      = PC_COMPOSITE_MULTIPLICATIVE;

898:   /* We attach functions so that they can be called on another PC without crashing the program */
899:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockSetFields_C", SNESMultiblockSetFields_Default));
900:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockSetIS_C", SNESMultiblockSetIS_Default));
901:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockSetType_C", SNESMultiblockSetType_Default));
902:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockSetBlockSize_C", SNESMultiblockSetBlockSize_Default));
903:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMultiblockGetSubSNES_C", SNESMultiblockGetSubSNES_Default));
904:   PetscFunctionReturn(PETSC_SUCCESS);
905: }