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