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 = PetscTypeCompare((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 { 150 mb->bs = 1; 151 } 152 } 153 154 ierr = PetscOptionsGetBool(((PetscObject) snes)->prefix, "-snes_multiblock_default", &flg, PETSC_NULL);CHKERRQ(ierr); 155 ierr = PetscOptionsGetBool(((PetscObject) snes)->prefix, "-snes_multiblock_detect_saddle_point", &stokes, PETSC_NULL);CHKERRQ(ierr); 156 if (stokes) { 157 IS zerodiags, rest; 158 PetscInt nmin, nmax; 159 160 ierr = MatGetOwnershipRange(snes->jacobian_pre, &nmin, &nmax);CHKERRQ(ierr); 161 ierr = MatFindZeroDiagonals(snes->jacobian_pre, &zerodiags);CHKERRQ(ierr); 162 ierr = ISComplement(zerodiags, nmin, nmax, &rest);CHKERRQ(ierr); 163 ierr = SNESMultiblockSetIS(snes, "0", rest);CHKERRQ(ierr); 164 ierr = SNESMultiblockSetIS(snes, "1", zerodiags);CHKERRQ(ierr); 165 ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 166 ierr = ISDestroy(&rest);CHKERRQ(ierr); 167 } else { 168 if (!flg) { 169 /* Allow user to set fields from command line, if bs was known at the time of SNESSetFromOptions_Multiblock() 170 then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 171 ierr = SNESMultiblockSetFieldsRuntime_Private(snes);CHKERRQ(ierr); 172 if (mb->defined) {ierr = PetscInfo(snes, "Blocks defined using the options database\n");CHKERRQ(ierr);} 173 } 174 if (flg || !mb->defined) { 175 ierr = PetscInfo(snes, "Using default splitting of fields\n");CHKERRQ(ierr); 176 for(i = 0; i < mb->bs; ++i) { 177 char name[8]; 178 179 ierr = PetscSNPrintf(name, sizeof name, "%D", i);CHKERRQ(ierr); 180 ierr = SNESMultiblockSetFields(snes, name, 1, &i);CHKERRQ(ierr); 181 } 182 mb->defaultblocks = PETSC_TRUE; 183 } 184 } 185 } 186 } else if (mb->numBlocks == 1) { 187 if (blocks->is) { 188 IS is2; 189 PetscInt nmin, nmax; 190 191 ierr = MatGetOwnershipRange(snes->jacobian_pre, &nmin, &nmax);CHKERRQ(ierr); 192 ierr = ISComplement(blocks->is, nmin, nmax, &is2);CHKERRQ(ierr); 193 ierr = SNESMultiblockSetIS(snes, "1", is2);CHKERRQ(ierr); 194 ierr = ISDestroy(&is2);CHKERRQ(ierr); 195 } else { 196 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must provide at least two sets of fields to SNES multiblock"); 197 } 198 } 199 if (mb->numBlocks < 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unhandled case, must have at least two blocks"); 200 PetscFunctionReturn(0); 201 } 202 203 /* 204 SNESSetUp_Multiblock - Sets up the internal data structures for the later use 205 of the SNESMULTIBLOCK nonlinear solver. 206 207 Input Parameters: 208 + snes - the SNES context 209 - x - the solution vector 210 211 Application Interface Routine: SNESSetUp() 212 */ 213 #undef __FUNCT__ 214 #define __FUNCT__ "SNESSetUp_Multiblock" 215 PetscErrorCode SNESSetUp_Multiblock(SNES snes) 216 { 217 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 218 BlockDesc blocks; 219 PetscInt i, numBlocks; 220 PetscErrorCode ierr; 221 222 PetscFunctionBegin; 223 /* ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); */ 224 ierr = SNESMultiblockSetDefaults(snes);CHKERRQ(ierr); 225 numBlocks = mb->numBlocks; 226 blocks = mb->blocks; 227 228 /* Create ISs */ 229 if (!mb->issetup) { 230 PetscInt ccsize, rstart, rend, nslots, bs; 231 PetscBool sorted; 232 233 mb->issetup = PETSC_TRUE; 234 bs = mb->bs; 235 ierr = MatGetOwnershipRange(snes->jacobian_pre, &rstart, &rend);CHKERRQ(ierr); 236 ierr = MatGetLocalSize(snes->jacobian_pre, PETSC_NULL, &ccsize);CHKERRQ(ierr); 237 nslots = (rend - rstart)/bs; 238 for(i = 0; i < numBlocks; ++i) { 239 if (mb->defaultblocks) { 240 ierr = ISCreateStride(((PetscObject) snes)->comm, nslots, rstart+i, numBlocks, &blocks->is);CHKERRQ(ierr); 241 } else if (!blocks->is) { 242 if (blocks->nfields > 1) { 243 PetscInt *ii, j, k, nfields = blocks->nfields, *fields = blocks->fields; 244 245 ierr = PetscMalloc(nfields*nslots*sizeof(PetscInt), &ii);CHKERRQ(ierr); 246 for(j = 0; j < nslots; ++j) { 247 for(k = 0; k < nfields; ++k) { 248 ii[nfields*j + k] = rstart + bs*j + fields[k]; 249 } 250 } 251 ierr = ISCreateGeneral(((PetscObject) snes)->comm, nslots*nfields, ii, PETSC_OWN_POINTER, &blocks->is);CHKERRQ(ierr); 252 } else { 253 ierr = ISCreateStride(((PetscObject) snes)->comm, nslots, rstart+blocks->fields[0], bs, &blocks->is);CHKERRQ(ierr); 254 } 255 } 256 ierr = ISSorted(blocks->is, &sorted);CHKERRQ(ierr); 257 if (!sorted) {SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Fields must be sorted when creating split");} 258 blocks = blocks->next; 259 } 260 } 261 262 #if 0 263 /* Create matrices */ 264 ilink = jac->head; 265 if (!jac->pmat) { 266 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 267 for (i=0; i<nsplit; i++) { 268 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 269 ilink = ilink->next; 270 } 271 } else { 272 for (i=0; i<nsplit; i++) { 273 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 274 ilink = ilink->next; 275 } 276 } 277 if (jac->realdiagonal) { 278 ilink = jac->head; 279 if (!jac->mat) { 280 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 281 for (i=0; i<nsplit; i++) { 282 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 283 ilink = ilink->next; 284 } 285 } else { 286 for (i=0; i<nsplit; i++) { 287 if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 288 ilink = ilink->next; 289 } 290 } 291 } else { 292 jac->mat = jac->pmat; 293 } 294 #endif 295 296 #if 0 297 if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 298 /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 299 ilink = jac->head; 300 if (!jac->Afield) { 301 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 302 for (i=0; i<nsplit; i++) { 303 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 304 ilink = ilink->next; 305 } 306 } else { 307 for (i=0; i<nsplit; i++) { 308 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 309 ilink = ilink->next; 310 } 311 } 312 } 313 #endif 314 315 if (mb->type == PC_COMPOSITE_SCHUR) { 316 #if 0 317 IS ccis; 318 PetscInt rstart,rend; 319 if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 320 321 /* When extracting off-diagonal submatrices, we take complements from this range */ 322 ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 323 324 /* need to handle case when one is resetting up the preconditioner */ 325 if (jac->schur) { 326 ilink = jac->head; 327 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 328 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 329 ierr = ISDestroy(&ccis);CHKERRQ(ierr); 330 ilink = ilink->next; 331 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 332 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 333 ierr = ISDestroy(&ccis);CHKERRQ(ierr); 334 ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 335 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 336 337 } else { 338 KSP ksp; 339 char schurprefix[256]; 340 341 /* extract the A01 and A10 matrices */ 342 ilink = jac->head; 343 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 344 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 345 ierr = ISDestroy(&ccis);CHKERRQ(ierr); 346 ilink = ilink->next; 347 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 348 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 349 ierr = ISDestroy(&ccis);CHKERRQ(ierr); 350 /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 351 ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 352 /* set tabbing and options prefix of KSP inside the MatSchur */ 353 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 354 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 355 ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 356 ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 357 ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 358 359 ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 360 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 361 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 362 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 363 if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 364 PC pc; 365 ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 366 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 367 /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 368 } 369 ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 370 ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 371 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 372 ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 373 374 ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 375 ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 376 ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 377 ilink = jac->head; 378 ilink->x = jac->x[0]; ilink->y = jac->y[0]; 379 ilink = ilink->next; 380 ilink->x = jac->x[1]; ilink->y = jac->y[1]; 381 } 382 #endif 383 } else { 384 /* Set up the individual SNESs */ 385 blocks = mb->blocks; 386 i = 0; 387 while (blocks) { 388 /*TODO: Set these correctly */ 389 /*ierr = SNESSetFunction(blocks->snes, blocks->x, func);CHKERRQ(ierr);*/ 390 /*ierr = SNESSetJacobian(blocks->snes, blocks->x, jac);CHKERRQ(ierr);*/ 391 ierr = VecDuplicate(blocks->snes->vec_sol, &blocks->x);CHKERRQ(ierr); 392 /* really want setfromoptions called in SNESSetFromOptions_Multiblock(), but it is not ready yet */ 393 ierr = SNESSetFromOptions(blocks->snes);CHKERRQ(ierr); 394 ierr = SNESSetUp(blocks->snes);CHKERRQ(ierr); 395 blocks = blocks->next; 396 i++; 397 } 398 } 399 400 /* Compute scatter contexts needed by multiplicative versions and non-default splits */ 401 if (!mb->blocks->sctx) { 402 Vec xtmp; 403 404 blocks = mb->blocks; 405 ierr = MatGetVecs(snes->jacobian_pre, &xtmp, PETSC_NULL);CHKERRQ(ierr); 406 while(blocks) { 407 ierr = VecScatterCreate(xtmp, blocks->is, blocks->x, PETSC_NULL, &blocks->sctx);CHKERRQ(ierr); 408 blocks = blocks->next; 409 } 410 ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 411 } 412 PetscFunctionReturn(0); 413 } 414 415 /* 416 SNESSetFromOptions_Multiblock - Sets various parameters for the SNESMULTIBLOCK method. 417 418 Input Parameter: 419 . snes - the SNES context 420 421 Application Interface Routine: SNESSetFromOptions() 422 */ 423 #undef __FUNCT__ 424 #define __FUNCT__ "SNESSetFromOptions_Multiblock" 425 static PetscErrorCode SNESSetFromOptions_Multiblock(SNES snes) 426 { 427 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 428 PCCompositeType ctype; 429 PetscInt bs; 430 PetscBool flg; 431 PetscErrorCode ierr; 432 433 PetscFunctionBegin; 434 ierr = PetscOptionsHead("SNES Multiblock options");CHKERRQ(ierr); 435 ierr = PetscOptionsInt("-snes_multiblock_block_size", "Blocksize that defines number of fields", "PCFieldSplitSetBlockSize", mb->bs, &bs, &flg);CHKERRQ(ierr); 436 if (flg) {ierr = SNESMultiblockSetBlockSize(snes, bs);CHKERRQ(ierr);} 437 ierr = PetscOptionsEnum("-snes_multiblock_type", "Type of composition", "PCFieldSplitSetType", PCCompositeTypes, (PetscEnum) mb->type, (PetscEnum *) &ctype, &flg);CHKERRQ(ierr); 438 if (flg) { 439 ierr = SNESMultiblockSetType(snes,ctype);CHKERRQ(ierr); 440 } 441 /* Only setup fields once */ 442 if ((mb->bs > 0) && (mb->numBlocks == 0)) { 443 /* only allow user to set fields from command line if bs is already known, otherwise user can set them in SNESMultiblockSetDefaults() */ 444 ierr = SNESMultiblockSetFieldsRuntime_Private(snes);CHKERRQ(ierr); 445 if (mb->defined) {ierr = PetscInfo(snes, "Blocks defined using the options database\n");CHKERRQ(ierr);} 446 } 447 ierr = PetscOptionsTail();CHKERRQ(ierr); 448 PetscFunctionReturn(0); 449 } 450 451 /* 452 SNESView_Multiblock - Prints info from the SNESMULTIBLOCK data structure. 453 454 Input Parameters: 455 + SNES - the SNES context 456 - viewer - visualization context 457 458 Application Interface Routine: SNESView() 459 */ 460 #undef __FUNCT__ 461 #define __FUNCT__ "SNESView_Multiblock" 462 static PetscErrorCode SNESView_Multiblock(SNES snes, PetscViewer viewer) 463 { 464 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 465 BlockDesc blocks = mb->blocks; 466 PetscBool iascii; 467 PetscErrorCode ierr; 468 469 PetscFunctionBegin; 470 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 471 if (iascii) { 472 ierr = PetscViewerASCIIPrintf(viewer," Multiblock with %s composition: total blocks = %D, blocksize = %D\n", PCCompositeTypes[mb->type], mb->numBlocks, mb->bs);CHKERRQ(ierr); 473 ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following SNES objects:\n");CHKERRQ(ierr); 474 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 475 while (blocks) { 476 if (blocks->fields) { 477 PetscInt j; 478 479 ierr = PetscViewerASCIIPrintf(viewer, "Block %s Fields ", blocks->name);CHKERRQ(ierr); 480 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr); 481 for(j = 0; j < blocks->nfields; ++j) { 482 if (j > 0) { 483 ierr = PetscViewerASCIIPrintf(viewer, ",");CHKERRQ(ierr); 484 } 485 ierr = PetscViewerASCIIPrintf(viewer, " %D", blocks->fields[j]);CHKERRQ(ierr); 486 } 487 ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 488 ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr); 489 } else { 490 ierr = PetscViewerASCIIPrintf(viewer, "Block %s Defined by IS\n", blocks->name);CHKERRQ(ierr); 491 } 492 ierr = SNESView(blocks->snes, viewer);CHKERRQ(ierr); 493 blocks = blocks->next; 494 } 495 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 496 } 497 PetscFunctionReturn(0); 498 } 499 500 /* 501 SNESSolve_Multiblock - Solves a nonlinear system with the Multiblock method. 502 503 Input Parameters: 504 . snes - the SNES context 505 506 Output Parameter: 507 . outits - number of iterations until termination 508 509 Application Interface Routine: SNESSolve() 510 */ 511 #undef __FUNCT__ 512 #define __FUNCT__ "SNESSolve_Multiblock" 513 PetscErrorCode SNESSolve_Multiblock(SNES snes) 514 { 515 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 516 Vec X, Y, F; 517 PetscReal fnorm; 518 PetscInt maxits, i; 519 PetscErrorCode ierr; 520 521 PetscFunctionBegin; 522 snes->reason = SNES_CONVERGED_ITERATING; 523 524 maxits = snes->max_its; /* maximum number of iterations */ 525 X = snes->vec_sol; /* X^n */ 526 Y = snes->vec_sol_update; /* \tilde X */ 527 F = snes->vec_func; /* residual vector */ 528 529 ierr = VecSetBlockSize(X, mb->bs);CHKERRQ(ierr); 530 ierr = VecSetBlockSize(Y, mb->bs);CHKERRQ(ierr); 531 ierr = VecSetBlockSize(F, mb->bs);CHKERRQ(ierr); 532 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 533 snes->iter = 0; 534 snes->norm = 0.; 535 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 536 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 537 if (snes->domainerror) { 538 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 539 PetscFunctionReturn(0); 540 } 541 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 542 if (PetscIsInfOrNanReal(fnorm)) {SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in norm");} 543 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 544 snes->norm = fnorm; 545 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 546 SNESLogConvHistory(snes,fnorm,0); 547 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 548 549 /* set parameter for default relative tolerance convergence test */ 550 snes->ttol = fnorm*snes->rtol; 551 /* test convergence */ 552 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 553 if (snes->reason) PetscFunctionReturn(0); 554 555 for(i = 0; i < maxits; i++) { 556 /* Call general purpose update function */ 557 if (snes->ops->update) { 558 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 559 } 560 /* Compute X^{new} from subsolves */ 561 if (mb->type == PC_COMPOSITE_ADDITIVE) { 562 BlockDesc blocks = mb->blocks; 563 564 if (mb->defaultblocks) { 565 /*TODO: Make an array of Vecs for this */ 566 /*ierr = VecStrideGatherAll(X, mb->x, INSERT_VALUES);CHKERRQ(ierr);*/ 567 while (blocks) { 568 ierr = SNESSolve(blocks->snes, PETSC_NULL, blocks->x);CHKERRQ(ierr); 569 blocks = blocks->next; 570 } 571 /*ierr = VecStrideScatterAll(mb->x, X, INSERT_VALUES);CHKERRQ(ierr);*/ 572 } else { 573 while (blocks) { 574 ierr = VecScatterBegin(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 575 ierr = VecScatterEnd(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 576 ierr = SNESSolve(blocks->snes, PETSC_NULL, blocks->x);CHKERRQ(ierr); 577 ierr = VecScatterBegin(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 578 ierr = VecScatterEnd(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 579 blocks = blocks->next; 580 } 581 } 582 } else { 583 SETERRQ1(((PetscObject) snes)->comm, PETSC_ERR_SUP, "Unsupported or unknown composition", (int) mb->type); 584 } 585 CHKMEMQ; 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)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm"); 590 591 if (snes->nfuncs >= snes->max_funcs) { 592 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 593 break; 594 } 595 if (snes->domainerror) { 596 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 597 PetscFunctionReturn(0); 598 } 599 /* Monitor convergence */ 600 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 601 snes->iter = i+1; 602 snes->norm = fnorm; 603 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 604 SNESLogConvHistory(snes,snes->norm,0); 605 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 606 /* Test for convergence */ 607 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 608 if (snes->reason) break; 609 } 610 if (i == maxits) { 611 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 612 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 613 } 614 PetscFunctionReturn(0); 615 } 616 617 EXTERN_C_BEGIN 618 #undef __FUNCT__ 619 #define __FUNCT__ "SNESMultiblockSetFields_Default" 620 PetscErrorCode SNESMultiblockSetFields_Default(SNES snes, const char name[], PetscInt n, const PetscInt fields[]) 621 { 622 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 623 BlockDesc newblock, next = mb->blocks; 624 char prefix[128]; 625 PetscInt i; 626 PetscErrorCode ierr; 627 628 PetscFunctionBegin; 629 if (mb->defined) { 630 ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr); 631 PetscFunctionReturn(0); 632 } 633 for(i = 0; i < n; ++i) { 634 if (fields[i] >= mb->bs) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %D requested but only %D exist", fields[i], mb->bs);} 635 if (fields[i] < 0) {SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %D requested", fields[i]);} 636 } 637 ierr = PetscNew(struct _BlockDesc, &newblock);CHKERRQ(ierr); 638 if (name) { 639 ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr); 640 } else { 641 PetscInt len = floor(log10(mb->numBlocks))+1; 642 643 ierr = PetscMalloc((len+1)*sizeof(char), &newblock->name);CHKERRQ(ierr); 644 ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr); 645 } 646 newblock->nfields = n; 647 ierr = PetscMalloc(n*sizeof(PetscInt), &newblock->fields);CHKERRQ(ierr); 648 ierr = PetscMemcpy(newblock->fields, fields, n*sizeof(PetscInt));CHKERRQ(ierr); 649 newblock->next = PETSC_NULL; 650 ierr = SNESCreate(((PetscObject) snes)->comm, &newblock->snes);CHKERRQ(ierr); 651 ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr); 652 ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr); 653 ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr); 654 ierr = PetscSNPrintf(prefix, sizeof prefix, "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr); 655 ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr); 656 657 if (!next) { 658 mb->blocks = newblock; 659 newblock->previous = PETSC_NULL; 660 } else { 661 while (next->next) { 662 next = next->next; 663 } 664 next->next = newblock; 665 newblock->previous = next; 666 } 667 mb->numBlocks++; 668 PetscFunctionReturn(0); 669 } 670 EXTERN_C_END 671 672 EXTERN_C_BEGIN 673 #undef __FUNCT__ 674 #define __FUNCT__ "SNESMultiblockSetIS_Default" 675 PetscErrorCode SNESMultiblockSetIS_Default(SNES snes, const char name[], IS is) 676 { 677 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 678 BlockDesc newblock, next = mb->blocks; 679 char prefix[128]; 680 PetscErrorCode ierr; 681 682 PetscFunctionBegin; 683 if (mb->defined) { 684 ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr); 685 PetscFunctionReturn(0); 686 } 687 ierr = PetscNew(struct _BlockDesc, &newblock);CHKERRQ(ierr); 688 if (name) { 689 ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr); 690 } else { 691 PetscInt len = floor(log10(mb->numBlocks))+1; 692 693 ierr = PetscMalloc((len+1)*sizeof(char), &newblock->name);CHKERRQ(ierr); 694 ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr); 695 } 696 newblock->is = is; 697 ierr = PetscObjectReference((PetscObject) is);CHKERRQ(ierr); 698 newblock->next = PETSC_NULL; 699 ierr = SNESCreate(((PetscObject) snes)->comm, &newblock->snes);CHKERRQ(ierr); 700 ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr); 701 ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr); 702 ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr); 703 ierr = PetscSNPrintf(prefix, sizeof prefix, "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr); 704 ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr); 705 706 if (!next) { 707 mb->blocks = newblock; 708 newblock->previous = PETSC_NULL; 709 } else { 710 while (next->next) { 711 next = next->next; 712 } 713 next->next = newblock; 714 newblock->previous = next; 715 } 716 mb->numBlocks++; 717 PetscFunctionReturn(0); 718 } 719 EXTERN_C_END 720 721 EXTERN_C_BEGIN 722 #undef __FUNCT__ 723 #define __FUNCT__ "SNESMultiblockSetBlockSize_Default" 724 PetscErrorCode SNESMultiblockSetBlockSize_Default(SNES snes, PetscInt bs) 725 { 726 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 727 728 PetscFunctionBegin; 729 if (bs < 1) {SETERRQ1(((PetscObject) snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %D", bs);} 730 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);} 731 mb->bs = bs; 732 PetscFunctionReturn(0); 733 } 734 EXTERN_C_END 735 736 EXTERN_C_BEGIN 737 #undef __FUNCT__ 738 #define __FUNCT__ "SNESMultiblockGetSubSNES_Default" 739 PetscErrorCode SNESMultiblockGetSubSNES_Default(SNES snes, PetscInt *n, SNES **subsnes) 740 { 741 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 742 BlockDesc blocks = mb->blocks; 743 PetscInt cnt = 0; 744 PetscErrorCode ierr; 745 746 PetscFunctionBegin; 747 ierr = PetscMalloc(mb->numBlocks * sizeof(SNES), subsnes);CHKERRQ(ierr); 748 while (blocks) { 749 (*subsnes)[cnt++] = blocks->snes; 750 blocks = blocks->next; 751 } 752 if (cnt != mb->numBlocks) { 753 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); 754 } 755 if (n) {*n = mb->numBlocks;} 756 PetscFunctionReturn(0); 757 } 758 EXTERN_C_END 759 760 EXTERN_C_BEGIN 761 #undef __FUNCT__ 762 #define __FUNCT__ "SNESMultiblockSetType_Default" 763 PetscErrorCode SNESMultiblockSetType_Default(SNES snes, PCCompositeType type) 764 { 765 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 766 PetscErrorCode ierr; 767 768 PetscFunctionBegin; 769 mb->type = type; 770 if (type == PC_COMPOSITE_SCHUR) { 771 #if 1 772 SETERRQ(((PetscObject) snes)->comm, PETSC_ERR_SUP, "The Schur composite type is not yet supported"); 773 #else 774 snes->ops->solve = SNESSolve_Multiblock_Schur; 775 snes->ops->view = SNESView_Multiblock_Schur; 776 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Schur", SNESMultiblockGetSubSNES_Schur);CHKERRQ(ierr); 777 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", "SNESMultiblockSchurPrecondition_Default", SNESMultiblockSchurPrecondition_Default);CHKERRQ(ierr); 778 #endif 779 } else { 780 snes->ops->solve = SNESSolve_Multiblock; 781 snes->ops->view = SNESView_Multiblock; 782 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Default", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr); 783 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", "", 0);CHKERRQ(ierr); 784 } 785 PetscFunctionReturn(0); 786 } 787 EXTERN_C_END 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 PETSC_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(((PetscObject) snes)->comm, 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 PETSC_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(), SNESLS, SNESTR, SNESNRICHARDSON 961 M*/ 962 EXTERN_C_BEGIN 963 #undef __FUNCT__ 964 #define __FUNCT__ "SNESCreate_Multiblock" 965 PetscErrorCode SNESCreate_Multiblock(SNES snes) 966 { 967 SNES_Multiblock *mb; 968 PetscErrorCode ierr; 969 970 PetscFunctionBegin; 971 snes->ops->destroy = SNESDestroy_Multiblock; 972 snes->ops->setup = SNESSetUp_Multiblock; 973 snes->ops->setfromoptions = SNESSetFromOptions_Multiblock; 974 snes->ops->view = SNESView_Multiblock; 975 snes->ops->solve = SNESSolve_Multiblock; 976 snes->ops->reset = SNESReset_Multiblock; 977 978 snes->usesksp = PETSC_FALSE; 979 980 ierr = PetscNewLog(snes, SNES_Multiblock, &mb);CHKERRQ(ierr); 981 snes->data = (void*) mb; 982 mb->defined = PETSC_FALSE; 983 mb->numBlocks = 0; 984 mb->bs = -1; 985 mb->type = PC_COMPOSITE_MULTIPLICATIVE; 986 987 /* We attach functions so that they can be called on another PC without crashing the program */ 988 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetFields_C", "SNESMultiblockSetFields_Default", SNESMultiblockSetFields_Default);CHKERRQ(ierr); 989 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetIS_C", "SNESMultiblockSetIS_Default", SNESMultiblockSetIS_Default);CHKERRQ(ierr); 990 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetType_C", "SNESMultiblockSetType_Default", SNESMultiblockSetType_Default);CHKERRQ(ierr); 991 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetBlockSize_C", "SNESMultiblockSetBlockSize_Default", SNESMultiblockSetBlockSize_Default);CHKERRQ(ierr); 992 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Default", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr); 993 PetscFunctionReturn(0); 994 } 995 EXTERN_C_END 996