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 { 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 = PetscObjectTypeCompare((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 537 if (!snes->vec_func_init_set){ 538 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 539 if (snes->domainerror) { 540 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 541 PetscFunctionReturn(0); 542 } 543 } else { 544 snes->vec_func_init_set = PETSC_FALSE; 545 } 546 if (!snes->norm_init_set) { 547 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 548 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in norm"); 549 } else { 550 fnorm = snes->norm_init; 551 snes->norm_init_set = PETSC_FALSE; 552 } 553 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 554 snes->norm = fnorm; 555 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 556 SNESLogConvHistory(snes,fnorm,0); 557 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 558 559 /* set parameter for default relative tolerance convergence test */ 560 snes->ttol = fnorm*snes->rtol; 561 /* test convergence */ 562 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 563 if (snes->reason) PetscFunctionReturn(0); 564 565 for(i = 0; i < maxits; i++) { 566 /* Call general purpose update function */ 567 if (snes->ops->update) { 568 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 569 } 570 /* Compute X^{new} from subsolves */ 571 if (mb->type == PC_COMPOSITE_ADDITIVE) { 572 BlockDesc blocks = mb->blocks; 573 574 if (mb->defaultblocks) { 575 /*TODO: Make an array of Vecs for this */ 576 /*ierr = VecStrideGatherAll(X, mb->x, INSERT_VALUES);CHKERRQ(ierr);*/ 577 while (blocks) { 578 ierr = SNESSolve(blocks->snes, PETSC_NULL, blocks->x);CHKERRQ(ierr); 579 blocks = blocks->next; 580 } 581 /*ierr = VecStrideScatterAll(mb->x, X, INSERT_VALUES);CHKERRQ(ierr);*/ 582 } else { 583 while (blocks) { 584 ierr = VecScatterBegin(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 585 ierr = VecScatterEnd(blocks->sctx, X, blocks->x, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 586 ierr = SNESSolve(blocks->snes, PETSC_NULL, blocks->x);CHKERRQ(ierr); 587 ierr = VecScatterBegin(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 588 ierr = VecScatterEnd(blocks->sctx, blocks->x, X, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 589 blocks = blocks->next; 590 } 591 } 592 } else { 593 SETERRQ1(((PetscObject) snes)->comm, PETSC_ERR_SUP, "Unsupported or unknown composition", (int) mb->type); 594 } 595 CHKMEMQ; 596 /* Compute F(X^{new}) */ 597 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 598 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 599 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm"); 600 601 if (snes->nfuncs >= snes->max_funcs) { 602 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 603 break; 604 } 605 if (snes->domainerror) { 606 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 607 PetscFunctionReturn(0); 608 } 609 /* Monitor convergence */ 610 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 611 snes->iter = i+1; 612 snes->norm = fnorm; 613 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 614 SNESLogConvHistory(snes,snes->norm,0); 615 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 616 /* Test for convergence */ 617 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 618 if (snes->reason) break; 619 } 620 if (i == maxits) { 621 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 622 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 623 } 624 PetscFunctionReturn(0); 625 } 626 627 EXTERN_C_BEGIN 628 #undef __FUNCT__ 629 #define __FUNCT__ "SNESMultiblockSetFields_Default" 630 PetscErrorCode SNESMultiblockSetFields_Default(SNES snes, const char name[], PetscInt n, const PetscInt fields[]) 631 { 632 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 633 BlockDesc newblock, next = mb->blocks; 634 char prefix[128]; 635 PetscInt i; 636 PetscErrorCode ierr; 637 638 PetscFunctionBegin; 639 if (mb->defined) { 640 ierr = PetscInfo1(snes, "Ignoring new block \"%s\" because the blocks have already been defined\n", name);CHKERRQ(ierr); 641 PetscFunctionReturn(0); 642 } 643 for(i = 0; i < n; ++i) { 644 if (fields[i] >= mb->bs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %D requested but only %D exist", fields[i], mb->bs); 645 if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative field %D requested", fields[i]); 646 } 647 ierr = PetscNew(struct _BlockDesc, &newblock);CHKERRQ(ierr); 648 if (name) { 649 ierr = PetscStrallocpy(name, &newblock->name);CHKERRQ(ierr); 650 } else { 651 PetscInt len = floor(log10(mb->numBlocks))+1; 652 653 ierr = PetscMalloc((len+1)*sizeof(char), &newblock->name);CHKERRQ(ierr); 654 ierr = PetscSNPrintf(newblock->name, len, "%s", mb->numBlocks);CHKERRQ(ierr); 655 } 656 newblock->nfields = n; 657 ierr = PetscMalloc(n*sizeof(PetscInt), &newblock->fields);CHKERRQ(ierr); 658 ierr = PetscMemcpy(newblock->fields, fields, n*sizeof(PetscInt));CHKERRQ(ierr); 659 newblock->next = PETSC_NULL; 660 ierr = SNESCreate(((PetscObject) snes)->comm, &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 = PETSC_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 ierr = PetscObjectReference((PetscObject) is);CHKERRQ(ierr); 708 newblock->next = PETSC_NULL; 709 ierr = SNESCreate(((PetscObject) snes)->comm, &newblock->snes);CHKERRQ(ierr); 710 ierr = PetscObjectIncrementTabLevel((PetscObject) newblock->snes, (PetscObject) snes, 1);CHKERRQ(ierr); 711 ierr = SNESSetType(newblock->snes, SNESNRICHARDSON);CHKERRQ(ierr); 712 ierr = PetscLogObjectParent((PetscObject) snes, (PetscObject) newblock->snes);CHKERRQ(ierr); 713 ierr = PetscSNPrintf(prefix, sizeof prefix, "%smultiblock_%s_", ((PetscObject) snes)->prefix ? ((PetscObject) snes)->prefix : "", newblock->name);CHKERRQ(ierr); 714 ierr = SNESSetOptionsPrefix(newblock->snes, prefix);CHKERRQ(ierr); 715 716 if (!next) { 717 mb->blocks = newblock; 718 newblock->previous = PETSC_NULL; 719 } else { 720 while (next->next) { 721 next = next->next; 722 } 723 next->next = newblock; 724 newblock->previous = next; 725 } 726 mb->numBlocks++; 727 PetscFunctionReturn(0); 728 } 729 EXTERN_C_END 730 731 EXTERN_C_BEGIN 732 #undef __FUNCT__ 733 #define __FUNCT__ "SNESMultiblockSetBlockSize_Default" 734 PetscErrorCode SNESMultiblockSetBlockSize_Default(SNES snes, PetscInt bs) 735 { 736 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 737 738 PetscFunctionBegin; 739 if (bs < 1) SETERRQ1(((PetscObject) snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Blocksize must be positive, you gave %D", bs); 740 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); 741 mb->bs = bs; 742 PetscFunctionReturn(0); 743 } 744 EXTERN_C_END 745 746 EXTERN_C_BEGIN 747 #undef __FUNCT__ 748 #define __FUNCT__ "SNESMultiblockGetSubSNES_Default" 749 PetscErrorCode SNESMultiblockGetSubSNES_Default(SNES snes, PetscInt *n, SNES **subsnes) 750 { 751 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 752 BlockDesc blocks = mb->blocks; 753 PetscInt cnt = 0; 754 PetscErrorCode ierr; 755 756 PetscFunctionBegin; 757 ierr = PetscMalloc(mb->numBlocks * sizeof(SNES), subsnes);CHKERRQ(ierr); 758 while (blocks) { 759 (*subsnes)[cnt++] = blocks->snes; 760 blocks = blocks->next; 761 } 762 if (cnt != mb->numBlocks) { 763 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); 764 } 765 if (n) {*n = mb->numBlocks;} 766 PetscFunctionReturn(0); 767 } 768 EXTERN_C_END 769 770 EXTERN_C_BEGIN 771 #undef __FUNCT__ 772 #define __FUNCT__ "SNESMultiblockSetType_Default" 773 PetscErrorCode SNESMultiblockSetType_Default(SNES snes, PCCompositeType type) 774 { 775 SNES_Multiblock *mb = (SNES_Multiblock *) snes->data; 776 PetscErrorCode ierr; 777 778 PetscFunctionBegin; 779 mb->type = type; 780 if (type == PC_COMPOSITE_SCHUR) { 781 #if 1 782 SETERRQ(((PetscObject) snes)->comm, PETSC_ERR_SUP, "The Schur composite type is not yet supported"); 783 #else 784 snes->ops->solve = SNESSolve_Multiblock_Schur; 785 snes->ops->view = SNESView_Multiblock_Schur; 786 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Schur", SNESMultiblockGetSubSNES_Schur);CHKERRQ(ierr); 787 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", "SNESMultiblockSchurPrecondition_Default", SNESMultiblockSchurPrecondition_Default);CHKERRQ(ierr); 788 #endif 789 } else { 790 snes->ops->solve = SNESSolve_Multiblock; 791 snes->ops->view = SNESView_Multiblock; 792 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Default", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr); 793 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSchurPrecondition_C", "", 0);CHKERRQ(ierr); 794 } 795 PetscFunctionReturn(0); 796 } 797 EXTERN_C_END 798 799 #undef __FUNCT__ 800 #define __FUNCT__ "SNESMultiblockSetFields" 801 /*@ 802 SNESMultiblockSetFields - Sets the fields for one particular block in the solver 803 804 Logically Collective on SNES 805 806 Input Parameters: 807 + snes - the solver 808 . name - name of this block, if PETSC_NULL the number of the block is used 809 . n - the number of fields in this block 810 - fields - the fields in this block 811 812 Level: intermediate 813 814 Notes: Use SNESMultiblockSetIS() to set a completely general set of row indices as a block. 815 816 The SNESMultiblockSetFields() is for defining blocks as a group of strided indices, or fields. 817 For example, if the vector block size is three then one can define a block as field 0, or 818 1 or 2, or field 0,1 or 0,2 or 1,2 which means 819 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 820 where the numbered entries indicate what is in the block. 821 822 This function is called once per block (it creates a new block each time). Solve options 823 for this block will be available under the prefix -multiblock_BLOCKNAME_. 824 825 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetBlockSize(), SNESMultiblockSetIS() 826 @*/ 827 PetscErrorCode SNESMultiblockSetFields(SNES snes, const char name[], PetscInt n, const PetscInt *fields) 828 { 829 PetscErrorCode ierr; 830 831 PetscFunctionBegin; 832 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 833 PetscValidCharPointer(name, 2); 834 if (n < 1) SETERRQ2(((PetscObject) snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Provided number of fields %D in split \"%s\" not positive", n, name); 835 PetscValidIntPointer(fields, 3); 836 ierr = PetscTryMethod(snes, "SNESMultiblockSetFields_C", (SNES, const char[], PetscInt, const PetscInt *), (snes, name, n, fields));CHKERRQ(ierr); 837 PetscFunctionReturn(0); 838 } 839 840 #undef __FUNCT__ 841 #define __FUNCT__ "SNESMultiblockSetIS" 842 /*@ 843 SNESMultiblockSetIS - Sets the global row indices for the block 844 845 Logically Collective on SNES 846 847 Input Parameters: 848 + snes - the solver context 849 . name - name of this block, if PETSC_NULL the number of the block is used 850 - is - the index set that defines the global row indices in this block 851 852 Notes: 853 Use SNESMultiblockSetFields(), for blocks defined by strides. 854 855 This function is called once per block (it creates a new block each time). Solve options 856 for this block will be available under the prefix -multiblock_BLOCKNAME_. 857 858 Level: intermediate 859 860 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetBlockSize() 861 @*/ 862 PetscErrorCode SNESMultiblockSetIS(SNES snes, const char name[], IS is) 863 { 864 PetscErrorCode ierr; 865 866 PetscFunctionBegin; 867 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 868 PetscValidCharPointer(name, 2); 869 PetscValidHeaderSpecific(is, IS_CLASSID, 3); 870 ierr = PetscTryMethod(snes, "SNESMultiblockSetIS_C", (SNES, const char[], IS), (snes, name, is));CHKERRQ(ierr); 871 PetscFunctionReturn(0); 872 } 873 874 #undef __FUNCT__ 875 #define __FUNCT__ "SNESMultiblockSetType" 876 /*@ 877 SNESMultiblockSetType - Sets the type of block combination. 878 879 Collective on SNES 880 881 Input Parameters: 882 + snes - the solver context 883 - type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 884 885 Options Database Key: 886 . -snes_multiblock_type <type: one of multiplicative, additive, symmetric_multiplicative> - Sets block combination type 887 888 Level: Developer 889 890 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 891 .seealso: PCCompositeSetType() 892 @*/ 893 PetscErrorCode SNESMultiblockSetType(SNES snes, PCCompositeType type) 894 { 895 PetscErrorCode ierr; 896 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 899 ierr = PetscTryMethod(snes, "SNESMultiblockSetType_C", (SNES, PCCompositeType), (snes, type));CHKERRQ(ierr); 900 PetscFunctionReturn(0); 901 } 902 903 #undef __FUNCT__ 904 #define __FUNCT__ "SNESMultiblockSetBlockSize" 905 /*@ 906 SNESMultiblockSetBlockSize - Sets the block size for structured mesh block division. If not set the matrix block size is used. 907 908 Logically Collective on SNES 909 910 Input Parameters: 911 + snes - the solver context 912 - bs - the block size 913 914 Level: intermediate 915 916 .seealso: SNESMultiblockGetSubSNES(), SNESMULTIBLOCK, SNESMultiblockSetFields() 917 @*/ 918 PetscErrorCode SNESMultiblockSetBlockSize(SNES snes, PetscInt bs) 919 { 920 PetscErrorCode ierr; 921 922 PetscFunctionBegin; 923 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 924 PetscValidLogicalCollectiveInt(snes, bs, 2); 925 ierr = PetscTryMethod(snes, "SNESMultiblockSetBlockSize_C", (SNES, PetscInt), (snes,bs));CHKERRQ(ierr); 926 PetscFunctionReturn(0); 927 } 928 929 #undef __FUNCT__ 930 #define __FUNCT__ "SNESMultiblockGetSubSNES" 931 /*@C 932 SNESMultiblockGetSubSNES - Gets the SNES contexts for all blocks 933 934 Collective on SNES 935 936 Input Parameter: 937 . snes - the solver context 938 939 Output Parameters: 940 + n - the number of blocks 941 - subsnes - the array of SNES contexts 942 943 Note: 944 After SNESMultiblockGetSubSNES() the array of SNESs MUST be freed by the user 945 (not each SNES, just the array that contains them). 946 947 You must call SNESSetUp() before calling SNESMultiblockGetSubSNES(). 948 949 Level: advanced 950 951 .seealso: SNESMULTIBLOCK 952 @*/ 953 PetscErrorCode SNESMultiblockGetSubSNES(SNES snes, PetscInt *n, SNES *subsnes[]) 954 { 955 PetscErrorCode ierr; 956 957 PetscFunctionBegin; 958 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 959 if (n) PetscValidIntPointer(n, 2); 960 ierr = PetscUseMethod(snes, "SNESMultiblockGetSubSNES_C", (SNES, PetscInt*, SNES **), (snes, n, subsnes));CHKERRQ(ierr); 961 PetscFunctionReturn(0); 962 } 963 964 /*MC 965 SNESMULTIBLOCK - Multiblock nonlinear solver that can use overlapping or nonoverlapping blocks, organized 966 additively (Jacobi) or multiplicatively (Gauss-Seidel). 967 968 Level: beginner 969 970 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNRICHARDSON 971 M*/ 972 EXTERN_C_BEGIN 973 #undef __FUNCT__ 974 #define __FUNCT__ "SNESCreate_Multiblock" 975 PetscErrorCode SNESCreate_Multiblock(SNES snes) 976 { 977 SNES_Multiblock *mb; 978 PetscErrorCode ierr; 979 980 PetscFunctionBegin; 981 snes->ops->destroy = SNESDestroy_Multiblock; 982 snes->ops->setup = SNESSetUp_Multiblock; 983 snes->ops->setfromoptions = SNESSetFromOptions_Multiblock; 984 snes->ops->view = SNESView_Multiblock; 985 snes->ops->solve = SNESSolve_Multiblock; 986 snes->ops->reset = SNESReset_Multiblock; 987 988 snes->usesksp = PETSC_FALSE; 989 990 ierr = PetscNewLog(snes, SNES_Multiblock, &mb);CHKERRQ(ierr); 991 snes->data = (void*) mb; 992 mb->defined = PETSC_FALSE; 993 mb->numBlocks = 0; 994 mb->bs = -1; 995 mb->type = PC_COMPOSITE_MULTIPLICATIVE; 996 997 /* We attach functions so that they can be called on another PC without crashing the program */ 998 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetFields_C", "SNESMultiblockSetFields_Default", SNESMultiblockSetFields_Default);CHKERRQ(ierr); 999 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetIS_C", "SNESMultiblockSetIS_Default", SNESMultiblockSetIS_Default);CHKERRQ(ierr); 1000 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetType_C", "SNESMultiblockSetType_Default", SNESMultiblockSetType_Default);CHKERRQ(ierr); 1001 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockSetBlockSize_C", "SNESMultiblockSetBlockSize_Default", SNESMultiblockSetBlockSize_Default);CHKERRQ(ierr); 1002 ierr = PetscObjectComposeFunctionDynamic((PetscObject) snes, "SNESMultiblockGetSubSNES_C", "SNESMultiblockGetSubSNES_Default", SNESMultiblockGetSubSNES_Default);CHKERRQ(ierr); 1003 PetscFunctionReturn(0); 1004 } 1005 EXTERN_C_END 1006