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