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