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