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; 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_ADDITIVE && 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 ilink = jac->head; 315 if (!jac->Afield) { 316 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 317 for (i=0; i<nsplit; i++) { 318 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 319 ilink = ilink->next; 320 } 321 } else { 322 for (i=0; i<nsplit; i++) { 323 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 324 ilink = ilink->next; 325 } 326 } 327 } 328 329 if (jac->type == PC_COMPOSITE_SCHUR) { 330 IS ccis; 331 PetscInt rstart,rend; 332 if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 333 334 /* When extracting off-diagonal submatrices, we take complements from this range */ 335 ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 336 337 /* need to handle case when one is resetting up the preconditioner */ 338 if (jac->schur) { 339 ilink = jac->head; 340 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 341 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 342 ierr = ISDestroy(ccis);CHKERRQ(ierr); 343 ilink = ilink->next; 344 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 345 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 346 ierr = ISDestroy(ccis);CHKERRQ(ierr); 347 ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 348 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 349 350 } else { 351 KSP ksp; 352 353 /* extract the B and C matrices */ 354 ilink = jac->head; 355 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 356 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 357 ierr = ISDestroy(ccis);CHKERRQ(ierr); 358 ilink = ilink->next; 359 ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 360 ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 361 ierr = ISDestroy(ccis);CHKERRQ(ierr); 362 /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 363 ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 364 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 365 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 366 ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 367 368 ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 369 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 370 ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 371 if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 372 PC pc; 373 ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 374 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 375 /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 376 } 377 ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 378 ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 379 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 380 ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 381 382 ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 383 ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 384 ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 385 ilink = jac->head; 386 ilink->x = jac->x[0]; ilink->y = jac->y[0]; 387 ilink = ilink->next; 388 ilink->x = jac->x[1]; ilink->y = jac->y[1]; 389 } 390 } else { 391 /* set up the individual PCs */ 392 i = 0; 393 ilink = jac->head; 394 while (ilink) { 395 ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 396 /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 397 ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 398 ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 399 i++; 400 ilink = ilink->next; 401 } 402 403 /* create work vectors for each split */ 404 if (!jac->x) { 405 ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 406 ilink = jac->head; 407 for (i=0; i<nsplit; i++) { 408 Vec *vl,*vr; 409 410 ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 411 ilink->x = *vr; 412 ilink->y = *vl; 413 ierr = PetscFree(vr);CHKERRQ(ierr); 414 ierr = PetscFree(vl);CHKERRQ(ierr); 415 jac->x[i] = ilink->x; 416 jac->y[i] = ilink->y; 417 ilink = ilink->next; 418 } 419 } 420 } 421 422 423 if (!jac->head->sctx) { 424 Vec xtmp; 425 426 /* compute scatter contexts needed by multiplicative versions and non-default splits */ 427 428 ilink = jac->head; 429 ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 430 for (i=0; i<nsplit; i++) { 431 ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 432 ilink = ilink->next; 433 } 434 ierr = VecDestroy(xtmp);CHKERRQ(ierr); 435 } 436 PetscFunctionReturn(0); 437 } 438 439 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 440 (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 441 VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 442 KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 443 VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 444 VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 445 446 #undef __FUNCT__ 447 #define __FUNCT__ "PCApply_FieldSplit_Schur" 448 static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 449 { 450 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 451 PetscErrorCode ierr; 452 KSP ksp; 453 PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 454 455 PetscFunctionBegin; 456 ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 457 458 ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 459 ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 460 ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 461 ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 462 ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 463 ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 464 ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 465 466 ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 467 ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 468 ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 469 470 ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 471 ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 472 ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 473 ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 474 ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 475 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "PCApply_FieldSplit" 481 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 482 { 483 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 484 PetscErrorCode ierr; 485 PC_FieldSplitLink ilink = jac->head; 486 PetscInt cnt; 487 488 PetscFunctionBegin; 489 CHKMEMQ; 490 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 491 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 492 493 if (jac->type == PC_COMPOSITE_ADDITIVE) { 494 if (jac->defaultsplit) { 495 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 496 while (ilink) { 497 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 498 ilink = ilink->next; 499 } 500 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 501 } else { 502 ierr = VecSet(y,0.0);CHKERRQ(ierr); 503 while (ilink) { 504 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 505 ilink = ilink->next; 506 } 507 } 508 } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 509 if (!jac->w1) { 510 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 511 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 512 } 513 ierr = VecSet(y,0.0);CHKERRQ(ierr); 514 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 515 cnt = 1; 516 while (ilink->next) { 517 ilink = ilink->next; 518 /* compute the residual only over the part of the vector needed */ 519 ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 520 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 521 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 522 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 523 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 524 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 525 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 526 } 527 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 528 cnt -= 2; 529 while (ilink->previous) { 530 ilink = ilink->previous; 531 /* compute the residual only over the part of the vector needed */ 532 ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 533 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 534 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 535 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 536 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 537 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 538 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 539 } 540 } 541 } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 542 CHKMEMQ; 543 PetscFunctionReturn(0); 544 } 545 546 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 547 (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 548 VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 549 KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 550 VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 551 VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 552 553 #undef __FUNCT__ 554 #define __FUNCT__ "PCApply_FieldSplit" 555 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 556 { 557 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 558 PetscErrorCode ierr; 559 PC_FieldSplitLink ilink = jac->head; 560 561 PetscFunctionBegin; 562 CHKMEMQ; 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->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 640 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 641 if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 642 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 643 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 644 if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 645 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 646 if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 647 if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 648 if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 649 ierr = PetscFree(jac);CHKERRQ(ierr); 650 PetscFunctionReturn(0); 651 } 652 653 #undef __FUNCT__ 654 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 655 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 656 { 657 PetscErrorCode ierr; 658 PetscInt i = 0,nfields,*fields,bs; 659 PetscTruth flg; 660 char optionname[128]; 661 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 662 PCCompositeType ctype; 663 664 PetscFunctionBegin; 665 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 666 ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 667 ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 668 if (flg) { 669 ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 670 } 671 672 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 673 if (flg) { 674 ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 675 } 676 677 /* Only setup fields once */ 678 if ((jac->bs > 0) && (jac->nsplits == 0)) { 679 /* only allow user to set fields from command line if bs is already known. 680 otherwise user can set them in PCFieldSplitSetDefaults() */ 681 ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 682 while (PETSC_TRUE) { 683 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 684 nfields = jac->bs; 685 ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 686 if (!flg) break; 687 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 688 ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 689 } 690 ierr = PetscFree(fields);CHKERRQ(ierr); 691 } 692 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); 693 ierr = PetscOptionsTail();CHKERRQ(ierr); 694 PetscFunctionReturn(0); 695 } 696 697 /*------------------------------------------------------------------------------------*/ 698 699 EXTERN_C_BEGIN 700 #undef __FUNCT__ 701 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 702 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 703 { 704 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 705 PetscErrorCode ierr; 706 PC_FieldSplitLink ilink,next = jac->head; 707 char prefix[128]; 708 PetscInt i; 709 710 PetscFunctionBegin; 711 if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 712 for (i=0; i<n; i++) { 713 if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 714 if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 715 } 716 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 717 ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 718 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 719 ilink->nfields = n; 720 ilink->next = PETSC_NULL; 721 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 722 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 723 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 724 725 if (((PetscObject)pc)->prefix) { 726 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 727 } else { 728 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 729 } 730 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 731 732 if (!next) { 733 jac->head = ilink; 734 ilink->previous = PETSC_NULL; 735 } else { 736 while (next->next) { 737 next = next->next; 738 } 739 next->next = ilink; 740 ilink->previous = next; 741 } 742 jac->nsplits++; 743 PetscFunctionReturn(0); 744 } 745 EXTERN_C_END 746 747 EXTERN_C_BEGIN 748 #undef __FUNCT__ 749 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 750 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 751 { 752 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 753 PetscErrorCode ierr; 754 755 PetscFunctionBegin; 756 ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 757 ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 758 (*subksp)[1] = jac->kspschur; 759 *n = jac->nsplits; 760 PetscFunctionReturn(0); 761 } 762 EXTERN_C_END 763 764 EXTERN_C_BEGIN 765 #undef __FUNCT__ 766 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 767 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 768 { 769 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 770 PetscErrorCode ierr; 771 PetscInt cnt = 0; 772 PC_FieldSplitLink ilink = jac->head; 773 774 PetscFunctionBegin; 775 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 776 while (ilink) { 777 (*subksp)[cnt++] = ilink->ksp; 778 ilink = ilink->next; 779 } 780 if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 781 *n = jac->nsplits; 782 PetscFunctionReturn(0); 783 } 784 EXTERN_C_END 785 786 EXTERN_C_BEGIN 787 #undef __FUNCT__ 788 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 789 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 790 { 791 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 792 PetscErrorCode ierr; 793 PC_FieldSplitLink ilink, next = jac->head; 794 char prefix[128]; 795 796 PetscFunctionBegin; 797 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 798 ilink->is = is; 799 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 800 ilink->next = PETSC_NULL; 801 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 802 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 803 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 804 805 if (((PetscObject)pc)->prefix) { 806 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 807 } else { 808 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 809 } 810 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 811 812 if (!next) { 813 jac->head = ilink; 814 ilink->previous = PETSC_NULL; 815 } else { 816 while (next->next) { 817 next = next->next; 818 } 819 next->next = ilink; 820 ilink->previous = next; 821 } 822 jac->nsplits++; 823 824 PetscFunctionReturn(0); 825 } 826 EXTERN_C_END 827 828 #undef __FUNCT__ 829 #define __FUNCT__ "PCFieldSplitSetFields" 830 /*@ 831 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 832 833 Collective on PC 834 835 Input Parameters: 836 + pc - the preconditioner context 837 . n - the number of fields in this split 838 . fields - the fields in this split 839 840 Level: intermediate 841 842 Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 843 844 The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 845 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 846 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 847 where the numbered entries indicate what is in the field. 848 849 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 850 851 @*/ 852 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 853 { 854 PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 855 856 PetscFunctionBegin; 857 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 858 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 859 if (f) { 860 ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 861 } 862 PetscFunctionReturn(0); 863 } 864 865 #undef __FUNCT__ 866 #define __FUNCT__ "PCFieldSplitSetIS" 867 /*@ 868 PCFieldSplitSetIS - Sets the exact elements for field 869 870 Collective on PC 871 872 Input Parameters: 873 + pc - the preconditioner context 874 . is - the index set that defines the vector elements in this field 875 876 877 Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 878 879 Level: intermediate 880 881 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 882 883 @*/ 884 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 885 { 886 PetscErrorCode ierr,(*f)(PC,IS); 887 888 PetscFunctionBegin; 889 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 890 PetscValidHeaderSpecific(is,IS_COOKIE,1); 891 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 892 if (f) { 893 ierr = (*f)(pc,is);CHKERRQ(ierr); 894 } 895 PetscFunctionReturn(0); 896 } 897 898 #undef __FUNCT__ 899 #define __FUNCT__ "PCFieldSplitSetBlockSize" 900 /*@ 901 PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 902 fieldsplit preconditioner. If not set the matrix block size is used. 903 904 Collective on PC 905 906 Input Parameters: 907 + pc - the preconditioner context 908 - bs - the block size 909 910 Level: intermediate 911 912 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 913 914 @*/ 915 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 916 { 917 PetscErrorCode ierr,(*f)(PC,PetscInt); 918 919 PetscFunctionBegin; 920 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 921 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 922 if (f) { 923 ierr = (*f)(pc,bs);CHKERRQ(ierr); 924 } 925 PetscFunctionReturn(0); 926 } 927 928 #undef __FUNCT__ 929 #define __FUNCT__ "PCFieldSplitGetSubKSP" 930 /*@C 931 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 932 933 Collective on KSP 934 935 Input Parameter: 936 . pc - the preconditioner context 937 938 Output Parameters: 939 + n - the number of split 940 - pc - the array of KSP contexts 941 942 Note: 943 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 944 (not the KSP just the array that contains them). 945 946 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 947 948 Level: advanced 949 950 .seealso: PCFIELDSPLIT 951 @*/ 952 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 953 { 954 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 955 956 PetscFunctionBegin; 957 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 958 PetscValidIntPointer(n,2); 959 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 960 if (f) { 961 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 962 } else { 963 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 964 } 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "PCFieldSplitSchurPrecondition" 970 /*@ 971 PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 972 D matrix. Otherwise no preconditioner is used. 973 974 Collective on PC 975 976 Input Parameters: 977 + pc - the preconditioner context 978 . ptype - which matrix to use for preconditioning the Schur complement 979 - userpre - matrix to use for preconditioning, or PETSC_NULL 980 981 Notes: 982 The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 983 There are currently no preconditioners that work directly with the Schur complement so setting 984 PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 985 986 Options Database: 987 . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 988 989 Level: intermediate 990 991 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 992 993 @*/ 994 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 995 { 996 PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 997 998 PetscFunctionBegin; 999 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1000 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 1001 if (f) { 1002 ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1003 } 1004 PetscFunctionReturn(0); 1005 } 1006 1007 EXTERN_C_BEGIN 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1010 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1011 { 1012 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1013 PetscErrorCode ierr; 1014 1015 PetscFunctionBegin; 1016 jac->schurpre = ptype; 1017 if (pre) { 1018 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1019 jac->schur_user = pre; 1020 ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1021 } 1022 PetscFunctionReturn(0); 1023 } 1024 EXTERN_C_END 1025 1026 #undef __FUNCT__ 1027 #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 1028 /*@C 1029 PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 1030 1031 Collective on KSP 1032 1033 Input Parameter: 1034 . pc - the preconditioner context 1035 1036 Output Parameters: 1037 + A - the (0,0) block 1038 . B - the (0,1) block 1039 . C - the (1,0) block 1040 - D - the (1,1) block 1041 1042 Level: advanced 1043 1044 .seealso: PCFIELDSPLIT 1045 @*/ 1046 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 1047 { 1048 PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1049 1050 PetscFunctionBegin; 1051 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1052 if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");} 1053 if (A) *A = jac->pmat[0]; 1054 if (B) *B = jac->B; 1055 if (C) *C = jac->C; 1056 if (D) *D = jac->pmat[1]; 1057 PetscFunctionReturn(0); 1058 } 1059 1060 EXTERN_C_BEGIN 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1063 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 1064 { 1065 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1066 PetscErrorCode ierr; 1067 1068 PetscFunctionBegin; 1069 jac->type = type; 1070 if (type == PC_COMPOSITE_SCHUR) { 1071 pc->ops->apply = PCApply_FieldSplit_Schur; 1072 pc->ops->view = PCView_FieldSplit_Schur; 1073 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1074 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1075 1076 } else { 1077 pc->ops->apply = PCApply_FieldSplit; 1078 pc->ops->view = PCView_FieldSplit; 1079 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1080 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1081 } 1082 PetscFunctionReturn(0); 1083 } 1084 EXTERN_C_END 1085 1086 EXTERN_C_BEGIN 1087 #undef __FUNCT__ 1088 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 1089 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 1090 { 1091 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1092 1093 PetscFunctionBegin; 1094 if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1095 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); 1096 jac->bs = bs; 1097 PetscFunctionReturn(0); 1098 } 1099 EXTERN_C_END 1100 1101 #undef __FUNCT__ 1102 #define __FUNCT__ "PCFieldSplitSetType" 1103 /*@ 1104 PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 1105 1106 Collective on PC 1107 1108 Input Parameter: 1109 . pc - the preconditioner context 1110 . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1111 1112 Options Database Key: 1113 . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 1114 1115 Level: Developer 1116 1117 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1118 1119 .seealso: PCCompositeSetType() 1120 1121 @*/ 1122 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 1123 { 1124 PetscErrorCode ierr,(*f)(PC,PCCompositeType); 1125 1126 PetscFunctionBegin; 1127 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1128 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 1129 if (f) { 1130 ierr = (*f)(pc,type);CHKERRQ(ierr); 1131 } 1132 PetscFunctionReturn(0); 1133 } 1134 1135 /* -------------------------------------------------------------------------------------*/ 1136 /*MC 1137 PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1138 fields or groups of fields 1139 1140 1141 To set options on the solvers for each block append -fieldsplit_ to all the PC 1142 options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 1143 1144 To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 1145 and set the options directly on the resulting KSP object 1146 1147 Level: intermediate 1148 1149 Options Database Keys: 1150 + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 1151 . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 1152 been supplied explicitly by -pc_fieldsplit_%d_fields 1153 . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 1154 . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1155 . -pc_fieldsplit_schur_precondition <true,false> default is true 1156 1157 - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 1158 for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 1159 1160 1161 Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1162 to define a field by an arbitrary collection of entries. 1163 1164 If no fields are set the default is used. The fields are defined by entries strided by bs, 1165 beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1166 if this is not called the block size defaults to the blocksize of the second matrix passed 1167 to KSPSetOperators()/PCSetOperators(). 1168 1169 Currently for the multiplicative version, the updated residual needed for the next field 1170 solve is computed via a matrix vector product over the entire array. An optimization would be 1171 to update the residual only for the part of the right hand side associated with the next field 1172 solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1173 part of the matrix needed to just update part of the residual). 1174 1175 For the Schur complement preconditioner if J = ( A B ) 1176 ( C D ) 1177 the preconditioner is 1178 (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1179 (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1180 where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1181 inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1182 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1183 D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1184 option. 1185 1186 If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1187 is used automatically for a second block. 1188 1189 Concepts: physics based preconditioners, block preconditioners 1190 1191 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1192 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1193 M*/ 1194 1195 EXTERN_C_BEGIN 1196 #undef __FUNCT__ 1197 #define __FUNCT__ "PCCreate_FieldSplit" 1198 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 1199 { 1200 PetscErrorCode ierr; 1201 PC_FieldSplit *jac; 1202 1203 PetscFunctionBegin; 1204 ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 1205 jac->bs = -1; 1206 jac->nsplits = 0; 1207 jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1208 jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1209 1210 pc->data = (void*)jac; 1211 1212 pc->ops->apply = PCApply_FieldSplit; 1213 pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 1214 pc->ops->setup = PCSetUp_FieldSplit; 1215 pc->ops->destroy = PCDestroy_FieldSplit; 1216 pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 1217 pc->ops->view = PCView_FieldSplit; 1218 pc->ops->applyrichardson = 0; 1219 1220 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 1221 PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1222 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 1223 PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1224 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1225 PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 1226 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 1227 PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 1228 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 1229 PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 1230 PetscFunctionReturn(0); 1231 } 1232 EXTERN_C_END 1233 1234 1235 /*MC 1236 Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1237 overview of these methods. 1238 1239 Consider the solution to ( A B ) (x_1) = (b_1) 1240 ( C D ) (x_2) (b_2) 1241 1242 Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1243 B' 0) (x_2) (b_2) 1244 1245 One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1246 for this block system. 1247 1248 Consider an additional matrix (Ap Bp) 1249 (Cp Dp) where some or all of the entries may be the same as 1250 in the original matrix (for example Ap == A). 1251 1252 In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1253 approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1254 1255 Block Jacobi: x_1 = A^ b_1 1256 x_2 = D^ b_2 1257 1258 Lower block Gauss-Seidel: x_1 = A^ b_1 1259 x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1260 1261 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) 1262 Interestingly this form is not actually a symmetric matrix, the symmetric version is 1263 x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1264 1265 Level: intermediate 1266 1267 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1268 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1269 M*/ 1270