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 EXTERN_C_BEGIN 625 #undef __FUNCT__ 626 #define __FUNCT__ "SNESMultiblockSetFields_Default" 627 PetscErrorCode SNESMultiblockSetFields_Default(SNES snes, const char name[], PetscInt n, const PetscInt fields[]) 628 { 629 SNES_Multiblock *mb = (SNES_Multiblock*) snes->data; 630 BlockDesc newblock, next = mb->blocks; 631 char prefix[128]; 632 PetscInt i; 633 PetscErrorCode ierr; 634 635 PetscFunctionBegin; 636 if (mb->defined) { 637 ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr); 638 PetscFunctionReturn(0); 639 } 640 for (i = 0; i < n; ++i) { 641 if (fields[i] >= mb->bs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %D requested but only %D exist", fields[i], mb->bs); 642 if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %D requested", fields[i]); 643 } 644 ierr = PetscNew(struct _BlockDesc, &newblock);CHKERRQ(ierr); 645 if (name) { 646 ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr); 647 } else { 648 PetscInt len = floor(log10(mb->numBlocks))+1; 649 650 ierr = PetscMalloc((len+1)*sizeof(char), &newblock->name);CHKERRQ(ierr); 651 ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr); 652 } 653 newblock->nfields = n; 654 655 ierr = PetscMalloc(n*sizeof(PetscInt), &newblock->fields);CHKERRQ(ierr); 656 ierr = PetscMemcpy(newblock->fields, fields, n*sizeof(PetscInt));CHKERRQ(ierr); 657 658 newblock->next = NULL; 659 660 ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &newblock->snes);CHKERRQ(ierr); 661 ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr); 662 ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr); 663 ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr); 664 ierr = PetscSNPrintf(prefix, sizeof(prefix), "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr); 665 ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr); 666 667 if (!next) { 668 mb->blocks = newblock; 669 newblock->previous = NULL; 670 } else { 671 while (next->next) { 672 next = next->next; 673 } 674 next->next = newblock; 675 newblock->previous = next; 676 } 677 mb->numBlocks++; 678 PetscFunctionReturn(0); 679 } 680 EXTERN_C_END 681 682 EXTERN_C_BEGIN 683 #undef __FUNCT__ 684 #define __FUNCT__ "SNESMultiblockSetIS_Default" 685 PetscErrorCode SNESMultiblockSetIS_Default(SNES snes, const char name[], IS is) 686 { 687 SNES_Multiblock *mb = (SNES_Multiblock*) snes->data; 688 BlockDesc newblock, next = mb->blocks; 689 char prefix[128]; 690 PetscErrorCode ierr; 691 692 PetscFunctionBegin; 693 if (mb->defined) { 694 ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr); 695 PetscFunctionReturn(0); 696 } 697 ierr = PetscNew(struct _BlockDesc, &newblock);CHKERRQ(ierr); 698 if (name) { 699 ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr); 700 } else { 701 PetscInt len = floor(log10(mb->numBlocks))+1; 702 703 ierr = PetscMalloc((len+1)*sizeof(char), &newblock->name);CHKERRQ(ierr); 704 ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr); 705 } 706 newblock->is = is; 707 708 ierr = PetscObjectReference((PetscObject) is);CHKERRQ(ierr); 709 710 newblock->next = NULL; 711 712 ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &newblock->snes);CHKERRQ(ierr); 713 ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr); 714 ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr); 715 ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr); 716 ierr = PetscSNPrintf(prefix, sizeof(prefix), "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr); 717 ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr); 718 719 if (!next) { 720 mb->blocks = newblock; 721 newblock->previous = NULL; 722 } else { 723 while (next->next) { 724 next = next->next; 725 } 726 next->next = newblock; 727 newblock->previous = next; 728 } 729 mb->numBlocks++; 730 PetscFunctionReturn(0); 731 } 732 EXTERN_C_END 733 734 EXTERN_C_BEGIN 735 #undef __FUNCT__ 736 #define __FUNCT__ "SNESMultiblockSetBlockSize_Default" 737 PetscErrorCode SNESMultiblockSetBlockSize_Default(SNES snes, PetscInt bs) 738 { 739 SNES_Multiblock *mb = (SNES_Multiblock*) snes->data; 740 741 PetscFunctionBegin; 742 if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %D", bs); 743 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); 744 mb->bs = bs; 745 PetscFunctionReturn(0); 746 } 747 EXTERN_C_END 748 749 EXTERN_C_BEGIN 750 #undef __FUNCT__ 751 #define __FUNCT__ "SNESMultiblockGetSubSNES_Default" 752 PetscErrorCode SNESMultiblockGetSubSNES_Default(SNES snes, PetscInt *n, SNES **subsnes) 753 { 754 SNES_Multiblock *mb = (SNES_Multiblock*) snes->data; 755 BlockDesc blocks = mb->blocks; 756 PetscInt cnt = 0; 757 PetscErrorCode ierr; 758 759 PetscFunctionBegin; 760 ierr = PetscMalloc(mb->numBlocks * sizeof(SNES), subsnes);CHKERRQ(ierr); 761 while (blocks) { 762 (*subsnes)[cnt++] = blocks->snes; 763 blocks = blocks->next; 764 } 765 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); 766 767 if (n) *n = mb->numBlocks; 768 PetscFunctionReturn(0); 769 } 770 EXTERN_C_END 771 772 EXTERN_C_BEGIN 773 #undef __FUNCT__ 774 #define __FUNCT__ "SNESMultiblockSetType_Default" 775 PetscErrorCode SNESMultiblockSetType_Default(SNES snes, PCCompositeType type) 776 { 777 SNES_Multiblock *mb = (SNES_Multiblock*) snes->data; 778 PetscErrorCode ierr; 779 780 PetscFunctionBegin; 781 mb->type = type; 782 if (type == PC_COMPOSITE_SCHUR) { 783 #if 1 784 SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "The Schur composite type is not yet supported"); 785 #else 786 snes->ops->solve = SNESSolve_Multiblock_Schur; 787 snes->ops->view = SNESView_Multiblock_Schur; 788 789 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Schur", SNESMultiblockGetSubSNES_Schur);CHKERRQ(ierr); 790 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", "SNESMultiblockSchurPrecondition_Default", SNESMultiblockSchurPrecondition_Default);CHKERRQ(ierr); 791 #endif 792 } else { 793 snes->ops->solve = SNESSolve_Multiblock; 794 snes->ops->view = SNESView_Multiblock; 795 796 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Default", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr); 797 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", "", 0);CHKERRQ(ierr); 798 } 799 PetscFunctionReturn(0); 800 } 801 EXTERN_C_END 802 803 #undef __FUNCT__ 804 #define __FUNCT__ "SNESMultiblockSetFields" 805 /*@ 806 SNESMultiblockSetFields - Sets the fields for one particular block in the solver 807 808 Logically Collective on SNES 809 810 Input Parameters: 811 + snes - the solver 812 . name - name of this block, if NULL the number of the block is used 813 . n - the number of fields in this block 814 - fields - the fields in this block 815 816 Level: intermediate 817 818 Notes: Use SNESMultiblockSetIS() to set a completely general set of row indices as a block. 819 820 The SNESMultiblockSetFields() is for defining blocks as a group of strided indices, or fields. 821 For example, if the vector block size is three then one can define a block as field 0, or 822 1 or 2, or field 0,1 or 0,2 or 1,2 which means 823 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 824 where the numbered entries indicate what is in the block. 825 826 This function is called once per block (it creates a new block each time). Solve options 827 for this block will be available under the prefix -multiblock_BLOCKNAME_. 828 829 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetBlockSize(), SNESMultiblockSetIS() 830 @*/ 831 PetscErrorCode SNESMultiblockSetFields(SNES snes, const char name[], PetscInt n, const PetscInt *fields) 832 { 833 PetscErrorCode ierr; 834 835 PetscFunctionBegin; 836 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 837 PetscValidCharPointer(name, 2); 838 if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Provided number of fields %D in split \"%s\" not positive", n, name); 839 PetscValidIntPointer(fields, 3); 840 ierr = PetscTryMethod(snes, "SNESMultiblockSetFields_C", (SNES, const char[], PetscInt, const PetscInt*), (snes, name, n, fields));CHKERRQ(ierr); 841 PetscFunctionReturn(0); 842 } 843 844 #undef __FUNCT__ 845 #define __FUNCT__ "SNESMultiblockSetIS" 846 /*@ 847 SNESMultiblockSetIS - Sets the global row indices for the block 848 849 Logically Collective on SNES 850 851 Input Parameters: 852 + snes - the solver context 853 . name - name of this block, if NULL the number of the block is used 854 - is - the index set that defines the global row indices in this block 855 856 Notes: 857 Use SNESMultiblockSetFields(), for blocks defined by strides. 858 859 This function is called once per block (it creates a new block each time). Solve options 860 for this block will be available under the prefix -multiblock_BLOCKNAME_. 861 862 Level: intermediate 863 864 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetBlockSize() 865 @*/ 866 PetscErrorCode SNESMultiblockSetIS(SNES snes, const char name[], IS is) 867 { 868 PetscErrorCode ierr; 869 870 PetscFunctionBegin; 871 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 872 PetscValidCharPointer(name, 2); 873 PetscValidHeaderSpecific(is, IS_CLASSID, 3); 874 ierr = PetscTryMethod(snes, "SNESMultiblockSetIS_C", (SNES, const char[], IS), (snes, name, is));CHKERRQ(ierr); 875 PetscFunctionReturn(0); 876 } 877 878 #undef __FUNCT__ 879 #define __FUNCT__ "SNESMultiblockSetType" 880 /*@ 881 SNESMultiblockSetType - Sets the type of block combination. 882 883 Collective on SNES 884 885 Input Parameters: 886 + snes - the solver context 887 - type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 888 889 Options Database Key: 890 . -snes_multiblock_type <type: one of multiplicative, additive, symmetric_multiplicative> - Sets block combination type 891 892 Level: Developer 893 894 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 895 .seealso: PCCompositeSetType() 896 @*/ 897 PetscErrorCode SNESMultiblockSetType(SNES snes, PCCompositeType type) 898 { 899 PetscErrorCode ierr; 900 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 903 ierr = PetscTryMethod(snes, "SNESMultiblockSetType_C", (SNES, PCCompositeType), (snes, type));CHKERRQ(ierr); 904 PetscFunctionReturn(0); 905 } 906 907 #undef __FUNCT__ 908 #define __FUNCT__ "SNESMultiblockSetBlockSize" 909 /*@ 910 SNESMultiblockSetBlockSize - Sets the block size for structured mesh block division. If not set the matrix block size is used. 911 912 Logically Collective on SNES 913 914 Input Parameters: 915 + snes - the solver context 916 - bs - the block size 917 918 Level: intermediate 919 920 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetFields() 921 @*/ 922 PetscErrorCode SNESMultiblockSetBlockSize(SNES snes, PetscInt bs) 923 { 924 PetscErrorCode ierr; 925 926 PetscFunctionBegin; 927 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 928 PetscValidLogicalCollectiveInt(snes, bs, 2); 929 ierr = PetscTryMethod(snes, "SNESMultiblockSetBlockSize_C", (SNES, PetscInt), (snes,bs));CHKERRQ(ierr); 930 PetscFunctionReturn(0); 931 } 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "SNESMultiblockGetSubSNES" 935 /*@C 936 SNESMultiblockGetSubSNES - Gets the SNES contexts for all blocks 937 938 Collective on SNES 939 940 Input Parameter: 941 . snes - the solver context 942 943 Output Parameters: 944 + n - the number of blocks 945 - subsnes - the array of SNES contexts 946 947 Note: 948 After SNESMultiblockGetSubSNES() the array of SNESs MUST be freed by the user 949 (not each SNES, just the array that contains them). 950 951 You must call SNESSetUp() before calling SNESMultiblockGetSubSNES(). 952 953 Level: advanced 954 955 .seealso: SNESMULTIBLOCK 956 @*/ 957 PetscErrorCode SNESMultiblockGetSubSNES(SNES snes, PetscInt *n, SNES *subsnes[]) 958 { 959 PetscErrorCode ierr; 960 961 PetscFunctionBegin; 962 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 963 if (n) PetscValidIntPointer(n, 2); 964 ierr = PetscUseMethod(snes, "SNESMultiblockGetSubSNES_C", (SNES, PetscInt*, SNES **), (snes, n, subsnes));CHKERRQ(ierr); 965 PetscFunctionReturn(0); 966 } 967 968 /*MC 969 SNESMULTIBLOCK - Multiblock nonlinear solver that can use overlapping or nonoverlapping blocks, organized 970 additively (Jacobi) or multiplicatively (Gauss-Seidel). 971 972 Level: beginner 973 974 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNRICHARDSON 975 M*/ 976 EXTERN_C_BEGIN 977 #undef __FUNCT__ 978 #define __FUNCT__ "SNESCreate_Multiblock" 979 PetscErrorCode SNESCreate_Multiblock(SNES snes) 980 { 981 SNES_Multiblock *mb; 982 PetscErrorCode ierr; 983 984 PetscFunctionBegin; 985 snes->ops->destroy = SNESDestroy_Multiblock; 986 snes->ops->setup = SNESSetUp_Multiblock; 987 snes->ops->setfromoptions = SNESSetFromOptions_Multiblock; 988 snes->ops->view = SNESView_Multiblock; 989 snes->ops->solve = SNESSolve_Multiblock; 990 snes->ops->reset = SNESReset_Multiblock; 991 992 snes->usesksp = PETSC_FALSE; 993 994 ierr = PetscNewLog(snes, SNES_Multiblock, &mb);CHKERRQ(ierr); 995 snes->data = (void*) mb; 996 mb->defined = PETSC_FALSE; 997 mb->numBlocks = 0; 998 mb->bs = -1; 999 mb->type = PC_COMPOSITE_MULTIPLICATIVE; 1000 1001 /* We attach functions so that they can be called on another PC without crashing the program */ 1002 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetFields_C", "SNESMultiblockSetFields_Default", SNESMultiblockSetFields_Default);CHKERRQ(ierr); 1003 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetIS_C", "SNESMultiblockSetIS_Default", SNESMultiblockSetIS_Default);CHKERRQ(ierr); 1004 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetType_C", "SNESMultiblockSetType_Default", SNESMultiblockSetType_Default);CHKERRQ(ierr); 1005 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockSetBlockSize_C", "SNESMultiblockSetBlockSize_Default", SNESMultiblockSetBlockSize_Default);CHKERRQ(ierr); 1006 ierr = PetscObjectComposeFunction((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Default", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr); 1007 PetscFunctionReturn(0); 1008 } 1009 EXTERN_C_END 1010