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