1 #define PETSCKSP_DLL 2 3 #include "private/pcimpl.h" /*I "petscpc.h" I*/ 4 5 const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6 7 typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 8 struct _PC_FieldSplitLink { 9 KSP ksp; 10 Vec x,y; 11 char *splitname; 12 PetscInt nfields; 13 PetscInt *fields; 14 VecScatter sctx; 15 IS is; 16 PC_FieldSplitLink next,previous; 17 }; 18 19 typedef struct { 20 PCCompositeType type; 21 PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 22 PetscTruth splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 23 PetscTruth realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 24 PetscInt bs; /* Block size for IS and Mat structures */ 25 PetscInt nsplits; /* Number of field divisions defined */ 26 Vec *x,*y,w1,w2; 27 Mat *mat; /* The diagonal block for each split */ 28 Mat *pmat; /* The preconditioning diagonal block for each split */ 29 Mat *Afield; /* The rows of the matrix associated with each split */ 30 PetscTruth issetup; 31 /* Only used when Schur complement preconditioning is used */ 32 Mat B; /* The (0,1) block */ 33 Mat C; /* The (1,0) block */ 34 Mat schur; /* The Schur complement S = D - C A^{-1} B */ 35 Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 36 PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 37 KSP kspschur; /* The solver for S */ 38 PC_FieldSplitLink head; 39 } PC_FieldSplit; 40 41 /* 42 Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 43 inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 44 PC you could change this. 45 */ 46 47 /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 48 * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 49 static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 50 { 51 switch (jac->schurpre) { 52 case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 53 case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 54 case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 55 default: 56 return jac->schur_user ? jac->schur_user : jac->pmat[1]; 57 } 58 } 59 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "PCView_FieldSplit" 63 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 64 { 65 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 66 PetscErrorCode ierr; 67 PetscTruth iascii; 68 PetscInt i,j; 69 PC_FieldSplitLink ilink = jac->head; 70 71 PetscFunctionBegin; 72 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 73 if (iascii) { 74 ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 75 ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 76 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 77 for (i=0; i<jac->nsplits; i++) { 78 if (ilink->fields) { 79 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 80 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 81 for (j=0; j<ilink->nfields; j++) { 82 if (j > 0) { 83 ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 84 } 85 ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 86 } 87 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 88 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 89 } else { 90 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 91 } 92 ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 93 ilink = ilink->next; 94 } 95 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 96 } else { 97 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 98 } 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "PCView_FieldSplit_Schur" 104 static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 105 { 106 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 107 PetscErrorCode ierr; 108 PetscTruth iascii; 109 PetscInt i,j; 110 PC_FieldSplitLink ilink = jac->head; 111 KSP ksp; 112 113 PetscFunctionBegin; 114 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 115 if (iascii) { 116 ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 117 ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 118 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 119 for (i=0; i<jac->nsplits; i++) { 120 if (ilink->fields) { 121 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 122 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 123 for (j=0; j<ilink->nfields; j++) { 124 if (j > 0) { 125 ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 126 } 127 ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 128 } 129 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 130 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 131 } else { 132 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 133 } 134 ilink = ilink->next; 135 } 136 ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 137 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 138 if (jac->schur) { 139 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 140 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 141 } else { 142 ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 143 } 144 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 145 ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 146 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 147 if (jac->kspschur) { 148 ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 149 } else { 150 ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 151 } 152 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 153 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 154 } else { 155 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 156 } 157 PetscFunctionReturn(0); 158 } 159 160 #undef __FUNCT__ 161 #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 162 /* Precondition: jac->bs is set to a meaningful value */ 163 static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 164 { 165 PetscErrorCode ierr; 166 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 167 PetscInt i,nfields,*ifields; 168 PetscTruth flg; 169 char optionname[128],splitname[8]; 170 171 PetscFunctionBegin; 172 ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 173 for (i=0,flg=PETSC_TRUE; ; i++) { 174 ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 175 ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 176 nfields = jac->bs; 177 ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 178 if (!flg) break; 179 if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 180 ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr); 181 } 182 if (i > 0) { 183 /* Makes command-line setting of splits take precedence over setting them in code. 184 Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 185 create new splits, which would probably not be what the user wanted. */ 186 jac->splitdefined = PETSC_TRUE; 187 } 188 ierr = PetscFree(ifields);CHKERRQ(ierr); 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "PCFieldSplitSetDefaults" 194 static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 195 { 196 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 197 PetscErrorCode ierr; 198 PC_FieldSplitLink ilink = jac->head; 199 PetscTruth flg = PETSC_FALSE; 200 PetscInt i; 201 202 PetscFunctionBegin; 203 if (!ilink) { 204 205 if (jac->bs <= 0) { 206 if (pc->pmat) { 207 ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 208 } else { 209 jac->bs = 1; 210 } 211 } 212 213 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 214 if (!flg) { 215 /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 216 then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 217 ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 218 if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 219 } 220 if (flg || !jac->splitdefined) { 221 ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 222 for (i=0; i<jac->bs; i++) { 223 char splitname[8]; 224 ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 225 ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr); 226 } 227 jac->defaultsplit = PETSC_TRUE; 228 } 229 } else if (jac->nsplits == 1) { 230 if (ilink->is) { 231 IS is2; 232 PetscInt nmin,nmax; 233 234 ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 235 ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 236 ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 237 ierr = ISDestroy(is2);CHKERRQ(ierr); 238 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 239 } 240 if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 241 PetscFunctionReturn(0); 242 } 243 244 245 #undef __FUNCT__ 246 #define __FUNCT__ "PCSetUp_FieldSplit" 247 static PetscErrorCode PCSetUp_FieldSplit(PC pc) 248 { 249 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 250 PetscErrorCode ierr; 251 PC_FieldSplitLink ilink; 252 PetscInt i,nsplit,ccsize; 253 MatStructure flag = pc->flag; 254 PetscTruth sorted; 255 256 PetscFunctionBegin; 257 ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 258 nsplit = jac->nsplits; 259 ilink = jac->head; 260 261 /* get the matrices for each split */ 262 if (!jac->issetup) { 263 PetscInt rstart,rend,nslots,bs; 264 265 jac->issetup = PETSC_TRUE; 266 267 /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 268 bs = jac->bs; 269 ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 270 ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 271 nslots = (rend - rstart)/bs; 272 for (i=0; i<nsplit; i++) { 273 if (jac->defaultsplit) { 274 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 275 } else if (!ilink->is) { 276 if (ilink->nfields > 1) { 277 PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 278 ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 279 for (j=0; j<nslots; j++) { 280 for (k=0; k<nfields; k++) { 281 ii[nfields*j + k] = rstart + bs*j + fields[k]; 282 } 283 } 284 ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 285 ierr = PetscFree(ii);CHKERRQ(ierr); 286 } else { 287 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 288 } 289 } 290 ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 291 if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 292 ilink = ilink->next; 293 } 294 } 295 296 ilink = jac->head; 297 if (!jac->pmat) { 298 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 299 for (i=0; i<nsplit; i++) { 300 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 301 ilink = ilink->next; 302 } 303 } else { 304 for (i=0; i<nsplit; i++) { 305 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 306 ilink = ilink->next; 307 } 308 } 309 if (jac->realdiagonal) { 310 ilink = jac->head; 311 if (!jac->mat) { 312 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 313 for (i=0; i<nsplit; i++) { 314 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 315 ilink = ilink->next; 316 } 317 } else { 318 for (i=0; i<nsplit; i++) { 319 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 320 ilink = ilink->next; 321 } 322 } 323 } else { 324 jac->mat = jac->pmat; 325 } 326 327 if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 328 /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 329 ilink = jac->head; 330 if (!jac->Afield) { 331 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 332 for (i=0; i<nsplit; i++) { 333 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 334 ilink = ilink->next; 335 } 336 } else { 337 for (i=0; i<nsplit; i++) { 338 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 339 ilink = ilink->next; 340 } 341 } 342 } 343 344 if (jac->type == PC_COMPOSITE_SCHUR) { 345 IS ccis; 346 PetscInt rstart,rend; 347 if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 348 349 /* When extracting off-diagonal submatrices, we take complements from this range */ 350 ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 351 352 /* need to handle case when one is resetting up the preconditioner */ 353 if (jac->schur) { 354 ilink = jac->head; 355 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 356 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 357 ierr = ISDestroy(ccis);CHKERRQ(ierr); 358 ilink = ilink->next; 359 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 360 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 361 ierr = ISDestroy(ccis);CHKERRQ(ierr); 362 ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 363 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 364 365 } else { 366 KSP ksp; 367 char schurprefix[256]; 368 369 /* extract the B and C matrices */ 370 ilink = jac->head; 371 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 372 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 373 ierr = ISDestroy(ccis);CHKERRQ(ierr); 374 ilink = ilink->next; 375 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 376 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 377 ierr = ISDestroy(ccis);CHKERRQ(ierr); 378 /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 379 ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 380 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 381 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 382 ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 383 384 ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 385 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 386 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 387 if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 388 PC pc; 389 ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 390 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 391 /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 392 } 393 ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 394 ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 395 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 396 ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 397 398 ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 399 ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 400 ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 401 ilink = jac->head; 402 ilink->x = jac->x[0]; ilink->y = jac->y[0]; 403 ilink = ilink->next; 404 ilink->x = jac->x[1]; ilink->y = jac->y[1]; 405 } 406 } else { 407 /* set up the individual PCs */ 408 i = 0; 409 ilink = jac->head; 410 while (ilink) { 411 ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 412 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 413 ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 414 ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 415 i++; 416 ilink = ilink->next; 417 } 418 419 /* create work vectors for each split */ 420 if (!jac->x) { 421 ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 422 ilink = jac->head; 423 for (i=0; i<nsplit; i++) { 424 Vec *vl,*vr; 425 426 ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 427 ilink->x = *vr; 428 ilink->y = *vl; 429 ierr = PetscFree(vr);CHKERRQ(ierr); 430 ierr = PetscFree(vl);CHKERRQ(ierr); 431 jac->x[i] = ilink->x; 432 jac->y[i] = ilink->y; 433 ilink = ilink->next; 434 } 435 } 436 } 437 438 439 if (!jac->head->sctx) { 440 Vec xtmp; 441 442 /* compute scatter contexts needed by multiplicative versions and non-default splits */ 443 444 ilink = jac->head; 445 ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 446 for (i=0; i<nsplit; i++) { 447 ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 448 ilink = ilink->next; 449 } 450 ierr = VecDestroy(xtmp);CHKERRQ(ierr); 451 } 452 PetscFunctionReturn(0); 453 } 454 455 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 456 (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 457 VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 458 KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 459 VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 460 VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 461 462 #undef __FUNCT__ 463 #define __FUNCT__ "PCApply_FieldSplit_Schur" 464 static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 465 { 466 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 467 PetscErrorCode ierr; 468 KSP ksp; 469 PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 470 471 PetscFunctionBegin; 472 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 473 474 ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 475 ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 476 ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 477 ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 478 ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 479 ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 480 ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 481 482 ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 483 ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 484 ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 485 486 ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 487 ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 488 ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 489 ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 490 ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 491 492 PetscFunctionReturn(0); 493 } 494 495 #undef __FUNCT__ 496 #define __FUNCT__ "PCApply_FieldSplit" 497 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 498 { 499 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 500 PetscErrorCode ierr; 501 PC_FieldSplitLink ilink = jac->head; 502 PetscInt cnt; 503 504 PetscFunctionBegin; 505 CHKMEMQ; 506 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 507 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 508 509 if (jac->type == PC_COMPOSITE_ADDITIVE) { 510 if (jac->defaultsplit) { 511 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 512 while (ilink) { 513 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 514 ilink = ilink->next; 515 } 516 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 517 } else { 518 ierr = VecSet(y,0.0);CHKERRQ(ierr); 519 while (ilink) { 520 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 521 ilink = ilink->next; 522 } 523 } 524 } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 525 if (!jac->w1) { 526 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 527 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 528 } 529 ierr = VecSet(y,0.0);CHKERRQ(ierr); 530 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 531 cnt = 1; 532 while (ilink->next) { 533 ilink = ilink->next; 534 /* compute the residual only over the part of the vector needed */ 535 ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 536 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 537 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 538 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 539 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 540 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 541 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 542 } 543 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 544 cnt -= 2; 545 while (ilink->previous) { 546 ilink = ilink->previous; 547 /* compute the residual only over the part of the vector needed */ 548 ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 549 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 550 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 551 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 552 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 553 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 554 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 555 } 556 } 557 } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 558 CHKMEMQ; 559 PetscFunctionReturn(0); 560 } 561 562 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 563 (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 564 VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 565 KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 566 VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 567 VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "PCApply_FieldSplit" 571 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 572 { 573 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 574 PetscErrorCode ierr; 575 PC_FieldSplitLink ilink = jac->head; 576 577 PetscFunctionBegin; 578 CHKMEMQ; 579 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 580 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 581 582 if (jac->type == PC_COMPOSITE_ADDITIVE) { 583 if (jac->defaultsplit) { 584 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 585 while (ilink) { 586 ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 587 ilink = ilink->next; 588 } 589 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 590 } else { 591 ierr = VecSet(y,0.0);CHKERRQ(ierr); 592 while (ilink) { 593 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 594 ilink = ilink->next; 595 } 596 } 597 } else { 598 if (!jac->w1) { 599 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 600 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 601 } 602 ierr = VecSet(y,0.0);CHKERRQ(ierr); 603 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 604 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 605 while (ilink->next) { 606 ilink = ilink->next; 607 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 608 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 609 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 610 } 611 while (ilink->previous) { 612 ilink = ilink->previous; 613 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 614 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 615 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 616 } 617 } else { 618 while (ilink->next) { /* get to last entry in linked list */ 619 ilink = ilink->next; 620 } 621 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 622 while (ilink->previous) { 623 ilink = ilink->previous; 624 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 625 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 626 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 627 } 628 } 629 } 630 CHKMEMQ; 631 PetscFunctionReturn(0); 632 } 633 634 #undef __FUNCT__ 635 #define __FUNCT__ "PCDestroy_FieldSplit" 636 static PetscErrorCode PCDestroy_FieldSplit(PC pc) 637 { 638 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 639 PetscErrorCode ierr; 640 PC_FieldSplitLink ilink = jac->head,next; 641 642 PetscFunctionBegin; 643 while (ilink) { 644 ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 645 if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 646 if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 647 if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 648 if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 649 next = ilink->next; 650 ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 651 ierr = PetscFree(ilink);CHKERRQ(ierr); 652 ilink = next; 653 } 654 ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 655 if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 656 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 657 if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 658 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 659 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 660 if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 661 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 662 if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 663 if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 664 if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 665 ierr = PetscFree(jac);CHKERRQ(ierr); 666 PetscFunctionReturn(0); 667 } 668 669 #undef __FUNCT__ 670 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 671 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 672 { 673 PetscErrorCode ierr; 674 PetscInt bs; 675 PetscTruth flg; 676 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 677 PCCompositeType ctype; 678 679 PetscFunctionBegin; 680 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 681 ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 682 ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 683 if (flg) { 684 ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 685 } 686 687 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 688 if (flg) { 689 ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 690 } 691 692 /* Only setup fields once */ 693 if ((jac->bs > 0) && (jac->nsplits == 0)) { 694 /* only allow user to set fields from command line if bs is already known. 695 otherwise user can set them in PCFieldSplitSetDefaults() */ 696 ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 697 if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 698 } 699 ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,PETSC_NULL);CHKERRQ(ierr); 700 ierr = PetscOptionsTail();CHKERRQ(ierr); 701 PetscFunctionReturn(0); 702 } 703 704 /*------------------------------------------------------------------------------------*/ 705 706 EXTERN_C_BEGIN 707 #undef __FUNCT__ 708 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 709 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 710 { 711 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 712 PetscErrorCode ierr; 713 PC_FieldSplitLink ilink,next = jac->head; 714 char prefix[128]; 715 PetscInt i; 716 717 PetscFunctionBegin; 718 if (jac->splitdefined) { 719 ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 720 PetscFunctionReturn(0); 721 } 722 for (i=0; i<n; i++) { 723 if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 724 if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 725 } 726 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 727 ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 728 ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 729 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 730 ilink->nfields = n; 731 ilink->next = PETSC_NULL; 732 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 733 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 734 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 735 736 ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr); 737 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 738 739 if (!next) { 740 jac->head = ilink; 741 ilink->previous = PETSC_NULL; 742 } else { 743 while (next->next) { 744 next = next->next; 745 } 746 next->next = ilink; 747 ilink->previous = next; 748 } 749 jac->nsplits++; 750 PetscFunctionReturn(0); 751 } 752 EXTERN_C_END 753 754 EXTERN_C_BEGIN 755 #undef __FUNCT__ 756 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 757 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 758 { 759 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 760 PetscErrorCode ierr; 761 762 PetscFunctionBegin; 763 ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 764 ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 765 (*subksp)[1] = jac->kspschur; 766 *n = jac->nsplits; 767 PetscFunctionReturn(0); 768 } 769 EXTERN_C_END 770 771 EXTERN_C_BEGIN 772 #undef __FUNCT__ 773 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 774 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 775 { 776 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 777 PetscErrorCode ierr; 778 PetscInt cnt = 0; 779 PC_FieldSplitLink ilink = jac->head; 780 781 PetscFunctionBegin; 782 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 783 while (ilink) { 784 (*subksp)[cnt++] = ilink->ksp; 785 ilink = ilink->next; 786 } 787 if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 788 *n = jac->nsplits; 789 PetscFunctionReturn(0); 790 } 791 EXTERN_C_END 792 793 EXTERN_C_BEGIN 794 #undef __FUNCT__ 795 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 796 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 797 { 798 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 799 PetscErrorCode ierr; 800 PC_FieldSplitLink ilink, next = jac->head; 801 char prefix[128]; 802 803 PetscFunctionBegin; 804 if (jac->splitdefined) { 805 ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 806 PetscFunctionReturn(0); 807 } 808 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 809 ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 810 ilink->is = is; 811 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 812 ilink->next = PETSC_NULL; 813 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 814 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 815 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 816 817 ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr); 818 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 819 820 if (!next) { 821 jac->head = ilink; 822 ilink->previous = PETSC_NULL; 823 } else { 824 while (next->next) { 825 next = next->next; 826 } 827 next->next = ilink; 828 ilink->previous = next; 829 } 830 jac->nsplits++; 831 832 PetscFunctionReturn(0); 833 } 834 EXTERN_C_END 835 836 #undef __FUNCT__ 837 #define __FUNCT__ "PCFieldSplitSetFields" 838 /*@ 839 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 840 841 Collective on PC 842 843 Input Parameters: 844 + pc - the preconditioner context 845 . splitname - name of this split 846 . n - the number of fields in this split 847 - fields - the fields in this split 848 849 Level: intermediate 850 851 Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 852 853 The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 854 size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean 855 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 856 where the numbered entries indicate what is in the field. 857 858 This function is called once per split (it creates a new split each time). Solve options 859 for this split will be available under the prefix -fieldsplit_SPLITNAME_. 860 861 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 862 863 @*/ 864 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 865 { 866 PetscErrorCode ierr,(*f)(PC,const char[],PetscInt,const PetscInt *); 867 868 PetscFunctionBegin; 869 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 870 PetscValidCharPointer(splitname,2); 871 if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 872 PetscValidIntPointer(fields,3); 873 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 874 if (f) { 875 ierr = (*f)(pc,splitname,n,fields);CHKERRQ(ierr); 876 } 877 PetscFunctionReturn(0); 878 } 879 880 #undef __FUNCT__ 881 #define __FUNCT__ "PCFieldSplitSetIS" 882 /*@ 883 PCFieldSplitSetIS - Sets the exact elements for field 884 885 Collective on PC 886 887 Input Parameters: 888 + pc - the preconditioner context 889 . splitname - name of this split 890 - is - the index set that defines the vector elements in this field 891 892 893 Notes: 894 Use PCFieldSplitSetFields(), for fields defined by strided types. 895 896 This function is called once per split (it creates a new split each time). Solve options 897 for this split will be available under the prefix -fieldsplit_SPLITNAME_. 898 899 Level: intermediate 900 901 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 902 903 @*/ 904 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 905 { 906 PetscErrorCode ierr,(*f)(PC,const char[],IS); 907 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 910 PetscValidCharPointer(splitname,2); 911 PetscValidHeaderSpecific(is,IS_CLASSID,3); 912 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 913 if (f) { 914 ierr = (*f)(pc,splitname,is);CHKERRQ(ierr); 915 } 916 PetscFunctionReturn(0); 917 } 918 919 #undef __FUNCT__ 920 #define __FUNCT__ "PCFieldSplitSetBlockSize" 921 /*@ 922 PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 923 fieldsplit preconditioner. If not set the matrix block size is used. 924 925 Collective on PC 926 927 Input Parameters: 928 + pc - the preconditioner context 929 - bs - the block size 930 931 Level: intermediate 932 933 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 934 935 @*/ 936 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 937 { 938 PetscErrorCode ierr,(*f)(PC,PetscInt); 939 940 PetscFunctionBegin; 941 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 942 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 943 if (f) { 944 ierr = (*f)(pc,bs);CHKERRQ(ierr); 945 } 946 PetscFunctionReturn(0); 947 } 948 949 #undef __FUNCT__ 950 #define __FUNCT__ "PCFieldSplitGetSubKSP" 951 /*@C 952 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 953 954 Collective on KSP 955 956 Input Parameter: 957 . pc - the preconditioner context 958 959 Output Parameters: 960 + n - the number of split 961 - pc - the array of KSP contexts 962 963 Note: 964 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 965 (not the KSP just the array that contains them). 966 967 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 968 969 Level: advanced 970 971 .seealso: PCFIELDSPLIT 972 @*/ 973 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 974 { 975 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 976 977 PetscFunctionBegin; 978 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 979 PetscValidIntPointer(n,2); 980 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 981 if (f) { 982 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 983 } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 984 PetscFunctionReturn(0); 985 } 986 987 #undef __FUNCT__ 988 #define __FUNCT__ "PCFieldSplitSchurPrecondition" 989 /*@ 990 PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 991 D matrix. Otherwise no preconditioner is used. 992 993 Collective on PC 994 995 Input Parameters: 996 + pc - the preconditioner context 997 . ptype - which matrix to use for preconditioning the Schur complement 998 - userpre - matrix to use for preconditioning, or PETSC_NULL 999 1000 Notes: 1001 The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 1002 There are currently no preconditioners that work directly with the Schur complement so setting 1003 PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 1004 1005 Options Database: 1006 . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1007 1008 Level: intermediate 1009 1010 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1011 1012 @*/ 1013 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1014 { 1015 PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 1016 1017 PetscFunctionBegin; 1018 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1019 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 1020 if (f) { 1021 ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1022 } 1023 PetscFunctionReturn(0); 1024 } 1025 1026 EXTERN_C_BEGIN 1027 #undef __FUNCT__ 1028 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1029 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1030 { 1031 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1032 PetscErrorCode ierr; 1033 1034 PetscFunctionBegin; 1035 jac->schurpre = ptype; 1036 if (pre) { 1037 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1038 jac->schur_user = pre; 1039 ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1040 } 1041 PetscFunctionReturn(0); 1042 } 1043 EXTERN_C_END 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 1047 /*@C 1048 PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 1049 1050 Collective on KSP 1051 1052 Input Parameter: 1053 . pc - the preconditioner context 1054 1055 Output Parameters: 1056 + A - the (0,0) block 1057 . B - the (0,1) block 1058 . C - the (1,0) block 1059 - D - the (1,1) block 1060 1061 Level: advanced 1062 1063 .seealso: PCFIELDSPLIT 1064 @*/ 1065 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 1066 { 1067 PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1068 1069 PetscFunctionBegin; 1070 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1071 if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1072 if (A) *A = jac->pmat[0]; 1073 if (B) *B = jac->B; 1074 if (C) *C = jac->C; 1075 if (D) *D = jac->pmat[1]; 1076 PetscFunctionReturn(0); 1077 } 1078 1079 EXTERN_C_BEGIN 1080 #undef __FUNCT__ 1081 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1082 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 1083 { 1084 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1085 PetscErrorCode ierr; 1086 1087 PetscFunctionBegin; 1088 jac->type = type; 1089 if (type == PC_COMPOSITE_SCHUR) { 1090 pc->ops->apply = PCApply_FieldSplit_Schur; 1091 pc->ops->view = PCView_FieldSplit_Schur; 1092 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1093 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1094 1095 } else { 1096 pc->ops->apply = PCApply_FieldSplit; 1097 pc->ops->view = PCView_FieldSplit; 1098 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1099 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1100 } 1101 PetscFunctionReturn(0); 1102 } 1103 EXTERN_C_END 1104 1105 EXTERN_C_BEGIN 1106 #undef __FUNCT__ 1107 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 1108 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 1109 { 1110 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1111 1112 PetscFunctionBegin; 1113 if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1114 if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 1115 jac->bs = bs; 1116 PetscFunctionReturn(0); 1117 } 1118 EXTERN_C_END 1119 1120 #undef __FUNCT__ 1121 #define __FUNCT__ "PCFieldSplitSetType" 1122 /*@ 1123 PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 1124 1125 Collective on PC 1126 1127 Input Parameter: 1128 . pc - the preconditioner context 1129 . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1130 1131 Options Database Key: 1132 . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 1133 1134 Level: Developer 1135 1136 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1137 1138 .seealso: PCCompositeSetType() 1139 1140 @*/ 1141 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 1142 { 1143 PetscErrorCode ierr,(*f)(PC,PCCompositeType); 1144 1145 PetscFunctionBegin; 1146 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1147 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 1148 if (f) { 1149 ierr = (*f)(pc,type);CHKERRQ(ierr); 1150 } 1151 PetscFunctionReturn(0); 1152 } 1153 1154 /* -------------------------------------------------------------------------------------*/ 1155 /*MC 1156 PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1157 fields or groups of fields 1158 1159 1160 To set options on the solvers for each block append -fieldsplit_ to all the PC 1161 options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 1162 1163 To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 1164 and set the options directly on the resulting KSP object 1165 1166 Level: intermediate 1167 1168 Options Database Keys: 1169 + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 1170 . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 1171 been supplied explicitly by -pc_fieldsplit_%d_fields 1172 . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 1173 . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1174 . -pc_fieldsplit_schur_precondition <true,false> default is true 1175 1176 - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 1177 for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 1178 1179 1180 Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1181 to define a field by an arbitrary collection of entries. 1182 1183 If no fields are set the default is used. The fields are defined by entries strided by bs, 1184 beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1185 if this is not called the block size defaults to the blocksize of the second matrix passed 1186 to KSPSetOperators()/PCSetOperators(). 1187 1188 Currently for the multiplicative version, the updated residual needed for the next field 1189 solve is computed via a matrix vector product over the entire array. An optimization would be 1190 to update the residual only for the part of the right hand side associated with the next field 1191 solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1192 part of the matrix needed to just update part of the residual). 1193 1194 For the Schur complement preconditioner if J = ( A B ) 1195 ( C D ) 1196 the preconditioner is 1197 (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1198 (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1199 where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1200 inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1201 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1202 D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1203 option. 1204 1205 If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1206 is used automatically for a second block. 1207 1208 Concepts: physics based preconditioners, block preconditioners 1209 1210 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1211 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1212 M*/ 1213 1214 EXTERN_C_BEGIN 1215 #undef __FUNCT__ 1216 #define __FUNCT__ "PCCreate_FieldSplit" 1217 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 1218 { 1219 PetscErrorCode ierr; 1220 PC_FieldSplit *jac; 1221 1222 PetscFunctionBegin; 1223 ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 1224 jac->bs = -1; 1225 jac->nsplits = 0; 1226 jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1227 jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1228 1229 pc->data = (void*)jac; 1230 1231 pc->ops->apply = PCApply_FieldSplit; 1232 pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 1233 pc->ops->setup = PCSetUp_FieldSplit; 1234 pc->ops->destroy = PCDestroy_FieldSplit; 1235 pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 1236 pc->ops->view = PCView_FieldSplit; 1237 pc->ops->applyrichardson = 0; 1238 1239 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 1240 PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1241 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 1242 PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1243 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1244 PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 1245 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 1246 PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 1247 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 1248 PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 1249 PetscFunctionReturn(0); 1250 } 1251 EXTERN_C_END 1252 1253 1254 /*MC 1255 Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1256 overview of these methods. 1257 1258 Consider the solution to ( A B ) (x_1) = (b_1) 1259 ( C D ) (x_2) (b_2) 1260 1261 Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1262 B' 0) (x_2) (b_2) 1263 1264 One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1265 for this block system. 1266 1267 Consider an additional matrix (Ap Bp) 1268 (Cp Dp) where some or all of the entries may be the same as 1269 in the original matrix (for example Ap == A). 1270 1271 In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1272 approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1273 1274 Block Jacobi: x_1 = A^ b_1 1275 x_2 = D^ b_2 1276 1277 Lower block Gauss-Seidel: x_1 = A^ b_1 1278 x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1279 1280 Symmetric Gauss-Seidel: x_1 = x_1 + A^(b_1 - A x_1 - B x_2) variant x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2) 1281 Interestingly this form is not actually a symmetric matrix, the symmetric version is 1282 x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1283 1284 Level: intermediate 1285 1286 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1287 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1288 M*/ 1289