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