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