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