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 IS is,cis; 16 PetscInt csize; 17 PC_FieldSplitLink next,previous; 18 }; 19 20 typedef struct { 21 PCCompositeType type; /* additive or multiplicative */ 22 PetscTruth defaultsplit; 23 PetscInt bs; 24 PetscInt nsplits; 25 Vec *x,*y,w1,w2; 26 Mat *pmat; 27 PetscTruth issetup; 28 PC_FieldSplitLink head; 29 } PC_FieldSplit; 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "PCView_FieldSplit" 33 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 34 { 35 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 36 PetscErrorCode ierr; 37 PetscTruth iascii; 38 PetscInt i,j; 39 PC_FieldSplitLink ilink = jac->head; 40 41 PetscFunctionBegin; 42 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 43 if (iascii) { 44 ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 45 ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 46 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 47 for (i=0; i<jac->nsplits; i++) { 48 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 49 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 50 for (j=0; j<ilink->nfields; j++) { 51 if (j > 0) { 52 ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 53 } 54 ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 55 } 56 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 57 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 58 ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 59 ilink = ilink->next; 60 } 61 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 62 } else { 63 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 64 } 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "PCFieldSplitSetDefaults" 70 static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 71 { 72 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 73 PetscErrorCode ierr; 74 PC_FieldSplitLink ilink = jac->head; 75 PetscInt i; 76 PetscTruth flg = PETSC_FALSE,*fields; 77 78 PetscFunctionBegin; 79 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 80 if (!ilink || flg) { 81 ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 82 if (jac->bs <= 0) { 83 if (pc->pmat) { 84 ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 85 } else { 86 jac->bs = 1; 87 } 88 } 89 ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 90 ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 91 while (ilink) { 92 for (i=0; i<ilink->nfields; i++) { 93 fields[ilink->fields[i]] = PETSC_TRUE; 94 } 95 ilink = ilink->next; 96 } 97 jac->defaultsplit = PETSC_TRUE; 98 for (i=0; i<jac->bs; i++) { 99 if (!fields[i]) { 100 ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 101 } else { 102 jac->defaultsplit = PETSC_FALSE; 103 } 104 } 105 } 106 PetscFunctionReturn(0); 107 } 108 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "PCSetUp_FieldSplit" 112 static PetscErrorCode PCSetUp_FieldSplit(PC pc) 113 { 114 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 115 PetscErrorCode ierr; 116 PC_FieldSplitLink ilink; 117 PetscInt i,nsplit,ccsize; 118 MatStructure flag = pc->flag; 119 PetscTruth sorted; 120 121 PetscFunctionBegin; 122 ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 123 nsplit = jac->nsplits; 124 ilink = jac->head; 125 126 /* get the matrices for each split */ 127 if (!jac->issetup) { 128 PetscInt rstart,rend,nslots,bs; 129 130 jac->issetup = PETSC_TRUE; 131 132 /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 133 bs = jac->bs; 134 ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 135 ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 136 nslots = (rend - rstart)/bs; 137 for (i=0; i<nsplit; i++) { 138 if (jac->defaultsplit) { 139 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 140 ilink->csize = ccsize/nsplit; 141 } else if (!ilink->is) { 142 if (ilink->nfields > 1) { 143 PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 144 ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 145 for (j=0; j<nslots; j++) { 146 for (k=0; k<nfields; k++) { 147 ii[nfields*j + k] = rstart + bs*j + fields[k]; 148 } 149 } 150 ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 151 ierr = PetscFree(ii);CHKERRQ(ierr); 152 } else { 153 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 154 } 155 ilink->csize = (ccsize/bs)*ilink->nfields; 156 ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 157 if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 158 } 159 ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 160 ilink = ilink->next; 161 } 162 } 163 164 ilink = jac->head; 165 if (!jac->pmat) { 166 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 167 for (i=0; i<nsplit; i++) { 168 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 169 ilink = ilink->next; 170 } 171 } else { 172 for (i=0; i<nsplit; i++) { 173 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 174 ilink = ilink->next; 175 } 176 } 177 178 /* set up the individual PCs */ 179 i = 0; 180 ilink = jac->head; 181 while (ilink) { 182 ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 183 ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 184 ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 185 i++; 186 ilink = ilink->next; 187 } 188 189 /* create work vectors for each split */ 190 if (!jac->x) { 191 Vec xtmp; 192 ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 193 ilink = jac->head; 194 for (i=0; i<nsplit; i++) { 195 Vec *vl,*vr; 196 197 ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 198 ilink->x = *vr; 199 ilink->y = *vl; 200 ierr = PetscFree(vr);CHKERRQ(ierr); 201 ierr = PetscFree(vl);CHKERRQ(ierr); 202 jac->x[i] = ilink->x; 203 jac->y[i] = ilink->y; 204 ilink = ilink->next; 205 } 206 /* compute scatter contexts needed by multiplicative versions and non-default splits */ 207 208 ilink = jac->head; 209 ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 210 for (i=0; i<nsplit; i++) { 211 ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 212 ilink = ilink->next; 213 } 214 ierr = VecDestroy(xtmp);CHKERRQ(ierr); 215 } 216 PetscFunctionReturn(0); 217 } 218 219 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 220 (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 221 VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 222 KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 223 VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 224 VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 225 226 #undef __FUNCT__ 227 #define __FUNCT__ "PCApply_FieldSplit" 228 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 229 { 230 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 231 PetscErrorCode ierr; 232 PC_FieldSplitLink ilink = jac->head; 233 PetscInt bs; 234 235 PetscFunctionBegin; 236 CHKMEMQ; 237 ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 238 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 239 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 240 241 if (jac->type == PC_COMPOSITE_ADDITIVE) { 242 if (jac->defaultsplit) { 243 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 244 while (ilink) { 245 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 246 ilink = ilink->next; 247 } 248 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 249 } else { 250 ierr = VecSet(y,0.0);CHKERRQ(ierr); 251 while (ilink) { 252 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 253 ilink = ilink->next; 254 } 255 } 256 } else { 257 if (!jac->w1) { 258 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 259 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 260 } 261 ierr = VecSet(y,0.0);CHKERRQ(ierr); 262 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 263 while (ilink->next) { 264 ilink = ilink->next; 265 ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 266 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 267 ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 268 } 269 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 270 while (ilink->previous) { 271 ilink = ilink->previous; 272 ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 273 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 274 ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 275 } 276 } 277 } 278 CHKMEMQ; 279 PetscFunctionReturn(0); 280 } 281 282 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 283 (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 284 VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 285 KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 286 VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 287 VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "PCApply_FieldSplit" 291 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 292 { 293 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 294 PetscErrorCode ierr; 295 PC_FieldSplitLink ilink = jac->head; 296 PetscInt bs; 297 298 PetscFunctionBegin; 299 CHKMEMQ; 300 ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 301 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 302 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 303 304 if (jac->type == PC_COMPOSITE_ADDITIVE) { 305 if (jac->defaultsplit) { 306 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 307 while (ilink) { 308 ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 309 ilink = ilink->next; 310 } 311 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 312 } else { 313 ierr = VecSet(y,0.0);CHKERRQ(ierr); 314 while (ilink) { 315 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 316 ilink = ilink->next; 317 } 318 } 319 } else { 320 if (!jac->w1) { 321 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 322 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 323 } 324 ierr = VecSet(y,0.0);CHKERRQ(ierr); 325 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 326 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 327 while (ilink->next) { 328 ilink = ilink->next; 329 ierr = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr); 330 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 331 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 332 } 333 while (ilink->previous) { 334 ilink = ilink->previous; 335 ierr = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr); 336 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 337 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 338 } 339 } else { 340 while (ilink->next) { /* get to last entry in linked list */ 341 ilink = ilink->next; 342 } 343 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 344 while (ilink->previous) { 345 ilink = ilink->previous; 346 ierr = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr); 347 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 348 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 349 } 350 } 351 } 352 CHKMEMQ; 353 PetscFunctionReturn(0); 354 } 355 356 #undef __FUNCT__ 357 #define __FUNCT__ "PCDestroy_FieldSplit" 358 static PetscErrorCode PCDestroy_FieldSplit(PC pc) 359 { 360 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 361 PetscErrorCode ierr; 362 PC_FieldSplitLink ilink = jac->head,next; 363 364 PetscFunctionBegin; 365 while (ilink) { 366 ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 367 if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 368 if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 369 if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 370 if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 371 if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 372 next = ilink->next; 373 ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 374 ierr = PetscFree(ilink);CHKERRQ(ierr); 375 ilink = next; 376 } 377 ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 378 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 379 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 380 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 381 ierr = PetscFree(jac);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 387 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 388 { 389 PetscErrorCode ierr; 390 PetscInt i = 0,nfields,*fields,bs; 391 PetscTruth flg; 392 char optionname[128]; 393 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 394 395 PetscFunctionBegin; 396 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 397 ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 398 if (flg) { 399 ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 400 } 401 if (jac->bs <= 0) { 402 403 ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr); 404 } 405 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 406 ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 407 while (PETSC_TRUE) { 408 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 409 nfields = jac->bs; 410 ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 411 if (!flg) break; 412 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 413 ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 414 } 415 ierr = PetscFree(fields);CHKERRQ(ierr); 416 ierr = PetscOptionsTail();CHKERRQ(ierr); 417 PetscFunctionReturn(0); 418 } 419 420 /*------------------------------------------------------------------------------------*/ 421 422 EXTERN_C_BEGIN 423 #undef __FUNCT__ 424 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 425 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 426 { 427 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 428 PetscErrorCode ierr; 429 PC_FieldSplitLink ilink,next = jac->head; 430 char prefix[128]; 431 PetscInt i; 432 433 PetscFunctionBegin; 434 if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 435 for (i=0; i<n; i++) { 436 if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 437 if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 438 } 439 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 440 ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 441 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 442 ilink->nfields = n; 443 ilink->next = PETSC_NULL; 444 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 445 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 446 447 if (((PetscObject)pc)->prefix) { 448 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 449 } else { 450 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 451 } 452 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 453 454 if (!next) { 455 jac->head = ilink; 456 ilink->previous = PETSC_NULL; 457 } else { 458 while (next->next) { 459 next = next->next; 460 } 461 next->next = ilink; 462 ilink->previous = next; 463 } 464 jac->nsplits++; 465 PetscFunctionReturn(0); 466 } 467 EXTERN_C_END 468 469 470 EXTERN_C_BEGIN 471 #undef __FUNCT__ 472 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 473 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 474 { 475 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 476 PetscErrorCode ierr; 477 PetscInt cnt = 0; 478 PC_FieldSplitLink ilink = jac->head; 479 480 PetscFunctionBegin; 481 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 482 while (ilink) { 483 (*subksp)[cnt++] = ilink->ksp; 484 ilink = ilink->next; 485 } 486 if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 487 *n = jac->nsplits; 488 PetscFunctionReturn(0); 489 } 490 EXTERN_C_END 491 492 EXTERN_C_BEGIN 493 #undef __FUNCT__ 494 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 495 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 496 { 497 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 498 PetscErrorCode ierr; 499 PC_FieldSplitLink ilink, next = jac->head; 500 char prefix[128]; 501 502 PetscFunctionBegin; 503 ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,0,PetscInt,&ilink->fields);CHKERRQ(ierr); 504 ilink->nfields = -1; 505 ilink->next = PETSC_NULL; 506 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 507 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 508 509 if (((PetscObject)pc)->prefix) { 510 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 511 } else { 512 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 513 } 514 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 515 516 if (!next) { 517 jac->head = ilink; 518 ilink->previous = PETSC_NULL; 519 } else { 520 while (next->next) { 521 next = next->next; 522 } 523 next->next = ilink; 524 ilink->previous = next; 525 } 526 jac->nsplits++; 527 528 PetscFunctionReturn(0); 529 } 530 EXTERN_C_END 531 532 #undef __FUNCT__ 533 #define __FUNCT__ "PCFieldSplitSetFields" 534 /*@ 535 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 536 537 Collective on PC 538 539 Input Parameters: 540 + pc - the preconditioner context 541 . n - the number of fields in this split 542 . fields - the fields in this split 543 544 Level: intermediate 545 546 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 547 548 @*/ 549 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 550 { 551 PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 552 553 PetscFunctionBegin; 554 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 555 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 556 if (f) { 557 ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 558 } 559 PetscFunctionReturn(0); 560 } 561 562 #undef __FUNCT__ 563 #define __FUNCT__ "PCFieldSplitSetIS" 564 /*@ 565 PCFieldSplitSetIS - Sets the exact elements for field 566 567 Collective on PC 568 569 Input Parameters: 570 + pc - the preconditioner context 571 . is - the index set that defines the vector elements in this field 572 573 Level: intermediate 574 575 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 576 577 @*/ 578 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 579 { 580 PetscErrorCode ierr,(*f)(PC,IS); 581 582 PetscFunctionBegin; 583 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 584 PetscValidHeaderSpecific(is,IS_COOKIE,1); 585 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 586 if (f) { 587 ierr = (*f)(pc,is);CHKERRQ(ierr); 588 } 589 PetscFunctionReturn(0); 590 } 591 592 #undef __FUNCT__ 593 #define __FUNCT__ "PCFieldSplitSetBlockSize" 594 /*@ 595 PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 596 fieldsplit preconditioner. If not set the matrix block size is used. 597 598 Collective on PC 599 600 Input Parameters: 601 + pc - the preconditioner context 602 - bs - the block size 603 604 Level: intermediate 605 606 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 607 608 @*/ 609 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 610 { 611 PetscErrorCode ierr,(*f)(PC,PetscInt); 612 613 PetscFunctionBegin; 614 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 615 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 616 if (f) { 617 ierr = (*f)(pc,bs);CHKERRQ(ierr); 618 } 619 PetscFunctionReturn(0); 620 } 621 622 #undef __FUNCT__ 623 #define __FUNCT__ "PCFieldSplitGetSubKSP" 624 /*@C 625 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 626 627 Collective on KSP 628 629 Input Parameter: 630 . pc - the preconditioner context 631 632 Output Parameters: 633 + n - the number of split 634 - pc - the array of KSP contexts 635 636 Note: 637 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 638 639 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 640 641 Level: advanced 642 643 .seealso: PCFIELDSPLIT 644 @*/ 645 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 646 { 647 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 648 649 PetscFunctionBegin; 650 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 651 PetscValidIntPointer(n,2); 652 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 653 if (f) { 654 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 655 } else { 656 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 657 } 658 PetscFunctionReturn(0); 659 } 660 661 EXTERN_C_BEGIN 662 #undef __FUNCT__ 663 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 664 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 665 { 666 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 667 668 PetscFunctionBegin; 669 jac->type = type; 670 PetscFunctionReturn(0); 671 } 672 EXTERN_C_END 673 674 EXTERN_C_BEGIN 675 #undef __FUNCT__ 676 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 677 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 678 { 679 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 680 681 PetscFunctionBegin; 682 if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 683 if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 684 jac->bs = bs; 685 PetscFunctionReturn(0); 686 } 687 EXTERN_C_END 688 689 #undef __FUNCT__ 690 #define __FUNCT__ "PCFieldSplitSetType" 691 /*@ 692 PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 693 694 Collective on PC 695 696 Input Parameter: 697 . pc - the preconditioner context 698 . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 699 700 Options Database Key: 701 . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type 702 703 Level: Developer 704 705 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 706 707 .seealso: PCCompositeSetType() 708 709 @*/ 710 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 711 { 712 PetscErrorCode ierr,(*f)(PC,PCCompositeType); 713 714 PetscFunctionBegin; 715 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 716 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 717 if (f) { 718 ierr = (*f)(pc,type);CHKERRQ(ierr); 719 } 720 PetscFunctionReturn(0); 721 } 722 723 /* -------------------------------------------------------------------------------------*/ 724 /*MC 725 PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 726 fields or groups of fields 727 728 729 To set options on the solvers for each block append -sub_ to all the PC 730 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 731 732 To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 733 and set the options directly on the resulting KSP object 734 735 Level: intermediate 736 737 Options Database Keys: 738 + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 739 . -pc_splitfield_default - automatically add any fields to additional splits that have not 740 been supplied explicitly by -pc_splitfield_%d_fields 741 . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 742 - -pc_splitfield_type <additive,multiplicative> 743 744 Concepts: physics based preconditioners 745 746 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 747 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 748 M*/ 749 750 EXTERN_C_BEGIN 751 #undef __FUNCT__ 752 #define __FUNCT__ "PCCreate_FieldSplit" 753 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 754 { 755 PetscErrorCode ierr; 756 PC_FieldSplit *jac; 757 758 PetscFunctionBegin; 759 ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 760 jac->bs = -1; 761 jac->nsplits = 0; 762 jac->type = PC_COMPOSITE_ADDITIVE; 763 764 pc->data = (void*)jac; 765 766 pc->ops->apply = PCApply_FieldSplit; 767 pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 768 pc->ops->setup = PCSetUp_FieldSplit; 769 pc->ops->destroy = PCDestroy_FieldSplit; 770 pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 771 pc->ops->view = PCView_FieldSplit; 772 pc->ops->applyrichardson = 0; 773 774 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 775 PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 776 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 777 PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 778 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 779 PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 780 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 781 PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 782 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 783 PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 784 PetscFunctionReturn(0); 785 } 786 EXTERN_C_END 787 788 789