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