1 #define PETSCKSP_DLL 2 3 /* 4 5 */ 6 #include "src/ksp/pc/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 = PetscLogInfo((pc,"PCFieldSplitSetDefaults: 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 206 PetscFunctionBegin; 207 if (jac->type == PC_COMPOSITE_ADDITIVE) { 208 if (jac->defaultsplit) { 209 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 210 while (ilink) { 211 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 212 ilink = ilink->next; 213 } 214 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 215 } else { 216 PetscInt i = 0; 217 218 ierr = VecSet(y,0.0);CHKERRQ(ierr); 219 while (ilink) { 220 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 221 ilink = ilink->next; 222 i++; 223 } 224 } 225 } else { 226 if (!jac->w1) { 227 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 228 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 229 } 230 ierr = VecSet(y,0.0);CHKERRQ(ierr); 231 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 232 while (ilink->next) { 233 ilink = ilink->next; 234 ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 235 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 236 ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 237 } 238 } 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "PCDestroy_FieldSplit" 244 static PetscErrorCode PCDestroy_FieldSplit(PC pc) 245 { 246 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 247 PetscErrorCode ierr; 248 PC_FieldSplitLink ilink = jac->head,next; 249 250 PetscFunctionBegin; 251 while (ilink) { 252 ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 253 if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 254 if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 255 if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 256 next = ilink->next; 257 ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr); 258 ilink = next; 259 } 260 if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);} 261 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 262 if (jac->is) { 263 PetscInt i; 264 for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);} 265 ierr = PetscFree(jac->is);CHKERRQ(ierr); 266 } 267 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 268 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 269 ierr = PetscFree(jac);CHKERRQ(ierr); 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 275 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 276 /* This does not call KSPSetFromOptions() on the subksp's, see PCSetFromOptionsBJacobi/ASM() */ 277 { 278 PetscErrorCode ierr; 279 PetscInt i = 0,nfields,fields[12]; 280 PetscTruth flg; 281 char optionname[128]; 282 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 283 284 PetscFunctionBegin; 285 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 286 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 287 while (PETSC_TRUE) { 288 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 289 nfields = 12; 290 ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 291 if (!flg) break; 292 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 293 ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 294 } 295 ierr = PetscOptionsTail();CHKERRQ(ierr); 296 PetscFunctionReturn(0); 297 } 298 299 /*------------------------------------------------------------------------------------*/ 300 301 EXTERN_C_BEGIN 302 #undef __FUNCT__ 303 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 304 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 305 { 306 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 307 PetscErrorCode ierr; 308 PC_FieldSplitLink ilink,next = jac->head; 309 char prefix[128]; 310 311 PetscFunctionBegin; 312 if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 313 ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr); 314 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 315 ilink->nfields = n; 316 ilink->next = PETSC_NULL; 317 ierr = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr); 318 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 319 320 if (pc->prefix) { 321 sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits); 322 } else { 323 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 324 } 325 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 326 327 if (!next) { 328 jac->head = ilink; 329 } else { 330 while (next->next) { 331 next = next->next; 332 } 333 next->next = ilink; 334 } 335 jac->nsplits++; 336 PetscFunctionReturn(0); 337 } 338 EXTERN_C_END 339 340 341 EXTERN_C_BEGIN 342 #undef __FUNCT__ 343 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 344 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 345 { 346 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 347 PetscErrorCode ierr; 348 PetscInt cnt = 0; 349 PC_FieldSplitLink ilink = jac->head; 350 351 PetscFunctionBegin; 352 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 353 while (ilink) { 354 (*subksp)[cnt++] = ilink->ksp; 355 ilink = ilink->next; 356 } 357 if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 358 *n = jac->nsplits; 359 PetscFunctionReturn(0); 360 } 361 EXTERN_C_END 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "PCFieldSplitSetFields" 365 /*@ 366 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 367 368 Collective on PC 369 370 Input Parameters: 371 + pc - the preconditioner context 372 . n - the number of fields in this split 373 . fields - the fields in this split 374 375 Level: intermediate 376 377 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT 378 379 @*/ 380 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 381 { 382 PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 383 384 PetscFunctionBegin; 385 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 386 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 387 if (f) { 388 ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 389 } 390 PetscFunctionReturn(0); 391 } 392 393 #undef __FUNCT__ 394 #define __FUNCT__ "PCFieldSplitGetSubKSP" 395 /*@C 396 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 397 398 Collective on KSP 399 400 Input Parameter: 401 . pc - the preconditioner context 402 403 Output Parameters: 404 + n - the number of split 405 - pc - the array of KSP contexts 406 407 Note: 408 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 409 410 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 411 412 Level: advanced 413 414 .seealso: PCFIELDSPLIT 415 @*/ 416 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 417 { 418 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 419 420 PetscFunctionBegin; 421 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 422 PetscValidIntPointer(n,2); 423 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 424 if (f) { 425 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 426 } else { 427 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 428 } 429 PetscFunctionReturn(0); 430 } 431 432 EXTERN_C_BEGIN 433 #undef __FUNCT__ 434 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 435 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 436 { 437 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 438 439 PetscFunctionBegin; 440 jac->type = type; 441 PetscFunctionReturn(0); 442 } 443 EXTERN_C_END 444 445 #undef __FUNCT__ 446 #define __FUNCT__ "PCFieldSplitSetType" 447 /*@C 448 PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 449 450 Collective on PC 451 452 Input Parameter: 453 . pc - the preconditioner context 454 . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE 455 456 Options Database Key: 457 . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type 458 459 Level: Developer 460 461 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 462 463 .seealso: PCCompositeSetType() 464 465 @*/ 466 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 467 { 468 PetscErrorCode ierr,(*f)(PC,PCCompositeType); 469 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 472 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 473 if (f) { 474 ierr = (*f)(pc,type);CHKERRQ(ierr); 475 } 476 PetscFunctionReturn(0); 477 } 478 479 /* -------------------------------------------------------------------------------------*/ 480 /*MC 481 PCFIELDSPLIT - Preconditioner created by combining seperate preconditioners for individual 482 fields or groups of fields 483 484 485 To set options on the solvers for each block append -sub_ to all the PC 486 options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1 487 488 To set the options on the solvers seperate for each block call PCFieldSplitGetSubKSP() 489 and set the options directly on the resulting KSP object 490 491 Level: intermediate 492 493 Options Database Keys: 494 + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 495 . -pc_splitfield_default - automatically add any fields to additional splits that have not 496 been supplied explicitly by -pc_splitfield_%d_fields 497 - -pc_splitfield_type <additive,multiplicative> 498 499 Concepts: physics based preconditioners 500 501 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 502 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 503 M*/ 504 505 EXTERN_C_BEGIN 506 #undef __FUNCT__ 507 #define __FUNCT__ "PCCreate_FieldSplit" 508 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 509 { 510 PetscErrorCode ierr; 511 PC_FieldSplit *jac; 512 513 PetscFunctionBegin; 514 ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr); 515 ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr); 516 jac->bs = -1; 517 jac->nsplits = 0; 518 jac->type = PC_COMPOSITE_ADDITIVE; 519 pc->data = (void*)jac; 520 521 pc->ops->apply = PCApply_FieldSplit; 522 pc->ops->setup = PCSetUp_FieldSplit; 523 pc->ops->destroy = PCDestroy_FieldSplit; 524 pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 525 pc->ops->view = PCView_FieldSplit; 526 pc->ops->applyrichardson = 0; 527 528 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 529 PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 530 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 531 PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 532 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 533 PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 534 PetscFunctionReturn(0); 535 } 536 EXTERN_C_END 537 538 539