1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 30971522cSBarry Smith /* 40971522cSBarry Smith 50971522cSBarry Smith */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 70971522cSBarry Smith 80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 90971522cSBarry Smith struct _PC_FieldSplitLink { 1069a612a9SBarry Smith KSP ksp; 110971522cSBarry Smith Vec x,y; 120971522cSBarry Smith PetscInt nfields; 130971522cSBarry Smith PetscInt *fields; 141b9fc7fcSBarry Smith VecScatter sctx; 15704ba839SBarry Smith IS is,cis; 16704ba839SBarry Smith PetscInt csize; 1751f519a2SBarry Smith PC_FieldSplitLink next,previous; 180971522cSBarry Smith }; 190971522cSBarry Smith 200971522cSBarry Smith typedef struct { 2168dd23aaSBarry Smith PCCompositeType type; 2297bbdb24SBarry Smith PetscTruth defaultsplit; 230971522cSBarry Smith PetscInt bs; 24704ba839SBarry Smith PetscInt nsplits; 2579416396SBarry Smith Vec *x,*y,w1,w2; 2697bbdb24SBarry Smith Mat *pmat; 2768dd23aaSBarry Smith Mat *Afield; /* the rows of the matrix associated with each field */ 28704ba839SBarry Smith PetscTruth issetup; 293b224e63SBarry Smith Mat B,C,schur; /* only used when Schur complement preconditioning is used */ 303b224e63SBarry Smith KSP kspschur; 31*e69d4d44SBarry Smith PetscTruth schurpre; /* preconditioner for the Schur complement is built from pmat[1] == D */ 3297bbdb24SBarry Smith PC_FieldSplitLink head; 330971522cSBarry Smith } PC_FieldSplit; 340971522cSBarry Smith 3516913363SBarry Smith /* 3616913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 3716913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 3816913363SBarry Smith PC you could change this. 3916913363SBarry Smith */ 400971522cSBarry Smith #undef __FUNCT__ 410971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 420971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 430971522cSBarry Smith { 440971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 450971522cSBarry Smith PetscErrorCode ierr; 460971522cSBarry Smith PetscTruth iascii; 470971522cSBarry Smith PetscInt i,j; 485a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 490971522cSBarry Smith 500971522cSBarry Smith PetscFunctionBegin; 510971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 520971522cSBarry Smith if (iascii) { 5351f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 5469a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 550971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 560971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 571ab39975SBarry Smith if (ilink->fields) { 580971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 5979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 605a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 6179416396SBarry Smith if (j > 0) { 6279416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 6379416396SBarry Smith } 645a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 650971522cSBarry Smith } 660971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6779416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 681ab39975SBarry Smith } else { 691ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 701ab39975SBarry Smith } 715a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 725a9f2f41SSatish Balay ilink = ilink->next; 730971522cSBarry Smith } 740971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 750971522cSBarry Smith } else { 760971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 770971522cSBarry Smith } 780971522cSBarry Smith PetscFunctionReturn(0); 790971522cSBarry Smith } 800971522cSBarry Smith 810971522cSBarry Smith #undef __FUNCT__ 823b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 833b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 843b224e63SBarry Smith { 853b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 863b224e63SBarry Smith PetscErrorCode ierr; 873b224e63SBarry Smith PetscTruth iascii; 883b224e63SBarry Smith PetscInt i,j; 893b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 903b224e63SBarry Smith KSP ksp; 913b224e63SBarry Smith 923b224e63SBarry Smith PetscFunctionBegin; 933b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 943b224e63SBarry Smith if (iascii) { 953b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 963b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 973b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 983b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 993b224e63SBarry Smith if (ilink->fields) { 1003b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1013b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1023b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1033b224e63SBarry Smith if (j > 0) { 1043b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1053b224e63SBarry Smith } 1063b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1073b224e63SBarry Smith } 1083b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1093b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1103b224e63SBarry Smith } else { 1113b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1123b224e63SBarry Smith } 1133b224e63SBarry Smith ilink = ilink->next; 1143b224e63SBarry Smith } 1153b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1163b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1173b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1183b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 1193b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1203b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1213b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1223b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 1233b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1243b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1253b224e63SBarry Smith } else { 1263b224e63SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1273b224e63SBarry Smith } 1283b224e63SBarry Smith PetscFunctionReturn(0); 1293b224e63SBarry Smith } 1303b224e63SBarry Smith 1313b224e63SBarry Smith #undef __FUNCT__ 13269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 13369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1340971522cSBarry Smith { 1350971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1360971522cSBarry Smith PetscErrorCode ierr; 1375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 138d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 139d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 140d32f9abdSBarry Smith char optionname[128]; 1410971522cSBarry Smith 1420971522cSBarry Smith PetscFunctionBegin; 143d32f9abdSBarry Smith if (!ilink) { 144d32f9abdSBarry Smith 145521d7252SBarry Smith if (jac->bs <= 0) { 146704ba839SBarry Smith if (pc->pmat) { 147521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 148704ba839SBarry Smith } else { 149704ba839SBarry Smith jac->bs = 1; 150704ba839SBarry Smith } 151521d7252SBarry Smith } 152d32f9abdSBarry Smith 153d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 154d32f9abdSBarry Smith if (!flg) { 155d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 156d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 157d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 158d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 159d32f9abdSBarry Smith while (PETSC_TRUE) { 160d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 161d32f9abdSBarry Smith nfields = jac->bs; 162d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 163d32f9abdSBarry Smith if (!flg2) break; 164d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 165d32f9abdSBarry Smith flg = PETSC_FALSE; 166d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 167d32f9abdSBarry Smith } 168d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 169d32f9abdSBarry Smith } 170d32f9abdSBarry Smith 171d32f9abdSBarry Smith if (flg) { 172d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 17379416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 17479416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1755a9f2f41SSatish Balay while (ilink) { 1765a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1775a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 178521d7252SBarry Smith } 1795a9f2f41SSatish Balay ilink = ilink->next; 18079416396SBarry Smith } 18197bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 18279416396SBarry Smith for (i=0; i<jac->bs; i++) { 18379416396SBarry Smith if (!fields[i]) { 18479416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 18579416396SBarry Smith } else { 18679416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 18779416396SBarry Smith } 18879416396SBarry Smith } 18968dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 190521d7252SBarry Smith } 191d32f9abdSBarry Smith } 19269a612a9SBarry Smith PetscFunctionReturn(0); 19369a612a9SBarry Smith } 19469a612a9SBarry Smith 19569a612a9SBarry Smith 19669a612a9SBarry Smith #undef __FUNCT__ 19769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 19869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 19969a612a9SBarry Smith { 20069a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 20169a612a9SBarry Smith PetscErrorCode ierr; 2025a9f2f41SSatish Balay PC_FieldSplitLink ilink; 203cf502942SBarry Smith PetscInt i,nsplit,ccsize; 20469a612a9SBarry Smith MatStructure flag = pc->flag; 20568dd23aaSBarry Smith PetscTruth sorted,getsub; 20669a612a9SBarry Smith 20769a612a9SBarry Smith PetscFunctionBegin; 20869a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 20997bbdb24SBarry Smith nsplit = jac->nsplits; 2105a9f2f41SSatish Balay ilink = jac->head; 21197bbdb24SBarry Smith 21297bbdb24SBarry Smith /* get the matrices for each split */ 213704ba839SBarry Smith if (!jac->issetup) { 2141b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 21597bbdb24SBarry Smith 216704ba839SBarry Smith jac->issetup = PETSC_TRUE; 217704ba839SBarry Smith 218704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 21951f519a2SBarry Smith bs = jac->bs; 22097bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 221cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2221b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2231b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2241b9fc7fcSBarry Smith if (jac->defaultsplit) { 225704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 226704ba839SBarry Smith } else if (!ilink->is) { 227ccb205f8SBarry Smith if (ilink->nfields > 1) { 2285a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2295a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2301b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2311b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2321b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 23397bbdb24SBarry Smith } 23497bbdb24SBarry Smith } 235704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 236ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 237ccb205f8SBarry Smith } else { 238704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 239ccb205f8SBarry Smith } 2403e197d65SBarry Smith } 2413e197d65SBarry Smith ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr); 242704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 2431b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 244704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 245704ba839SBarry Smith ilink = ilink->next; 2461b9fc7fcSBarry Smith } 2471b9fc7fcSBarry Smith } 2481b9fc7fcSBarry Smith 249704ba839SBarry Smith ilink = jac->head; 25097bbdb24SBarry Smith if (!jac->pmat) { 251cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 252cf502942SBarry Smith for (i=0; i<nsplit; i++) { 253704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 254704ba839SBarry Smith ilink = ilink->next; 255cf502942SBarry Smith } 25697bbdb24SBarry Smith } else { 257cf502942SBarry Smith for (i=0; i<nsplit; i++) { 258704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 259704ba839SBarry Smith ilink = ilink->next; 260cf502942SBarry Smith } 26197bbdb24SBarry Smith } 26297bbdb24SBarry Smith 26368dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 26468dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 2653b224e63SBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 26668dd23aaSBarry Smith ilink = jac->head; 26768dd23aaSBarry Smith if (!jac->Afield) { 26868dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 26968dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 27011755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 27168dd23aaSBarry Smith ilink = ilink->next; 27268dd23aaSBarry Smith } 27368dd23aaSBarry Smith } else { 27468dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 27511755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 27668dd23aaSBarry Smith ilink = ilink->next; 27768dd23aaSBarry Smith } 27868dd23aaSBarry Smith } 27968dd23aaSBarry Smith } 28068dd23aaSBarry Smith 2813b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 2823b224e63SBarry Smith IS ccis; 283*e69d4d44SBarry Smith PetscInt N,nlocal,nis; 2843b224e63SBarry Smith if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 28568dd23aaSBarry Smith 2863b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 2873b224e63SBarry Smith if (jac->schur) { 2883b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 2893b224e63SBarry Smith ilink = jac->head; 2903b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 291*e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 292*e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 293*e69d4d44SBarry Smith nlocal = nlocal - nis; 294*e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 2953b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 2963b224e63SBarry Smith ilink = ilink->next; 2973b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 298*e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 299*e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 300*e69d4d44SBarry Smith nlocal = nlocal - nis; 301*e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3023b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3033b224e63SBarry Smith ierr = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 304*e69d4d44SBarry Smith if (jac->schurpre) { 305*e69d4d44SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],pc->flag);CHKERRQ(ierr); 306*e69d4d44SBarry Smith } else { 3073b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,pc->flag);CHKERRQ(ierr); 308*e69d4d44SBarry Smith } 3093b224e63SBarry Smith 3103b224e63SBarry Smith } else { 3111cee3971SBarry Smith KSP ksp; 3123b224e63SBarry Smith 3133b224e63SBarry Smith /* extract the B and C matrices */ 3143b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3153b224e63SBarry Smith ilink = jac->head; 3163b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 317*e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 318*e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 319*e69d4d44SBarry Smith nlocal = nlocal - nis; 320*e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3213b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3223b224e63SBarry Smith ilink = ilink->next; 3233b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 324*e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 325*e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 326*e69d4d44SBarry Smith nlocal = nlocal - nis; 327*e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3283b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3293b224e63SBarry Smith ierr = MatCreateSchurComplement(jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr); 3301cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 331*e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3323b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3333b224e63SBarry Smith 3343b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3351cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 336*e69d4d44SBarry Smith if (jac->schurpre) { 337*e69d4d44SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 338*e69d4d44SBarry Smith } else { 3393b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 340*e69d4d44SBarry Smith } 3413b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 3423b224e63SBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"schur_");CHKERRQ(ierr); 3433b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3443b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3453b224e63SBarry Smith 3463b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3473b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3483b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3493b224e63SBarry Smith ilink = jac->head; 3503b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3513b224e63SBarry Smith ilink = ilink->next; 3523b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3533b224e63SBarry Smith } 3543b224e63SBarry Smith } else { 35597bbdb24SBarry Smith /* set up the individual PCs */ 35697bbdb24SBarry Smith i = 0; 3575a9f2f41SSatish Balay ilink = jac->head; 3585a9f2f41SSatish Balay while (ilink) { 3595a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 3603b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3615a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 3625a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 36397bbdb24SBarry Smith i++; 3645a9f2f41SSatish Balay ilink = ilink->next; 3650971522cSBarry Smith } 36697bbdb24SBarry Smith 36797bbdb24SBarry Smith /* create work vectors for each split */ 3681b9fc7fcSBarry Smith if (!jac->x) { 36979416396SBarry Smith Vec xtmp; 37097bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 3715a9f2f41SSatish Balay ilink = jac->head; 37297bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 373906ed7ccSBarry Smith Vec *vl,*vr; 374906ed7ccSBarry Smith 375906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 376906ed7ccSBarry Smith ilink->x = *vr; 377906ed7ccSBarry Smith ilink->y = *vl; 378906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 379906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 3805a9f2f41SSatish Balay jac->x[i] = ilink->x; 3815a9f2f41SSatish Balay jac->y[i] = ilink->y; 3825a9f2f41SSatish Balay ilink = ilink->next; 38397bbdb24SBarry Smith } 3843b224e63SBarry Smith } 3853b224e63SBarry Smith } 3863b224e63SBarry Smith 3873b224e63SBarry Smith 3883b224e63SBarry Smith if (1) { 3893b224e63SBarry Smith Vec xtmp; 3903b224e63SBarry Smith 39179416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 3921b9fc7fcSBarry Smith 3935a9f2f41SSatish Balay ilink = jac->head; 3941b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 3951b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 396704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 3975a9f2f41SSatish Balay ilink = ilink->next; 3981b9fc7fcSBarry Smith } 3991b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4001b9fc7fcSBarry Smith } 4010971522cSBarry Smith PetscFunctionReturn(0); 4020971522cSBarry Smith } 4030971522cSBarry Smith 4045a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 405ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 406ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4075a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 408ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 409ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 41079416396SBarry Smith 4110971522cSBarry Smith #undef __FUNCT__ 4123b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4133b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4143b224e63SBarry Smith { 4153b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4163b224e63SBarry Smith PetscErrorCode ierr; 4173b224e63SBarry Smith KSP ksp; 4183b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4193b224e63SBarry Smith 4203b224e63SBarry Smith PetscFunctionBegin; 4213b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4223b224e63SBarry Smith 4233b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4243b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4253b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4263b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4273b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4283b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4293b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4303b224e63SBarry Smith 4313b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4323b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4333b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4343b224e63SBarry Smith 4353b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4363b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4373b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4383b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4393b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4403b224e63SBarry Smith 4413b224e63SBarry Smith PetscFunctionReturn(0); 4423b224e63SBarry Smith } 4433b224e63SBarry Smith 4443b224e63SBarry Smith #undef __FUNCT__ 4450971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4460971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4470971522cSBarry Smith { 4480971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4490971522cSBarry Smith PetscErrorCode ierr; 4505a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4513e197d65SBarry Smith PetscInt bs,cnt; 4520971522cSBarry Smith 4530971522cSBarry Smith PetscFunctionBegin; 45451f519a2SBarry Smith CHKMEMQ; 455e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 45651f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 45751f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 45851f519a2SBarry Smith 45979416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4601b9fc7fcSBarry Smith if (jac->defaultsplit) { 4610971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 4625a9f2f41SSatish Balay while (ilink) { 4635a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4645a9f2f41SSatish Balay ilink = ilink->next; 4650971522cSBarry Smith } 4660971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 4671b9fc7fcSBarry Smith } else { 468efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4695a9f2f41SSatish Balay while (ilink) { 4705a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4715a9f2f41SSatish Balay ilink = ilink->next; 4721b9fc7fcSBarry Smith } 4731b9fc7fcSBarry Smith } 47416913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 47579416396SBarry Smith if (!jac->w1) { 47679416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 47779416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 47879416396SBarry Smith } 479efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4805a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4813e197d65SBarry Smith cnt = 1; 4825a9f2f41SSatish Balay while (ilink->next) { 4835a9f2f41SSatish Balay ilink = ilink->next; 4843e197d65SBarry Smith if (jac->Afield) { 4853e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 4863e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 4873e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 4883e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4893e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4903e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4913e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4923e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4933e197d65SBarry Smith } else { 4943e197d65SBarry Smith /* compute the residual over the entire vector */ 4951ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 496efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 4975a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 49879416396SBarry Smith } 4993e197d65SBarry Smith } 50051f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 50111755939SBarry Smith cnt -= 2; 50251f519a2SBarry Smith while (ilink->previous) { 50351f519a2SBarry Smith ilink = ilink->previous; 50411755939SBarry Smith if (jac->Afield) { 50511755939SBarry Smith /* compute the residual only over the part of the vector needed */ 50611755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 50711755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 50811755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 50911755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 51011755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 51111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 51211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 51311755939SBarry Smith } else { 5141ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 51551f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 51651f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 51779416396SBarry Smith } 51851f519a2SBarry Smith } 51911755939SBarry Smith } 52016913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 52151f519a2SBarry Smith CHKMEMQ; 5220971522cSBarry Smith PetscFunctionReturn(0); 5230971522cSBarry Smith } 5240971522cSBarry Smith 525421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 526ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 527ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 528421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 529ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 530ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 531421e10b8SBarry Smith 532421e10b8SBarry Smith #undef __FUNCT__ 533421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 534421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 535421e10b8SBarry Smith { 536421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 537421e10b8SBarry Smith PetscErrorCode ierr; 538421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 539421e10b8SBarry Smith PetscInt bs; 540421e10b8SBarry Smith 541421e10b8SBarry Smith PetscFunctionBegin; 542421e10b8SBarry Smith CHKMEMQ; 543421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 544421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 545421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 546421e10b8SBarry Smith 547421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 548421e10b8SBarry Smith if (jac->defaultsplit) { 549421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 550421e10b8SBarry Smith while (ilink) { 551421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 552421e10b8SBarry Smith ilink = ilink->next; 553421e10b8SBarry Smith } 554421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 555421e10b8SBarry Smith } else { 556421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 557421e10b8SBarry Smith while (ilink) { 558421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 559421e10b8SBarry Smith ilink = ilink->next; 560421e10b8SBarry Smith } 561421e10b8SBarry Smith } 562421e10b8SBarry Smith } else { 563421e10b8SBarry Smith if (!jac->w1) { 564421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 565421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 566421e10b8SBarry Smith } 567421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 568421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 569421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 570421e10b8SBarry Smith while (ilink->next) { 571421e10b8SBarry Smith ilink = ilink->next; 5729989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 573421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 574421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 575421e10b8SBarry Smith } 576421e10b8SBarry Smith while (ilink->previous) { 577421e10b8SBarry Smith ilink = ilink->previous; 5789989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 579421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 580421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 581421e10b8SBarry Smith } 582421e10b8SBarry Smith } else { 583421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 584421e10b8SBarry Smith ilink = ilink->next; 585421e10b8SBarry Smith } 586421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 587421e10b8SBarry Smith while (ilink->previous) { 588421e10b8SBarry Smith ilink = ilink->previous; 5899989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 590421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 591421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 592421e10b8SBarry Smith } 593421e10b8SBarry Smith } 594421e10b8SBarry Smith } 595421e10b8SBarry Smith CHKMEMQ; 596421e10b8SBarry Smith PetscFunctionReturn(0); 597421e10b8SBarry Smith } 598421e10b8SBarry Smith 5990971522cSBarry Smith #undef __FUNCT__ 6000971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6010971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6020971522cSBarry Smith { 6030971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6040971522cSBarry Smith PetscErrorCode ierr; 6055a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6060971522cSBarry Smith 6070971522cSBarry Smith PetscFunctionBegin; 6085a9f2f41SSatish Balay while (ilink) { 6095a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6105a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6115a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6125a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 613704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 614704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 6155a9f2f41SSatish Balay next = ilink->next; 616704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 617704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6185a9f2f41SSatish Balay ilink = next; 6190971522cSBarry Smith } 62005b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 62197bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 62268dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 62379416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 62479416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6253b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 6263b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6273b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6280971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6290971522cSBarry Smith PetscFunctionReturn(0); 6300971522cSBarry Smith } 6310971522cSBarry Smith 6320971522cSBarry Smith #undef __FUNCT__ 6330971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6340971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6350971522cSBarry Smith { 6361b9fc7fcSBarry Smith PetscErrorCode ierr; 63751f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 638*e69d4d44SBarry Smith PetscTruth flg,set; 6391b9fc7fcSBarry Smith char optionname[128]; 6409dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6413b224e63SBarry Smith PCCompositeType ctype; 6421b9fc7fcSBarry Smith 6430971522cSBarry Smith PetscFunctionBegin; 6441b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 64551f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 64651f519a2SBarry Smith if (flg) { 64751f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 64851f519a2SBarry Smith } 649704ba839SBarry Smith 6503b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6513b224e63SBarry Smith if (flg) { 6523b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6533b224e63SBarry Smith } 654d32f9abdSBarry Smith 655d32f9abdSBarry Smith if (jac->bs > 0) { 656d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 657d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 65851f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6591b9fc7fcSBarry Smith while (PETSC_TRUE) { 66013f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 66151f519a2SBarry Smith nfields = jac->bs; 6621b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6631b9fc7fcSBarry Smith if (!flg) break; 6641b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 6651b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 6661b9fc7fcSBarry Smith } 66751f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 668d32f9abdSBarry Smith } 669*e69d4d44SBarry Smith ierr = PetscOptionsTruth("-pc_fieldsplit_schur_precondition","Build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",jac->schurpre,&flg,&set);CHKERRQ(ierr); 670*e69d4d44SBarry Smith if (set) { 671*e69d4d44SBarry Smith ierr = PCFieldSplitSchurPrecondition(pc,flg);CHKERRQ(ierr); 672*e69d4d44SBarry Smith } 6731b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6740971522cSBarry Smith PetscFunctionReturn(0); 6750971522cSBarry Smith } 6760971522cSBarry Smith 6770971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 6780971522cSBarry Smith 6790971522cSBarry Smith EXTERN_C_BEGIN 6800971522cSBarry Smith #undef __FUNCT__ 6810971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 682dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 6830971522cSBarry Smith { 68497bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6850971522cSBarry Smith PetscErrorCode ierr; 6865a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 68769a612a9SBarry Smith char prefix[128]; 68851f519a2SBarry Smith PetscInt i; 6890971522cSBarry Smith 6900971522cSBarry Smith PetscFunctionBegin; 6910971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 69251f519a2SBarry Smith for (i=0; i<n; i++) { 69351f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 69451f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 69551f519a2SBarry Smith } 696704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 697704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 6985a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 6995a9f2f41SSatish Balay ilink->nfields = n; 7005a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7017adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7021cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7035a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 70469a612a9SBarry Smith 7057adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 7067adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 70769a612a9SBarry Smith } else { 70813f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 70969a612a9SBarry Smith } 7105a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7110971522cSBarry Smith 7120971522cSBarry Smith if (!next) { 7135a9f2f41SSatish Balay jac->head = ilink; 71451f519a2SBarry Smith ilink->previous = PETSC_NULL; 7150971522cSBarry Smith } else { 7160971522cSBarry Smith while (next->next) { 7170971522cSBarry Smith next = next->next; 7180971522cSBarry Smith } 7195a9f2f41SSatish Balay next->next = ilink; 72051f519a2SBarry Smith ilink->previous = next; 7210971522cSBarry Smith } 7220971522cSBarry Smith jac->nsplits++; 7230971522cSBarry Smith PetscFunctionReturn(0); 7240971522cSBarry Smith } 7250971522cSBarry Smith EXTERN_C_END 7260971522cSBarry Smith 727*e69d4d44SBarry Smith EXTERN_C_BEGIN 728*e69d4d44SBarry Smith #undef __FUNCT__ 729*e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 730*e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 731*e69d4d44SBarry Smith { 732*e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 733*e69d4d44SBarry Smith PetscErrorCode ierr; 734*e69d4d44SBarry Smith PetscInt cnt = 0; 735*e69d4d44SBarry Smith PC_FieldSplitLink ilink = jac->head; 736*e69d4d44SBarry Smith 737*e69d4d44SBarry Smith PetscFunctionBegin; 738*e69d4d44SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 739*e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 740*e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 741*e69d4d44SBarry Smith PetscFunctionReturn(0); 742*e69d4d44SBarry Smith } 743*e69d4d44SBarry Smith EXTERN_C_END 7440971522cSBarry Smith 7450971522cSBarry Smith EXTERN_C_BEGIN 7460971522cSBarry Smith #undef __FUNCT__ 74769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 748dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7490971522cSBarry Smith { 7500971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7510971522cSBarry Smith PetscErrorCode ierr; 7520971522cSBarry Smith PetscInt cnt = 0; 7535a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7540971522cSBarry Smith 7550971522cSBarry Smith PetscFunctionBegin; 75669a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7575a9f2f41SSatish Balay while (ilink) { 7585a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7595a9f2f41SSatish Balay ilink = ilink->next; 7600971522cSBarry Smith } 7610971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7620971522cSBarry Smith *n = jac->nsplits; 7630971522cSBarry Smith PetscFunctionReturn(0); 7640971522cSBarry Smith } 7650971522cSBarry Smith EXTERN_C_END 7660971522cSBarry Smith 767704ba839SBarry Smith EXTERN_C_BEGIN 768704ba839SBarry Smith #undef __FUNCT__ 769704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 770704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 771704ba839SBarry Smith { 772704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 773704ba839SBarry Smith PetscErrorCode ierr; 774704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 775704ba839SBarry Smith char prefix[128]; 776704ba839SBarry Smith 777704ba839SBarry Smith PetscFunctionBegin; 77816913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 7791ab39975SBarry Smith ilink->is = is; 7801ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 781704ba839SBarry Smith ilink->next = PETSC_NULL; 782704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7831cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 784704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 785704ba839SBarry Smith 786704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 787704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 788704ba839SBarry Smith } else { 789704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 790704ba839SBarry Smith } 791704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 792704ba839SBarry Smith 793704ba839SBarry Smith if (!next) { 794704ba839SBarry Smith jac->head = ilink; 795704ba839SBarry Smith ilink->previous = PETSC_NULL; 796704ba839SBarry Smith } else { 797704ba839SBarry Smith while (next->next) { 798704ba839SBarry Smith next = next->next; 799704ba839SBarry Smith } 800704ba839SBarry Smith next->next = ilink; 801704ba839SBarry Smith ilink->previous = next; 802704ba839SBarry Smith } 803704ba839SBarry Smith jac->nsplits++; 804704ba839SBarry Smith 805704ba839SBarry Smith PetscFunctionReturn(0); 806704ba839SBarry Smith } 807704ba839SBarry Smith EXTERN_C_END 808704ba839SBarry Smith 8090971522cSBarry Smith #undef __FUNCT__ 8100971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8110971522cSBarry Smith /*@ 8120971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8130971522cSBarry Smith 8140971522cSBarry Smith Collective on PC 8150971522cSBarry Smith 8160971522cSBarry Smith Input Parameters: 8170971522cSBarry Smith + pc - the preconditioner context 8180971522cSBarry Smith . n - the number of fields in this split 8190971522cSBarry Smith . fields - the fields in this split 8200971522cSBarry Smith 8210971522cSBarry Smith Level: intermediate 8220971522cSBarry Smith 823d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 824d32f9abdSBarry Smith 825d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 826d32f9abdSBarry Smith 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 827d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 828d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 829d32f9abdSBarry Smith 830d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8310971522cSBarry Smith 8320971522cSBarry Smith @*/ 833dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 8340971522cSBarry Smith { 8350971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 8360971522cSBarry Smith 8370971522cSBarry Smith PetscFunctionBegin; 8380971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8390971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8400971522cSBarry Smith if (f) { 8410971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8420971522cSBarry Smith } 8430971522cSBarry Smith PetscFunctionReturn(0); 8440971522cSBarry Smith } 8450971522cSBarry Smith 8460971522cSBarry Smith #undef __FUNCT__ 847704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 848704ba839SBarry Smith /*@ 849704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 850704ba839SBarry Smith 851704ba839SBarry Smith Collective on PC 852704ba839SBarry Smith 853704ba839SBarry Smith Input Parameters: 854704ba839SBarry Smith + pc - the preconditioner context 855704ba839SBarry Smith . is - the index set that defines the vector elements in this field 856704ba839SBarry Smith 857d32f9abdSBarry Smith 858d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 859d32f9abdSBarry Smith 860704ba839SBarry Smith Level: intermediate 861704ba839SBarry Smith 862704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 863704ba839SBarry Smith 864704ba839SBarry Smith @*/ 865704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 866704ba839SBarry Smith { 867704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 868704ba839SBarry Smith 869704ba839SBarry Smith PetscFunctionBegin; 870704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 871704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 872704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 873704ba839SBarry Smith if (f) { 874704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 875704ba839SBarry Smith } 876704ba839SBarry Smith PetscFunctionReturn(0); 877704ba839SBarry Smith } 878704ba839SBarry Smith 879704ba839SBarry Smith #undef __FUNCT__ 88051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 88151f519a2SBarry Smith /*@ 88251f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 88351f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 88451f519a2SBarry Smith 88551f519a2SBarry Smith Collective on PC 88651f519a2SBarry Smith 88751f519a2SBarry Smith Input Parameters: 88851f519a2SBarry Smith + pc - the preconditioner context 88951f519a2SBarry Smith - bs - the block size 89051f519a2SBarry Smith 89151f519a2SBarry Smith Level: intermediate 89251f519a2SBarry Smith 89351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 89451f519a2SBarry Smith 89551f519a2SBarry Smith @*/ 89651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 89751f519a2SBarry Smith { 89851f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 89951f519a2SBarry Smith 90051f519a2SBarry Smith PetscFunctionBegin; 90151f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 90251f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 90351f519a2SBarry Smith if (f) { 90451f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 90551f519a2SBarry Smith } 90651f519a2SBarry Smith PetscFunctionReturn(0); 90751f519a2SBarry Smith } 90851f519a2SBarry Smith 90951f519a2SBarry Smith #undef __FUNCT__ 91069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9110971522cSBarry Smith /*@C 91269a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9130971522cSBarry Smith 91469a612a9SBarry Smith Collective on KSP 9150971522cSBarry Smith 9160971522cSBarry Smith Input Parameter: 9170971522cSBarry Smith . pc - the preconditioner context 9180971522cSBarry Smith 9190971522cSBarry Smith Output Parameters: 9200971522cSBarry Smith + n - the number of split 92169a612a9SBarry Smith - pc - the array of KSP contexts 9220971522cSBarry Smith 9230971522cSBarry Smith Note: 924d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 925d32f9abdSBarry Smith (not the KSP just the array that contains them). 9260971522cSBarry Smith 92769a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9280971522cSBarry Smith 9290971522cSBarry Smith Level: advanced 9300971522cSBarry Smith 9310971522cSBarry Smith .seealso: PCFIELDSPLIT 9320971522cSBarry Smith @*/ 933dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9340971522cSBarry Smith { 93569a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9360971522cSBarry Smith 9370971522cSBarry Smith PetscFunctionBegin; 9380971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9390971522cSBarry Smith PetscValidIntPointer(n,2); 94069a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9410971522cSBarry Smith if (f) { 94269a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 9430971522cSBarry Smith } else { 94469a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9450971522cSBarry Smith } 9460971522cSBarry Smith PetscFunctionReturn(0); 9470971522cSBarry Smith } 9480971522cSBarry Smith 949*e69d4d44SBarry Smith #undef __FUNCT__ 950*e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 951*e69d4d44SBarry Smith /*@ 952*e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 953*e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 954*e69d4d44SBarry Smith 955*e69d4d44SBarry Smith Collective on PC 956*e69d4d44SBarry Smith 957*e69d4d44SBarry Smith Input Parameters: 958*e69d4d44SBarry Smith + pc - the preconditioner context 959*e69d4d44SBarry Smith - flg - build the preconditioner 960*e69d4d44SBarry Smith 961*e69d4d44SBarry Smith Options Database: 962*e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 963*e69d4d44SBarry Smith 964*e69d4d44SBarry Smith Level: intermediate 965*e69d4d44SBarry Smith 966*e69d4d44SBarry Smith Notes: What should we do if the user wants to provide a different matrix (like a mass matrix) to use to build the preconditioner 967*e69d4d44SBarry Smith 968*e69d4d44SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFIELDSPLIT 969*e69d4d44SBarry Smith 970*e69d4d44SBarry Smith @*/ 971*e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PetscTruth flg) 972*e69d4d44SBarry Smith { 973*e69d4d44SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 974*e69d4d44SBarry Smith 975*e69d4d44SBarry Smith PetscFunctionBegin; 976*e69d4d44SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 977*e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 978*e69d4d44SBarry Smith if (f) { 979*e69d4d44SBarry Smith ierr = (*f)(pc,flg);CHKERRQ(ierr); 980*e69d4d44SBarry Smith } 981*e69d4d44SBarry Smith PetscFunctionReturn(0); 982*e69d4d44SBarry Smith } 983*e69d4d44SBarry Smith 984*e69d4d44SBarry Smith EXTERN_C_BEGIN 985*e69d4d44SBarry Smith #undef __FUNCT__ 986*e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 987*e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PetscTruth flg) 988*e69d4d44SBarry Smith { 989*e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 990*e69d4d44SBarry Smith 991*e69d4d44SBarry Smith PetscFunctionBegin; 992*e69d4d44SBarry Smith jac->schurpre = flg; 993*e69d4d44SBarry Smith PetscFunctionReturn(0); 994*e69d4d44SBarry Smith } 995*e69d4d44SBarry Smith EXTERN_C_END 996*e69d4d44SBarry Smith 99779416396SBarry Smith EXTERN_C_BEGIN 99879416396SBarry Smith #undef __FUNCT__ 99979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1000dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 100179416396SBarry Smith { 100279416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1003*e69d4d44SBarry Smith PetscErrorCode ierr; 100479416396SBarry Smith 100579416396SBarry Smith PetscFunctionBegin; 100679416396SBarry Smith jac->type = type; 10073b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10083b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10093b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1010*e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1011*e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1012*e69d4d44SBarry Smith 10133b224e63SBarry Smith } else { 10143b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10153b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1016*e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1017*e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",0,0);CHKERRQ(ierr); 10183b224e63SBarry Smith } 101979416396SBarry Smith PetscFunctionReturn(0); 102079416396SBarry Smith } 102179416396SBarry Smith EXTERN_C_END 102279416396SBarry Smith 102351f519a2SBarry Smith EXTERN_C_BEGIN 102451f519a2SBarry Smith #undef __FUNCT__ 102551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 102651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 102751f519a2SBarry Smith { 102851f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 102951f519a2SBarry Smith 103051f519a2SBarry Smith PetscFunctionBegin; 103151f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 103251f519a2SBarry Smith 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); 103351f519a2SBarry Smith jac->bs = bs; 103451f519a2SBarry Smith PetscFunctionReturn(0); 103551f519a2SBarry Smith } 103651f519a2SBarry Smith EXTERN_C_END 103751f519a2SBarry Smith 103879416396SBarry Smith #undef __FUNCT__ 103979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1040bc08b0f1SBarry Smith /*@ 104179416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 104279416396SBarry Smith 104379416396SBarry Smith Collective on PC 104479416396SBarry Smith 104579416396SBarry Smith Input Parameter: 104679416396SBarry Smith . pc - the preconditioner context 1047ce9499c7SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 104879416396SBarry Smith 104979416396SBarry Smith Options Database Key: 1050ce9499c7SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type 105179416396SBarry Smith 105279416396SBarry Smith Level: Developer 105379416396SBarry Smith 105479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 105579416396SBarry Smith 105679416396SBarry Smith .seealso: PCCompositeSetType() 105779416396SBarry Smith 105879416396SBarry Smith @*/ 1059dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 106079416396SBarry Smith { 106179416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 106279416396SBarry Smith 106379416396SBarry Smith PetscFunctionBegin; 106479416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 106579416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 106679416396SBarry Smith if (f) { 106779416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 106879416396SBarry Smith } 106979416396SBarry Smith PetscFunctionReturn(0); 107079416396SBarry Smith } 107179416396SBarry Smith 10720971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 10730971522cSBarry Smith /*MC 1074a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 10750971522cSBarry Smith fields or groups of fields 10760971522cSBarry Smith 10770971522cSBarry Smith 10780971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 1079d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 10800971522cSBarry Smith 1081a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 108269a612a9SBarry Smith and set the options directly on the resulting KSP object 10830971522cSBarry Smith 10840971522cSBarry Smith Level: intermediate 10850971522cSBarry Smith 108679416396SBarry Smith Options Database Keys: 10872e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 10882e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 10892e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 109051f519a2SBarry Smith . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 1091*e69d4d44SBarry Smith . -pc_splitfield_type <additive,multiplicative,schur,symmetric_multiplicative> 1092*e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 109379416396SBarry Smith 10943b224e63SBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -schurAblock_ and -schur, 10953b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 10963b224e63SBarry Smith 10973b224e63SBarry Smith 1098d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1099d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1100d32f9abdSBarry Smith 1101d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1102d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1103d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1104d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1105d32f9abdSBarry Smith 1106d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1107d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1108d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1109d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1110d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1111d32f9abdSBarry Smith 1112*e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1113*e69d4d44SBarry Smith ( C D ) 1114*e69d4d44SBarry Smith the preconditioner is 1115*e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1116*e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1117*e69d4d44SBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -schurAblock_. The action of 1118*e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1119*e69d4d44SBarry Smith 0 it returns the KSP associated with -schurAblock_ while field number 1 gives -schur_ KSP. By default 1120*e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1121*e69d4d44SBarry Smith option. 1122*e69d4d44SBarry Smith 11230971522cSBarry Smith Concepts: physics based preconditioners 11240971522cSBarry Smith 11250971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1126*e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 11270971522cSBarry Smith M*/ 11280971522cSBarry Smith 11290971522cSBarry Smith EXTERN_C_BEGIN 11300971522cSBarry Smith #undef __FUNCT__ 11310971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1132dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 11330971522cSBarry Smith { 11340971522cSBarry Smith PetscErrorCode ierr; 11350971522cSBarry Smith PC_FieldSplit *jac; 11360971522cSBarry Smith 11370971522cSBarry Smith PetscFunctionBegin; 113838f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 11390971522cSBarry Smith jac->bs = -1; 11400971522cSBarry Smith jac->nsplits = 0; 11413e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1142*e69d4d44SBarry Smith jac->schurpre = PETSC_TRUE; 114351f519a2SBarry Smith 11440971522cSBarry Smith pc->data = (void*)jac; 11450971522cSBarry Smith 11460971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1147421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 11480971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 11490971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 11500971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 11510971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 11520971522cSBarry Smith pc->ops->applyrichardson = 0; 11530971522cSBarry Smith 115469a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 115569a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 11560971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 11570971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1158704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1159704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 116079416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 116179416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 116251f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 116351f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 11640971522cSBarry Smith PetscFunctionReturn(0); 11650971522cSBarry Smith } 11660971522cSBarry Smith EXTERN_C_END 11670971522cSBarry Smith 11680971522cSBarry Smith 1169