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