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