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