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 PetscTruth realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 25 PetscInt bs; /* Block size for IS and Mat structures */ 26 PetscInt nsplits; /* Number of field divisions defined */ 27 Vec *x,*y,w1,w2; 28 Mat *mat; /* The diagonal block for each split */ 29 Mat *pmat; /* The preconditioning diagonal block for each split */ 30 Mat *Afield; /* The rows of the matrix associated with each split */ 31 PetscTruth issetup; 32 /* Only used when Schur complement preconditioning is used */ 33 Mat B; /* The (0,1) block */ 34 Mat C; /* The (1,0) block */ 35 Mat schur; /* The Schur complement S = D - C A^{-1} B */ 36 Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 37 PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 38 KSP kspschur; /* The solver for S */ 39 PC_FieldSplitLink head; 40 } PC_FieldSplit; 41 42 /* 43 Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 44 inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 45 PC you could change this. 46 */ 47 48 /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 49 * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 50 static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 51 { 52 switch (jac->schurpre) { 53 case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 54 case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 55 case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 56 default: 57 return jac->schur_user ? jac->schur_user : jac->pmat[1]; 58 } 59 } 60 61 62 #undef __FUNCT__ 63 #define __FUNCT__ "PCView_FieldSplit" 64 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 65 { 66 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 67 PetscErrorCode ierr; 68 PetscTruth iascii; 69 PetscInt i,j; 70 PC_FieldSplitLink ilink = jac->head; 71 72 PetscFunctionBegin; 73 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 74 if (iascii) { 75 ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 76 ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 77 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 78 for (i=0; i<jac->nsplits; i++) { 79 if (ilink->fields) { 80 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 81 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 82 for (j=0; j<ilink->nfields; j++) { 83 if (j > 0) { 84 ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 85 } 86 ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 87 } 88 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 89 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 90 } else { 91 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 92 } 93 ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 94 ilink = ilink->next; 95 } 96 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 97 } else { 98 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 99 } 100 PetscFunctionReturn(0); 101 } 102 103 #undef __FUNCT__ 104 #define __FUNCT__ "PCView_FieldSplit_Schur" 105 static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 106 { 107 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 108 PetscErrorCode ierr; 109 PetscTruth iascii; 110 PetscInt i,j; 111 PC_FieldSplitLink ilink = jac->head; 112 KSP ksp; 113 114 PetscFunctionBegin; 115 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 116 if (iascii) { 117 ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 118 ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 119 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 120 for (i=0; i<jac->nsplits; i++) { 121 if (ilink->fields) { 122 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 124 for (j=0; j<ilink->nfields; j++) { 125 if (j > 0) { 126 ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 127 } 128 ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 129 } 130 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 131 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 132 } else { 133 ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 134 } 135 ilink = ilink->next; 136 } 137 ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 138 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 139 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 140 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 141 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 142 ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 143 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 144 ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 145 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 146 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 147 } else { 148 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 149 } 150 PetscFunctionReturn(0); 151 } 152 153 #undef __FUNCT__ 154 #define __FUNCT__ "PCFieldSplitSetDefaults" 155 static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 156 { 157 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 158 PetscErrorCode ierr; 159 PC_FieldSplitLink ilink = jac->head; 160 PetscInt i = 0,*ifields,nfields; 161 PetscTruth flg = PETSC_FALSE,*fields,flg2; 162 char optionname[128]; 163 164 PetscFunctionBegin; 165 if (!ilink) { 166 167 if (jac->bs <= 0) { 168 if (pc->pmat) { 169 ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 170 } else { 171 jac->bs = 1; 172 } 173 } 174 175 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 176 if (!flg) { 177 /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 178 then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 179 flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 180 ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 181 while (PETSC_TRUE) { 182 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 183 nfields = jac->bs; 184 ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 185 if (!flg2) break; 186 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 187 flg = PETSC_FALSE; 188 ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 189 } 190 ierr = PetscFree(ifields);CHKERRQ(ierr); 191 } 192 193 if (flg) { 194 ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 195 ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 196 ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 197 while (ilink) { 198 for (i=0; i<ilink->nfields; i++) { 199 fields[ilink->fields[i]] = PETSC_TRUE; 200 } 201 ilink = ilink->next; 202 } 203 jac->defaultsplit = PETSC_TRUE; 204 for (i=0; i<jac->bs; i++) { 205 if (!fields[i]) { 206 ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 207 } else { 208 jac->defaultsplit = PETSC_FALSE; 209 } 210 } 211 ierr = PetscFree(fields);CHKERRQ(ierr); 212 } 213 } else if (jac->nsplits == 1) { 214 if (ilink->is) { 215 IS is2; 216 PetscInt nmin,nmax; 217 218 ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 219 ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 220 ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr); 221 ierr = ISDestroy(is2);CHKERRQ(ierr); 222 } else { 223 SETERRQ(PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 224 } 225 } 226 PetscFunctionReturn(0); 227 } 228 229 230 #undef __FUNCT__ 231 #define __FUNCT__ "PCSetUp_FieldSplit" 232 static PetscErrorCode PCSetUp_FieldSplit(PC pc) 233 { 234 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 235 PetscErrorCode ierr; 236 PC_FieldSplitLink ilink; 237 PetscInt i,nsplit,ccsize; 238 MatStructure flag = pc->flag; 239 PetscTruth sorted,getsub; 240 241 PetscFunctionBegin; 242 ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 243 nsplit = jac->nsplits; 244 ilink = jac->head; 245 246 /* get the matrices for each split */ 247 if (!jac->issetup) { 248 PetscInt rstart,rend,nslots,bs; 249 250 jac->issetup = PETSC_TRUE; 251 252 /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 253 bs = jac->bs; 254 ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 255 ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 256 nslots = (rend - rstart)/bs; 257 for (i=0; i<nsplit; i++) { 258 if (jac->defaultsplit) { 259 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 260 } else if (!ilink->is) { 261 if (ilink->nfields > 1) { 262 PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 263 ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 264 for (j=0; j<nslots; j++) { 265 for (k=0; k<nfields; k++) { 266 ii[nfields*j + k] = rstart + bs*j + fields[k]; 267 } 268 } 269 ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 270 ierr = PetscFree(ii);CHKERRQ(ierr); 271 } else { 272 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 273 } 274 } 275 ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 276 if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 277 ilink = ilink->next; 278 } 279 } 280 281 ilink = jac->head; 282 if (!jac->pmat) { 283 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 284 for (i=0; i<nsplit; i++) { 285 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 286 ilink = ilink->next; 287 } 288 } else { 289 for (i=0; i<nsplit; i++) { 290 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 291 ilink = ilink->next; 292 } 293 } 294 if (jac->realdiagonal) { 295 ilink = jac->head; 296 if (!jac->mat) { 297 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 298 for (i=0; i<nsplit; i++) { 299 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 300 ilink = ilink->next; 301 } 302 } else { 303 for (i=0; i<nsplit; i++) { 304 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 305 ilink = ilink->next; 306 } 307 } 308 } else { 309 jac->mat = jac->pmat; 310 } 311 312 if (jac->type != PC_COMPOSITE_SCHUR) { 313 /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 314 ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 315 if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 316 ilink = jac->head; 317 if (!jac->Afield) { 318 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 319 for (i=0; i<nsplit; i++) { 320 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 321 ilink = ilink->next; 322 } 323 } else { 324 for (i=0; i<nsplit; i++) { 325 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 326 ilink = ilink->next; 327 } 328 } 329 } 330 } 331 332 if (jac->type == PC_COMPOSITE_SCHUR) { 333 IS ccis; 334 PetscInt rstart,rend; 335 if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 336 337 /* When extracting off-diagonal submatrices, we take complements from this range */ 338 ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 339 340 /* need to handle case when one is resetting up the preconditioner */ 341 if (jac->schur) { 342 ilink = jac->head; 343 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 344 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 345 ierr = ISDestroy(ccis);CHKERRQ(ierr); 346 ilink = ilink->next; 347 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 348 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 349 ierr = ISDestroy(ccis);CHKERRQ(ierr); 350 ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 351 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 352 353 } else { 354 KSP ksp; 355 356 /* extract the B and C matrices */ 357 ilink = jac->head; 358 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 359 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 360 ierr = ISDestroy(ccis);CHKERRQ(ierr); 361 ilink = ilink->next; 362 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 363 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 364 ierr = ISDestroy(ccis);CHKERRQ(ierr); 365 /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 366 ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 367 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 368 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 369 ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 370 371 ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 372 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 373 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 374 if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 375 PC pc; 376 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 377 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 378 /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 379 } 380 ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 381 ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 382 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 383 ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 384 385 ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 386 ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 387 ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 388 ilink = jac->head; 389 ilink->x = jac->x[0]; ilink->y = jac->y[0]; 390 ilink = ilink->next; 391 ilink->x = jac->x[1]; ilink->y = jac->y[1]; 392 } 393 } else { 394 /* set up the individual PCs */ 395 i = 0; 396 ilink = jac->head; 397 while (ilink) { 398 ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 399 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 400 ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 401 ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 402 i++; 403 ilink = ilink->next; 404 } 405 406 /* create work vectors for each split */ 407 if (!jac->x) { 408 ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 409 ilink = jac->head; 410 for (i=0; i<nsplit; i++) { 411 Vec *vl,*vr; 412 413 ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 414 ilink->x = *vr; 415 ilink->y = *vl; 416 ierr = PetscFree(vr);CHKERRQ(ierr); 417 ierr = PetscFree(vl);CHKERRQ(ierr); 418 jac->x[i] = ilink->x; 419 jac->y[i] = ilink->y; 420 ilink = ilink->next; 421 } 422 } 423 } 424 425 426 if (!jac->head->sctx) { 427 Vec xtmp; 428 429 /* compute scatter contexts needed by multiplicative versions and non-default splits */ 430 431 ilink = jac->head; 432 ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 433 for (i=0; i<nsplit; i++) { 434 ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 435 ilink = ilink->next; 436 } 437 ierr = VecDestroy(xtmp);CHKERRQ(ierr); 438 } 439 PetscFunctionReturn(0); 440 } 441 442 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 443 (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 444 VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 445 KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 446 VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 447 VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 448 449 #undef __FUNCT__ 450 #define __FUNCT__ "PCApply_FieldSplit_Schur" 451 static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 452 { 453 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 454 PetscErrorCode ierr; 455 KSP ksp; 456 PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 457 458 PetscFunctionBegin; 459 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 460 461 ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 462 ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 463 ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 464 ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 465 ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 466 ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 467 ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 468 469 ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 470 ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 471 ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 472 473 ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 474 ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 475 ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 476 ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 477 ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 478 479 PetscFunctionReturn(0); 480 } 481 482 #undef __FUNCT__ 483 #define __FUNCT__ "PCApply_FieldSplit" 484 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 485 { 486 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 487 PetscErrorCode ierr; 488 PC_FieldSplitLink ilink = jac->head; 489 PetscInt bs,cnt; 490 491 PetscFunctionBegin; 492 CHKMEMQ; 493 ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 494 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 495 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 496 497 if (jac->type == PC_COMPOSITE_ADDITIVE) { 498 if (jac->defaultsplit) { 499 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 500 while (ilink) { 501 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 502 ilink = ilink->next; 503 } 504 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 505 } else { 506 ierr = VecSet(y,0.0);CHKERRQ(ierr); 507 while (ilink) { 508 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 509 ilink = ilink->next; 510 } 511 } 512 } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 513 if (!jac->w1) { 514 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 515 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 516 } 517 ierr = VecSet(y,0.0);CHKERRQ(ierr); 518 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 519 cnt = 1; 520 while (ilink->next) { 521 ilink = ilink->next; 522 if (jac->Afield) { 523 /* compute the residual only over the part of the vector needed */ 524 ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 525 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 526 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 527 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 528 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 529 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 530 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 531 } else { 532 /* compute the residual over the entire vector */ 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 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 539 cnt -= 2; 540 while (ilink->previous) { 541 ilink = ilink->previous; 542 if (jac->Afield) { 543 /* compute the residual only over the part of the vector needed */ 544 ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 545 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 546 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 547 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 548 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 549 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 550 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 551 } else { 552 ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 553 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 554 ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 555 } 556 } 557 } 558 } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 559 CHKMEMQ; 560 PetscFunctionReturn(0); 561 } 562 563 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 564 (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 565 VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 566 KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 567 VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 568 VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 569 570 #undef __FUNCT__ 571 #define __FUNCT__ "PCApply_FieldSplit" 572 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 573 { 574 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 575 PetscErrorCode ierr; 576 PC_FieldSplitLink ilink = jac->head; 577 PetscInt bs; 578 579 PetscFunctionBegin; 580 CHKMEMQ; 581 ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 582 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 583 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 584 585 if (jac->type == PC_COMPOSITE_ADDITIVE) { 586 if (jac->defaultsplit) { 587 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 588 while (ilink) { 589 ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 590 ilink = ilink->next; 591 } 592 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 593 } else { 594 ierr = VecSet(y,0.0);CHKERRQ(ierr); 595 while (ilink) { 596 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 597 ilink = ilink->next; 598 } 599 } 600 } else { 601 if (!jac->w1) { 602 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 603 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 604 } 605 ierr = VecSet(y,0.0);CHKERRQ(ierr); 606 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 607 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 608 while (ilink->next) { 609 ilink = ilink->next; 610 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 611 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 612 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 613 } 614 while (ilink->previous) { 615 ilink = ilink->previous; 616 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 617 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 618 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 619 } 620 } else { 621 while (ilink->next) { /* get to last entry in linked list */ 622 ilink = ilink->next; 623 } 624 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 625 while (ilink->previous) { 626 ilink = ilink->previous; 627 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 628 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 629 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 630 } 631 } 632 } 633 CHKMEMQ; 634 PetscFunctionReturn(0); 635 } 636 637 #undef __FUNCT__ 638 #define __FUNCT__ "PCDestroy_FieldSplit" 639 static PetscErrorCode PCDestroy_FieldSplit(PC pc) 640 { 641 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 642 PetscErrorCode ierr; 643 PC_FieldSplitLink ilink = jac->head,next; 644 645 PetscFunctionBegin; 646 while (ilink) { 647 ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 648 if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 649 if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 650 if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 651 if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 652 next = ilink->next; 653 ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 654 ierr = PetscFree(ilink);CHKERRQ(ierr); 655 ilink = next; 656 } 657 ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 658 if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 659 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 660 if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 661 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 662 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 663 if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 664 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 665 if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 666 if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 667 if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 668 ierr = PetscFree(jac);CHKERRQ(ierr); 669 PetscFunctionReturn(0); 670 } 671 672 #undef __FUNCT__ 673 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 674 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 675 { 676 PetscErrorCode ierr; 677 PetscInt i = 0,nfields,*fields,bs; 678 PetscTruth flg; 679 char optionname[128]; 680 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 681 PCCompositeType ctype; 682 683 PetscFunctionBegin; 684 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 685 ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 686 ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 687 if (flg) { 688 ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 689 } 690 691 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 692 if (flg) { 693 ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 694 } 695 696 /* Only setup fields once */ 697 if ((jac->bs > 0) && (jac->nsplits == 0)) { 698 /* only allow user to set fields from command line if bs is already known. 699 otherwise user can set them in PCFieldSplitSetDefaults() */ 700 ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 701 while (PETSC_TRUE) { 702 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 703 nfields = jac->bs; 704 ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 705 if (!flg) break; 706 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 707 ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 708 } 709 ierr = PetscFree(fields);CHKERRQ(ierr); 710 } 711 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); 712 ierr = PetscOptionsTail();CHKERRQ(ierr); 713 PetscFunctionReturn(0); 714 } 715 716 /*------------------------------------------------------------------------------------*/ 717 718 EXTERN_C_BEGIN 719 #undef __FUNCT__ 720 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 721 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 722 { 723 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 724 PetscErrorCode ierr; 725 PC_FieldSplitLink ilink,next = jac->head; 726 char prefix[128]; 727 PetscInt i; 728 729 PetscFunctionBegin; 730 if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 731 for (i=0; i<n; i++) { 732 if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 733 if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 734 } 735 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 736 ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 737 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 738 ilink->nfields = n; 739 ilink->next = PETSC_NULL; 740 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 741 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 742 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 743 744 if (((PetscObject)pc)->prefix) { 745 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 746 } else { 747 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 748 } 749 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 750 751 if (!next) { 752 jac->head = ilink; 753 ilink->previous = PETSC_NULL; 754 } else { 755 while (next->next) { 756 next = next->next; 757 } 758 next->next = ilink; 759 ilink->previous = next; 760 } 761 jac->nsplits++; 762 PetscFunctionReturn(0); 763 } 764 EXTERN_C_END 765 766 EXTERN_C_BEGIN 767 #undef __FUNCT__ 768 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 769 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 770 { 771 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 772 PetscErrorCode ierr; 773 774 PetscFunctionBegin; 775 ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 776 ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 777 (*subksp)[1] = jac->kspschur; 778 *n = jac->nsplits; 779 PetscFunctionReturn(0); 780 } 781 EXTERN_C_END 782 783 EXTERN_C_BEGIN 784 #undef __FUNCT__ 785 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 786 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 787 { 788 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 789 PetscErrorCode ierr; 790 PetscInt cnt = 0; 791 PC_FieldSplitLink ilink = jac->head; 792 793 PetscFunctionBegin; 794 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 795 while (ilink) { 796 (*subksp)[cnt++] = ilink->ksp; 797 ilink = ilink->next; 798 } 799 if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 800 *n = jac->nsplits; 801 PetscFunctionReturn(0); 802 } 803 EXTERN_C_END 804 805 EXTERN_C_BEGIN 806 #undef __FUNCT__ 807 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 808 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 809 { 810 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 811 PetscErrorCode ierr; 812 PC_FieldSplitLink ilink, next = jac->head; 813 char prefix[128]; 814 815 PetscFunctionBegin; 816 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 817 ilink->is = is; 818 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 819 ilink->next = PETSC_NULL; 820 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 821 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 822 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 823 824 if (((PetscObject)pc)->prefix) { 825 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 826 } else { 827 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 828 } 829 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 830 831 if (!next) { 832 jac->head = ilink; 833 ilink->previous = PETSC_NULL; 834 } else { 835 while (next->next) { 836 next = next->next; 837 } 838 next->next = ilink; 839 ilink->previous = next; 840 } 841 jac->nsplits++; 842 843 PetscFunctionReturn(0); 844 } 845 EXTERN_C_END 846 847 #undef __FUNCT__ 848 #define __FUNCT__ "PCFieldSplitSetFields" 849 /*@ 850 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 851 852 Collective on PC 853 854 Input Parameters: 855 + pc - the preconditioner context 856 . n - the number of fields in this split 857 . fields - the fields in this split 858 859 Level: intermediate 860 861 Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 862 863 The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 864 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 865 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 866 where the numbered entries indicate what is in the field. 867 868 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 869 870 @*/ 871 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 872 { 873 PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 874 875 PetscFunctionBegin; 876 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 877 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 878 if (f) { 879 ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 880 } 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "PCFieldSplitSetIS" 886 /*@ 887 PCFieldSplitSetIS - Sets the exact elements for field 888 889 Collective on PC 890 891 Input Parameters: 892 + pc - the preconditioner context 893 . is - the index set that defines the vector elements in this field 894 895 896 Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 897 898 Level: intermediate 899 900 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 901 902 @*/ 903 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 904 { 905 PetscErrorCode ierr,(*f)(PC,IS); 906 907 PetscFunctionBegin; 908 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 909 PetscValidHeaderSpecific(is,IS_COOKIE,1); 910 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 911 if (f) { 912 ierr = (*f)(pc,is);CHKERRQ(ierr); 913 } 914 PetscFunctionReturn(0); 915 } 916 917 #undef __FUNCT__ 918 #define __FUNCT__ "PCFieldSplitSetBlockSize" 919 /*@ 920 PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 921 fieldsplit preconditioner. If not set the matrix block size is used. 922 923 Collective on PC 924 925 Input Parameters: 926 + pc - the preconditioner context 927 - bs - the block size 928 929 Level: intermediate 930 931 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 932 933 @*/ 934 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 935 { 936 PetscErrorCode ierr,(*f)(PC,PetscInt); 937 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 940 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 941 if (f) { 942 ierr = (*f)(pc,bs);CHKERRQ(ierr); 943 } 944 PetscFunctionReturn(0); 945 } 946 947 #undef __FUNCT__ 948 #define __FUNCT__ "PCFieldSplitGetSubKSP" 949 /*@C 950 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 951 952 Collective on KSP 953 954 Input Parameter: 955 . pc - the preconditioner context 956 957 Output Parameters: 958 + n - the number of split 959 - pc - the array of KSP contexts 960 961 Note: 962 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 963 (not the KSP just the array that contains them). 964 965 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 966 967 Level: advanced 968 969 .seealso: PCFIELDSPLIT 970 @*/ 971 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 972 { 973 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 974 975 PetscFunctionBegin; 976 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 977 PetscValidIntPointer(n,2); 978 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 979 if (f) { 980 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 981 } else { 982 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 983 } 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_COOKIE,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_COOKIE,1); 1071 if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(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(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1114 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); 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_COOKIE,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