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