1 #define PETSCKSP_DLL 2 3 /* 4 5 */ 6 #include "private/pcimpl.h" /*I "petscpc.h" I*/ 7 8 typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 9 struct _PC_FieldSplitLink { 10 KSP ksp; 11 Vec x,y; 12 PetscInt nfields; 13 PetscInt *fields; 14 VecScatter sctx; 15 PC_FieldSplitLink next; 16 }; 17 18 typedef struct { 19 PCCompositeType type; /* additive or multiplicative */ 20 PetscTruth defaultsplit; 21 PetscInt bs; 22 PetscInt nsplits; 23 Vec *x,*y,w1,w2; 24 Mat *pmat; 25 IS *is; 26 PC_FieldSplitLink head; 27 } PC_FieldSplit; 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "PCView_FieldSplit" 31 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 32 { 33 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 34 PetscErrorCode ierr; 35 PetscTruth iascii; 36 PetscInt i,j; 37 PC_FieldSplitLink ilink = jac->head; 38 39 PetscFunctionBegin; 40 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 41 if (iascii) { 42 ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 43 ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 44 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 45 for (i=0; i<jac->nsplits; i++) { 46 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 47 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 48 for (j=0; j<ilink->nfields; j++) { 49 if (j > 0) { 50 ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 51 } 52 ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 53 } 54 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 55 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 56 ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 57 ilink = ilink->next; 58 } 59 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 60 } else { 61 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 62 } 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "PCFieldSplitSetDefaults" 68 static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 69 { 70 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 71 PetscErrorCode ierr; 72 PC_FieldSplitLink ilink = jac->head; 73 PetscInt i; 74 PetscTruth flg = PETSC_FALSE,*fields; 75 76 PetscFunctionBegin; 77 ierr = PetscOptionsGetTruth(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 78 if (!ilink || flg) { 79 ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 80 if (jac->bs <= 0) { 81 ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 82 } 83 ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 84 ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 85 while (ilink) { 86 for (i=0; i<ilink->nfields; i++) { 87 fields[ilink->fields[i]] = PETSC_TRUE; 88 } 89 ilink = ilink->next; 90 } 91 jac->defaultsplit = PETSC_TRUE; 92 for (i=0; i<jac->bs; i++) { 93 if (!fields[i]) { 94 ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 95 } else { 96 jac->defaultsplit = PETSC_FALSE; 97 } 98 } 99 } 100 PetscFunctionReturn(0); 101 } 102 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "PCSetUp_FieldSplit" 106 static PetscErrorCode PCSetUp_FieldSplit(PC pc) 107 { 108 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 109 PetscErrorCode ierr; 110 PC_FieldSplitLink ilink; 111 PetscInt i,nsplit; 112 MatStructure flag = pc->flag; 113 114 PetscFunctionBegin; 115 ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 116 nsplit = jac->nsplits; 117 ilink = jac->head; 118 119 /* get the matrices for each split */ 120 if (!jac->is) { 121 PetscInt rstart,rend,nslots,bs; 122 123 ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 124 ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 125 nslots = (rend - rstart)/bs; 126 ierr = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr); 127 for (i=0; i<nsplit; i++) { 128 if (jac->defaultsplit) { 129 ierr = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr); 130 } else { 131 PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 132 PetscTruth sorted; 133 ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 134 for (j=0; j<nslots; j++) { 135 for (k=0; k<nfields; k++) { 136 ii[nfields*j + k] = rstart + bs*j + fields[k]; 137 } 138 } 139 ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr); 140 ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr); 141 if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 142 ierr = PetscFree(ii);CHKERRQ(ierr); 143 ilink = ilink->next; 144 } 145 } 146 } 147 148 if (!jac->pmat) { 149 ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);CHKERRQ(ierr); 150 } else { 151 ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);CHKERRQ(ierr); 152 } 153 154 /* set up the individual PCs */ 155 i = 0; 156 ilink = jac->head; 157 while (ilink) { 158 ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 159 ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 160 ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 161 i++; 162 ilink = ilink->next; 163 } 164 165 /* create work vectors for each split */ 166 if (!jac->x) { 167 Vec xtmp; 168 ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 169 ilink = jac->head; 170 for (i=0; i<nsplit; i++) { 171 Mat A; 172 ierr = KSPGetOperators(ilink->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr); 173 ierr = MatGetVecs(A,&ilink->x,&ilink->y);CHKERRQ(ierr); 174 jac->x[i] = ilink->x; 175 jac->y[i] = ilink->y; 176 ilink = ilink->next; 177 } 178 /* compute scatter contexts needed by multiplicative versions and non-default splits */ 179 180 ilink = jac->head; 181 ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 182 for (i=0; i<nsplit; i++) { 183 ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 184 ilink = ilink->next; 185 } 186 ierr = VecDestroy(xtmp);CHKERRQ(ierr); 187 } 188 PetscFunctionReturn(0); 189 } 190 191 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 192 (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 193 VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 194 KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 195 VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \ 196 VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx)) 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "PCApply_FieldSplit" 200 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 201 { 202 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 203 PetscErrorCode ierr; 204 PC_FieldSplitLink ilink = jac->head; 205 PetscInt bs; 206 207 PetscFunctionBegin; 208 ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 209 if (bs != jac->bs) { 210 SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector blocksize %D does not match Fieldsplit blocksize %D",bs,jac->bs); 211 } 212 if (jac->type == PC_COMPOSITE_ADDITIVE) { 213 if (jac->defaultsplit) { 214 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 215 while (ilink) { 216 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 217 ilink = ilink->next; 218 } 219 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 220 } else { 221 PetscInt i = 0; 222 223 ierr = VecSet(y,0.0);CHKERRQ(ierr); 224 while (ilink) { 225 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 226 ilink = ilink->next; 227 i++; 228 } 229 } 230 } else { 231 if (!jac->w1) { 232 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 233 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 234 } 235 ierr = VecSet(y,0.0);CHKERRQ(ierr); 236 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 237 while (ilink->next) { 238 ilink = ilink->next; 239 ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 240 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 241 ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 242 } 243 } 244 PetscFunctionReturn(0); 245 } 246 247 #undef __FUNCT__ 248 #define __FUNCT__ "PCDestroy_FieldSplit" 249 static PetscErrorCode PCDestroy_FieldSplit(PC pc) 250 { 251 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 252 PetscErrorCode ierr; 253 PC_FieldSplitLink ilink = jac->head,next; 254 255 PetscFunctionBegin; 256 while (ilink) { 257 ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 258 if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 259 if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 260 if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 261 next = ilink->next; 262 ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr); 263 ilink = next; 264 } 265 ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 266 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 267 if (jac->is) { 268 PetscInt i; 269 for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);} 270 ierr = PetscFree(jac->is);CHKERRQ(ierr); 271 } 272 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 273 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 274 ierr = PetscFree(jac);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 280 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 281 { 282 PetscErrorCode ierr; 283 PetscInt i = 0,nfields,fields[12]; 284 PetscTruth flg; 285 char optionname[128]; 286 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 287 288 PetscFunctionBegin; 289 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 290 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 291 while (PETSC_TRUE) { 292 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 293 nfields = 12; 294 ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 295 if (!flg) break; 296 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 297 ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 298 } 299 ierr = PetscOptionsTail();CHKERRQ(ierr); 300 PetscFunctionReturn(0); 301 } 302 303 /*------------------------------------------------------------------------------------*/ 304 305 EXTERN_C_BEGIN 306 #undef __FUNCT__ 307 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 308 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 309 { 310 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 311 PetscErrorCode ierr; 312 PC_FieldSplitLink ilink,next = jac->head; 313 char prefix[128]; 314 315 PetscFunctionBegin; 316 if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 317 ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr); 318 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 319 ilink->nfields = n; 320 ilink->next = PETSC_NULL; 321 ierr = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr); 322 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 323 324 if (pc->prefix) { 325 sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits); 326 } else { 327 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 328 } 329 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 330 331 if (!next) { 332 jac->head = ilink; 333 } else { 334 while (next->next) { 335 next = next->next; 336 } 337 next->next = ilink; 338 } 339 jac->nsplits++; 340 PetscFunctionReturn(0); 341 } 342 EXTERN_C_END 343 344 345 EXTERN_C_BEGIN 346 #undef __FUNCT__ 347 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 348 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 349 { 350 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 351 PetscErrorCode ierr; 352 PetscInt cnt = 0; 353 PC_FieldSplitLink ilink = jac->head; 354 355 PetscFunctionBegin; 356 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 357 while (ilink) { 358 (*subksp)[cnt++] = ilink->ksp; 359 ilink = ilink->next; 360 } 361 if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 362 *n = jac->nsplits; 363 PetscFunctionReturn(0); 364 } 365 EXTERN_C_END 366 367 #undef __FUNCT__ 368 #define __FUNCT__ "PCFieldSplitSetFields" 369 /*@ 370 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 371 372 Collective on PC 373 374 Input Parameters: 375 + pc - the preconditioner context 376 . n - the number of fields in this split 377 . fields - the fields in this split 378 379 Level: intermediate 380 381 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT 382 383 @*/ 384 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 385 { 386 PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 387 388 PetscFunctionBegin; 389 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 390 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 391 if (f) { 392 ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 393 } 394 PetscFunctionReturn(0); 395 } 396 397 #undef __FUNCT__ 398 #define __FUNCT__ "PCFieldSplitGetSubKSP" 399 /*@C 400 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 401 402 Collective on KSP 403 404 Input Parameter: 405 . pc - the preconditioner context 406 407 Output Parameters: 408 + n - the number of split 409 - pc - the array of KSP contexts 410 411 Note: 412 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 413 414 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 415 416 Level: advanced 417 418 .seealso: PCFIELDSPLIT 419 @*/ 420 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 421 { 422 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 423 424 PetscFunctionBegin; 425 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 426 PetscValidIntPointer(n,2); 427 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 428 if (f) { 429 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 430 } else { 431 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 432 } 433 PetscFunctionReturn(0); 434 } 435 436 EXTERN_C_BEGIN 437 #undef __FUNCT__ 438 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 439 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 440 { 441 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 442 443 PetscFunctionBegin; 444 jac->type = type; 445 PetscFunctionReturn(0); 446 } 447 EXTERN_C_END 448 449 #undef __FUNCT__ 450 #define __FUNCT__ "PCFieldSplitSetType" 451 /*@C 452 PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 453 454 Collective on PC 455 456 Input Parameter: 457 . pc - the preconditioner context 458 . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE 459 460 Options Database Key: 461 . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type 462 463 Level: Developer 464 465 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 466 467 .seealso: PCCompositeSetType() 468 469 @*/ 470 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 471 { 472 PetscErrorCode ierr,(*f)(PC,PCCompositeType); 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 476 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 477 if (f) { 478 ierr = (*f)(pc,type);CHKERRQ(ierr); 479 } 480 PetscFunctionReturn(0); 481 } 482 483 /* -------------------------------------------------------------------------------------*/ 484 /*MC 485 PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 486 fields or groups of fields 487 488 489 To set options on the solvers for each block append -sub_ to all the PC 490 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 491 492 To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 493 and set the options directly on the resulting KSP object 494 495 Level: intermediate 496 497 Options Database Keys: 498 + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 499 . -pc_splitfield_default - automatically add any fields to additional splits that have not 500 been supplied explicitly by -pc_splitfield_%d_fields 501 - -pc_splitfield_type <additive,multiplicative> 502 503 Concepts: physics based preconditioners 504 505 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 506 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 507 M*/ 508 509 EXTERN_C_BEGIN 510 #undef __FUNCT__ 511 #define __FUNCT__ "PCCreate_FieldSplit" 512 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 513 { 514 PetscErrorCode ierr; 515 PC_FieldSplit *jac; 516 517 PetscFunctionBegin; 518 ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr); 519 ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr); 520 jac->bs = -1; 521 jac->nsplits = 0; 522 jac->type = PC_COMPOSITE_ADDITIVE; 523 pc->data = (void*)jac; 524 525 pc->ops->apply = PCApply_FieldSplit; 526 pc->ops->setup = PCSetUp_FieldSplit; 527 pc->ops->destroy = PCDestroy_FieldSplit; 528 pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 529 pc->ops->view = PCView_FieldSplit; 530 pc->ops->applyrichardson = 0; 531 532 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 533 PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 534 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 535 PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 536 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 537 PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 538 PetscFunctionReturn(0); 539 } 540 EXTERN_C_END 541 542 543