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 if (jac->nsplits < 2) { 227 SETERRQ(PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 228 } 229 PetscFunctionReturn(0); 230 } 231 232 233 #undef __FUNCT__ 234 #define __FUNCT__ "PCSetUp_FieldSplit" 235 static PetscErrorCode PCSetUp_FieldSplit(PC pc) 236 { 237 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 238 PetscErrorCode ierr; 239 PC_FieldSplitLink ilink; 240 PetscInt i,nsplit,ccsize; 241 MatStructure flag = pc->flag; 242 PetscTruth sorted; 243 244 PetscFunctionBegin; 245 ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 246 nsplit = jac->nsplits; 247 ilink = jac->head; 248 249 /* get the matrices for each split */ 250 if (!jac->issetup) { 251 PetscInt rstart,rend,nslots,bs; 252 253 jac->issetup = PETSC_TRUE; 254 255 /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 256 bs = jac->bs; 257 ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 258 ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 259 nslots = (rend - rstart)/bs; 260 for (i=0; i<nsplit; i++) { 261 if (jac->defaultsplit) { 262 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 263 } else if (!ilink->is) { 264 if (ilink->nfields > 1) { 265 PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 266 ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 267 for (j=0; j<nslots; j++) { 268 for (k=0; k<nfields; k++) { 269 ii[nfields*j + k] = rstart + bs*j + fields[k]; 270 } 271 } 272 ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 273 ierr = PetscFree(ii);CHKERRQ(ierr); 274 } else { 275 ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 276 } 277 } 278 ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 279 if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 280 ilink = ilink->next; 281 } 282 } 283 284 ilink = jac->head; 285 if (!jac->pmat) { 286 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 287 for (i=0; i<nsplit; i++) { 288 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 289 ilink = ilink->next; 290 } 291 } else { 292 for (i=0; i<nsplit; i++) { 293 ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 294 ilink = ilink->next; 295 } 296 } 297 if (jac->realdiagonal) { 298 ilink = jac->head; 299 if (!jac->mat) { 300 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 301 for (i=0; i<nsplit; i++) { 302 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 303 ilink = ilink->next; 304 } 305 } else { 306 for (i=0; i<nsplit; i++) { 307 ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 308 ilink = ilink->next; 309 } 310 } 311 } else { 312 jac->mat = jac->pmat; 313 } 314 315 if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 316 /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 317 ilink = jac->head; 318 if (!jac->Afield) { 319 ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 320 for (i=0; i<nsplit; i++) { 321 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 322 ilink = ilink->next; 323 } 324 } else { 325 for (i=0; i<nsplit; i++) { 326 ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 327 ilink = ilink->next; 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(jac->kspschur,&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 cnt; 490 491 PetscFunctionBegin; 492 CHKMEMQ; 493 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 494 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 495 496 if (jac->type == PC_COMPOSITE_ADDITIVE) { 497 if (jac->defaultsplit) { 498 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 499 while (ilink) { 500 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 501 ilink = ilink->next; 502 } 503 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 504 } else { 505 ierr = VecSet(y,0.0);CHKERRQ(ierr); 506 while (ilink) { 507 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 508 ilink = ilink->next; 509 } 510 } 511 } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 512 if (!jac->w1) { 513 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 514 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 515 } 516 ierr = VecSet(y,0.0);CHKERRQ(ierr); 517 ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 518 cnt = 1; 519 while (ilink->next) { 520 ilink = ilink->next; 521 /* compute the residual only over the part of the vector needed */ 522 ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 523 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 524 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 525 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 526 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 527 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 528 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 529 } 530 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 531 cnt -= 2; 532 while (ilink->previous) { 533 ilink = ilink->previous; 534 /* compute the residual only over the part of the vector needed */ 535 ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 536 ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 537 ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 538 ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 539 ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 540 ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 541 ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 542 } 543 } 544 } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 545 CHKMEMQ; 546 PetscFunctionReturn(0); 547 } 548 549 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 550 (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 551 VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 552 KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 553 VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 554 VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 555 556 #undef __FUNCT__ 557 #define __FUNCT__ "PCApply_FieldSplit" 558 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 559 { 560 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 561 PetscErrorCode ierr; 562 PC_FieldSplitLink ilink = jac->head; 563 564 PetscFunctionBegin; 565 CHKMEMQ; 566 ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 567 ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 568 569 if (jac->type == PC_COMPOSITE_ADDITIVE) { 570 if (jac->defaultsplit) { 571 ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 572 while (ilink) { 573 ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 574 ilink = ilink->next; 575 } 576 ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 577 } else { 578 ierr = VecSet(y,0.0);CHKERRQ(ierr); 579 while (ilink) { 580 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 581 ilink = ilink->next; 582 } 583 } 584 } else { 585 if (!jac->w1) { 586 ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 587 ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 588 } 589 ierr = VecSet(y,0.0);CHKERRQ(ierr); 590 if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 591 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 592 while (ilink->next) { 593 ilink = ilink->next; 594 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 595 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 596 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 597 } 598 while (ilink->previous) { 599 ilink = ilink->previous; 600 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 601 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 602 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 603 } 604 } else { 605 while (ilink->next) { /* get to last entry in linked list */ 606 ilink = ilink->next; 607 } 608 ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 609 while (ilink->previous) { 610 ilink = ilink->previous; 611 ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 612 ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 613 ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 614 } 615 } 616 } 617 CHKMEMQ; 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "PCDestroy_FieldSplit" 623 static PetscErrorCode PCDestroy_FieldSplit(PC pc) 624 { 625 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 626 PetscErrorCode ierr; 627 PC_FieldSplitLink ilink = jac->head,next; 628 629 PetscFunctionBegin; 630 while (ilink) { 631 ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 632 if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 633 if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 634 if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 635 if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 636 next = ilink->next; 637 ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 638 ierr = PetscFree(ilink);CHKERRQ(ierr); 639 ilink = next; 640 } 641 ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 642 if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 643 if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 644 if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 645 if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 646 if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 647 if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 648 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 649 if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 650 if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 651 if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 652 ierr = PetscFree(jac);CHKERRQ(ierr); 653 PetscFunctionReturn(0); 654 } 655 656 #undef __FUNCT__ 657 #define __FUNCT__ "PCSetFromOptions_FieldSplit" 658 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 659 { 660 PetscErrorCode ierr; 661 PetscInt i = 0,nfields,*fields,bs; 662 PetscTruth flg; 663 char optionname[128]; 664 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 665 PCCompositeType ctype; 666 667 PetscFunctionBegin; 668 ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 669 ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 670 ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 671 if (flg) { 672 ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 673 } 674 675 ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 676 if (flg) { 677 ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 678 } 679 680 /* Only setup fields once */ 681 if ((jac->bs > 0) && (jac->nsplits == 0)) { 682 /* only allow user to set fields from command line if bs is already known. 683 otherwise user can set them in PCFieldSplitSetDefaults() */ 684 ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 685 while (PETSC_TRUE) { 686 sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 687 nfields = jac->bs; 688 ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 689 if (!flg) break; 690 if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 691 ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 692 } 693 ierr = PetscFree(fields);CHKERRQ(ierr); 694 } 695 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); 696 ierr = PetscOptionsTail();CHKERRQ(ierr); 697 PetscFunctionReturn(0); 698 } 699 700 /*------------------------------------------------------------------------------------*/ 701 702 EXTERN_C_BEGIN 703 #undef __FUNCT__ 704 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 705 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 706 { 707 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 708 PetscErrorCode ierr; 709 PC_FieldSplitLink ilink,next = jac->head; 710 char prefix[128]; 711 PetscInt i; 712 713 PetscFunctionBegin; 714 if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 715 for (i=0; i<n; i++) { 716 if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 717 if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 718 } 719 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 720 ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 721 ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 722 ilink->nfields = n; 723 ilink->next = PETSC_NULL; 724 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 725 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 726 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 727 728 if (((PetscObject)pc)->prefix) { 729 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 730 } else { 731 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 732 } 733 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 734 735 if (!next) { 736 jac->head = ilink; 737 ilink->previous = PETSC_NULL; 738 } else { 739 while (next->next) { 740 next = next->next; 741 } 742 next->next = ilink; 743 ilink->previous = next; 744 } 745 jac->nsplits++; 746 PetscFunctionReturn(0); 747 } 748 EXTERN_C_END 749 750 EXTERN_C_BEGIN 751 #undef __FUNCT__ 752 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 753 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 754 { 755 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 756 PetscErrorCode ierr; 757 758 PetscFunctionBegin; 759 ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 760 ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 761 (*subksp)[1] = jac->kspschur; 762 *n = jac->nsplits; 763 PetscFunctionReturn(0); 764 } 765 EXTERN_C_END 766 767 EXTERN_C_BEGIN 768 #undef __FUNCT__ 769 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 770 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 771 { 772 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 773 PetscErrorCode ierr; 774 PetscInt cnt = 0; 775 PC_FieldSplitLink ilink = jac->head; 776 777 PetscFunctionBegin; 778 ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 779 while (ilink) { 780 (*subksp)[cnt++] = ilink->ksp; 781 ilink = ilink->next; 782 } 783 if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 784 *n = jac->nsplits; 785 PetscFunctionReturn(0); 786 } 787 EXTERN_C_END 788 789 EXTERN_C_BEGIN 790 #undef __FUNCT__ 791 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 792 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 793 { 794 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 795 PetscErrorCode ierr; 796 PC_FieldSplitLink ilink, next = jac->head; 797 char prefix[128]; 798 799 PetscFunctionBegin; 800 ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 801 ilink->is = is; 802 ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 803 ilink->next = PETSC_NULL; 804 ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 805 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 806 ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 807 808 if (((PetscObject)pc)->prefix) { 809 sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 810 } else { 811 sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 812 } 813 ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 814 815 if (!next) { 816 jac->head = ilink; 817 ilink->previous = PETSC_NULL; 818 } else { 819 while (next->next) { 820 next = next->next; 821 } 822 next->next = ilink; 823 ilink->previous = next; 824 } 825 jac->nsplits++; 826 827 PetscFunctionReturn(0); 828 } 829 EXTERN_C_END 830 831 #undef __FUNCT__ 832 #define __FUNCT__ "PCFieldSplitSetFields" 833 /*@ 834 PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 835 836 Collective on PC 837 838 Input Parameters: 839 + pc - the preconditioner context 840 . n - the number of fields in this split 841 . fields - the fields in this split 842 843 Level: intermediate 844 845 Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 846 847 The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 848 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 849 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 850 where the numbered entries indicate what is in the field. 851 852 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 853 854 @*/ 855 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 856 { 857 PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 858 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 861 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 862 if (f) { 863 ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 864 } 865 PetscFunctionReturn(0); 866 } 867 868 #undef __FUNCT__ 869 #define __FUNCT__ "PCFieldSplitSetIS" 870 /*@ 871 PCFieldSplitSetIS - Sets the exact elements for field 872 873 Collective on PC 874 875 Input Parameters: 876 + pc - the preconditioner context 877 . is - the index set that defines the vector elements in this field 878 879 880 Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 881 882 Level: intermediate 883 884 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 885 886 @*/ 887 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 888 { 889 PetscErrorCode ierr,(*f)(PC,IS); 890 891 PetscFunctionBegin; 892 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 893 PetscValidHeaderSpecific(is,IS_COOKIE,1); 894 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 895 if (f) { 896 ierr = (*f)(pc,is);CHKERRQ(ierr); 897 } 898 PetscFunctionReturn(0); 899 } 900 901 #undef __FUNCT__ 902 #define __FUNCT__ "PCFieldSplitSetBlockSize" 903 /*@ 904 PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 905 fieldsplit preconditioner. If not set the matrix block size is used. 906 907 Collective on PC 908 909 Input Parameters: 910 + pc - the preconditioner context 911 - bs - the block size 912 913 Level: intermediate 914 915 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 916 917 @*/ 918 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 919 { 920 PetscErrorCode ierr,(*f)(PC,PetscInt); 921 922 PetscFunctionBegin; 923 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 924 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 925 if (f) { 926 ierr = (*f)(pc,bs);CHKERRQ(ierr); 927 } 928 PetscFunctionReturn(0); 929 } 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "PCFieldSplitGetSubKSP" 933 /*@C 934 PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 935 936 Collective on KSP 937 938 Input Parameter: 939 . pc - the preconditioner context 940 941 Output Parameters: 942 + n - the number of split 943 - pc - the array of KSP contexts 944 945 Note: 946 After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 947 (not the KSP just the array that contains them). 948 949 You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 950 951 Level: advanced 952 953 .seealso: PCFIELDSPLIT 954 @*/ 955 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 956 { 957 PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 958 959 PetscFunctionBegin; 960 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 961 PetscValidIntPointer(n,2); 962 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 963 if (f) { 964 ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 965 } else { 966 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 967 } 968 PetscFunctionReturn(0); 969 } 970 971 #undef __FUNCT__ 972 #define __FUNCT__ "PCFieldSplitSchurPrecondition" 973 /*@ 974 PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 975 D matrix. Otherwise no preconditioner is used. 976 977 Collective on PC 978 979 Input Parameters: 980 + pc - the preconditioner context 981 . ptype - which matrix to use for preconditioning the Schur complement 982 - userpre - matrix to use for preconditioning, or PETSC_NULL 983 984 Notes: 985 The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 986 There are currently no preconditioners that work directly with the Schur complement so setting 987 PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 988 989 Options Database: 990 . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 991 992 Level: intermediate 993 994 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 995 996 @*/ 997 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 998 { 999 PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 1000 1001 PetscFunctionBegin; 1002 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1003 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 1004 if (f) { 1005 ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1006 } 1007 PetscFunctionReturn(0); 1008 } 1009 1010 EXTERN_C_BEGIN 1011 #undef __FUNCT__ 1012 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1013 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1014 { 1015 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1016 PetscErrorCode ierr; 1017 1018 PetscFunctionBegin; 1019 jac->schurpre = ptype; 1020 if (pre) { 1021 if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1022 jac->schur_user = pre; 1023 ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1024 } 1025 PetscFunctionReturn(0); 1026 } 1027 EXTERN_C_END 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 1031 /*@C 1032 PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 1033 1034 Collective on KSP 1035 1036 Input Parameter: 1037 . pc - the preconditioner context 1038 1039 Output Parameters: 1040 + A - the (0,0) block 1041 . B - the (0,1) block 1042 . C - the (1,0) block 1043 - D - the (1,1) block 1044 1045 Level: advanced 1046 1047 .seealso: PCFIELDSPLIT 1048 @*/ 1049 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 1050 { 1051 PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1052 1053 PetscFunctionBegin; 1054 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1055 if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");} 1056 if (A) *A = jac->pmat[0]; 1057 if (B) *B = jac->B; 1058 if (C) *C = jac->C; 1059 if (D) *D = jac->pmat[1]; 1060 PetscFunctionReturn(0); 1061 } 1062 1063 EXTERN_C_BEGIN 1064 #undef __FUNCT__ 1065 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1066 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 1067 { 1068 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1069 PetscErrorCode ierr; 1070 1071 PetscFunctionBegin; 1072 jac->type = type; 1073 if (type == PC_COMPOSITE_SCHUR) { 1074 pc->ops->apply = PCApply_FieldSplit_Schur; 1075 pc->ops->view = PCView_FieldSplit_Schur; 1076 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1077 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1078 1079 } else { 1080 pc->ops->apply = PCApply_FieldSplit; 1081 pc->ops->view = PCView_FieldSplit; 1082 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1083 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1084 } 1085 PetscFunctionReturn(0); 1086 } 1087 EXTERN_C_END 1088 1089 EXTERN_C_BEGIN 1090 #undef __FUNCT__ 1091 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 1092 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 1093 { 1094 PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1095 1096 PetscFunctionBegin; 1097 if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1098 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); 1099 jac->bs = bs; 1100 PetscFunctionReturn(0); 1101 } 1102 EXTERN_C_END 1103 1104 #undef __FUNCT__ 1105 #define __FUNCT__ "PCFieldSplitSetType" 1106 /*@ 1107 PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 1108 1109 Collective on PC 1110 1111 Input Parameter: 1112 . pc - the preconditioner context 1113 . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1114 1115 Options Database Key: 1116 . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 1117 1118 Level: Developer 1119 1120 .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1121 1122 .seealso: PCCompositeSetType() 1123 1124 @*/ 1125 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 1126 { 1127 PetscErrorCode ierr,(*f)(PC,PCCompositeType); 1128 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1131 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 1132 if (f) { 1133 ierr = (*f)(pc,type);CHKERRQ(ierr); 1134 } 1135 PetscFunctionReturn(0); 1136 } 1137 1138 /* -------------------------------------------------------------------------------------*/ 1139 /*MC 1140 PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1141 fields or groups of fields 1142 1143 1144 To set options on the solvers for each block append -fieldsplit_ to all the PC 1145 options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 1146 1147 To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 1148 and set the options directly on the resulting KSP object 1149 1150 Level: intermediate 1151 1152 Options Database Keys: 1153 + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 1154 . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 1155 been supplied explicitly by -pc_fieldsplit_%d_fields 1156 . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 1157 . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1158 . -pc_fieldsplit_schur_precondition <true,false> default is true 1159 1160 - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 1161 for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 1162 1163 1164 Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1165 to define a field by an arbitrary collection of entries. 1166 1167 If no fields are set the default is used. The fields are defined by entries strided by bs, 1168 beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1169 if this is not called the block size defaults to the blocksize of the second matrix passed 1170 to KSPSetOperators()/PCSetOperators(). 1171 1172 Currently for the multiplicative version, the updated residual needed for the next field 1173 solve is computed via a matrix vector product over the entire array. An optimization would be 1174 to update the residual only for the part of the right hand side associated with the next field 1175 solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1176 part of the matrix needed to just update part of the residual). 1177 1178 For the Schur complement preconditioner if J = ( A B ) 1179 ( C D ) 1180 the preconditioner is 1181 (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1182 (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1183 where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1184 inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1185 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1186 D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1187 option. 1188 1189 If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1190 is used automatically for a second block. 1191 1192 Concepts: physics based preconditioners, block preconditioners 1193 1194 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1195 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1196 M*/ 1197 1198 EXTERN_C_BEGIN 1199 #undef __FUNCT__ 1200 #define __FUNCT__ "PCCreate_FieldSplit" 1201 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 1202 { 1203 PetscErrorCode ierr; 1204 PC_FieldSplit *jac; 1205 1206 PetscFunctionBegin; 1207 ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 1208 jac->bs = -1; 1209 jac->nsplits = 0; 1210 jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1211 jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1212 1213 pc->data = (void*)jac; 1214 1215 pc->ops->apply = PCApply_FieldSplit; 1216 pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 1217 pc->ops->setup = PCSetUp_FieldSplit; 1218 pc->ops->destroy = PCDestroy_FieldSplit; 1219 pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 1220 pc->ops->view = PCView_FieldSplit; 1221 pc->ops->applyrichardson = 0; 1222 1223 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 1224 PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1225 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 1226 PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1227 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1228 PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 1229 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 1230 PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 1231 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 1232 PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 1233 PetscFunctionReturn(0); 1234 } 1235 EXTERN_C_END 1236 1237 1238 /*MC 1239 Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1240 overview of these methods. 1241 1242 Consider the solution to ( A B ) (x_1) = (b_1) 1243 ( C D ) (x_2) (b_2) 1244 1245 Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1246 B' 0) (x_2) (b_2) 1247 1248 One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1249 for this block system. 1250 1251 Consider an additional matrix (Ap Bp) 1252 (Cp Dp) where some or all of the entries may be the same as 1253 in the original matrix (for example Ap == A). 1254 1255 In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1256 approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1257 1258 Block Jacobi: x_1 = A^ b_1 1259 x_2 = D^ b_2 1260 1261 Lower block Gauss-Seidel: x_1 = A^ b_1 1262 x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1263 1264 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) 1265 Interestingly this form is not actually a symmetric matrix, the symmetric version is 1266 x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1267 1268 Level: intermediate 1269 1270 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1271 PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1272 M*/ 1273