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