1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 33c48a1e8SJed Brown #include <petscdmcomposite.h> /*I "petscdmcomposite.h" I*/ 40971522cSBarry Smith 5443836d0SMatthew G Knepley /* 6443836d0SMatthew G Knepley There is a nice discussion of block preconditioners in 7443836d0SMatthew G Knepley 8443836d0SMatthew G Knepley [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier–Stokes equations 9443836d0SMatthew G Knepley Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808 101b8e4c5fSJed Brown http://chess.cs.umd.edu/~elman/papers/tax.pdf 11443836d0SMatthew G Knepley */ 12443836d0SMatthew G Knepley 13f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 14c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0}; 15c5d2311dSJed Brown 160971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 170971522cSBarry Smith struct _PC_FieldSplitLink { 1869a612a9SBarry Smith KSP ksp; 19443836d0SMatthew G Knepley Vec x,y,z; 20db4c96c1SJed Brown char *splitname; 210971522cSBarry Smith PetscInt nfields; 225d4c12cdSJungho Lee PetscInt *fields,*fields_col; 231b9fc7fcSBarry Smith VecScatter sctx; 245d4c12cdSJungho Lee IS is,is_col; 2551f519a2SBarry Smith PC_FieldSplitLink next,previous; 260971522cSBarry Smith }; 270971522cSBarry Smith 280971522cSBarry Smith typedef struct { 2968dd23aaSBarry Smith PCCompositeType type; 30ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 31ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 32ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 3330ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 3430ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 3579416396SBarry Smith Vec *x,*y,w1,w2; 36519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 37519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3830ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 39ace3abfcSBarry Smith PetscBool issetup; 4030ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 4130ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 4230ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 43443836d0SMatthew G Knepley Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */ 44084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 45084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 46c9c6ffaaSJed Brown PCFieldSplitSchurFactType schurfactorization; 4730ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 48443836d0SMatthew G Knepley KSP kspupper; /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */ 4997bbdb24SBarry Smith PC_FieldSplitLink head; 5063ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 51c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 520971522cSBarry Smith } PC_FieldSplit; 530971522cSBarry Smith 5416913363SBarry Smith /* 5516913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5616913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5716913363SBarry Smith PC you could change this. 5816913363SBarry Smith */ 59084e4875SJed Brown 60e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 61084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 62084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 63084e4875SJed Brown { 64084e4875SJed Brown switch (jac->schurpre) { 65084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 66084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 67084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 68084e4875SJed Brown default: 69084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 70084e4875SJed Brown } 71084e4875SJed Brown } 72084e4875SJed Brown 73084e4875SJed Brown 740971522cSBarry Smith #undef __FUNCT__ 750971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 760971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 770971522cSBarry Smith { 780971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 790971522cSBarry Smith PetscErrorCode ierr; 80ace3abfcSBarry Smith PetscBool iascii; 810971522cSBarry Smith PetscInt i,j; 825a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 830971522cSBarry Smith 840971522cSBarry Smith PetscFunctionBegin; 85251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 860971522cSBarry Smith if (iascii) { 872c9966d7SBarry Smith if (jac->bs > 0) { 8851f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 892c9966d7SBarry Smith } else { 902c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 912c9966d7SBarry Smith } 92a3df900dSMatthew G Knepley if (jac->realdiagonal) { 93a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 94a3df900dSMatthew G Knepley } 9569a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 960971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 970971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 981ab39975SBarry Smith if (ilink->fields) { 990971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 10079416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1015a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 10279416396SBarry Smith if (j > 0) { 10379416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 10479416396SBarry Smith } 1055a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1060971522cSBarry Smith } 1070971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 10879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1091ab39975SBarry Smith } else { 1101ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1111ab39975SBarry Smith } 1125a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1135a9f2f41SSatish Balay ilink = ilink->next; 1140971522cSBarry Smith } 1150971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1160971522cSBarry Smith } else { 11765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1180971522cSBarry Smith } 1190971522cSBarry Smith PetscFunctionReturn(0); 1200971522cSBarry Smith } 1210971522cSBarry Smith 1220971522cSBarry Smith #undef __FUNCT__ 1233b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1243b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1253b224e63SBarry Smith { 1263b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1273b224e63SBarry Smith PetscErrorCode ierr; 128ace3abfcSBarry Smith PetscBool iascii; 1293b224e63SBarry Smith PetscInt i,j; 1303b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1313b224e63SBarry Smith 1323b224e63SBarry Smith PetscFunctionBegin; 133251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1343b224e63SBarry Smith if (iascii) { 1352c9966d7SBarry Smith if (jac->bs > 0) { 136c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1372c9966d7SBarry Smith } else { 1382c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1392c9966d7SBarry Smith } 140a3df900dSMatthew G Knepley if (jac->realdiagonal) { 141a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 142a3df900dSMatthew G Knepley } 1433e8b8b31SMatthew G Knepley switch(jac->schurpre) { 1443e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1453e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 146b02e2d75SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_DIAG: 1473e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break; 1483e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1493e8b8b31SMatthew G Knepley if (jac->schur_user) { 1503e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1513e8b8b31SMatthew G Knepley } else { 1523e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr); 1533e8b8b31SMatthew G Knepley } 1543e8b8b31SMatthew G Knepley break; 1553e8b8b31SMatthew G Knepley default: 1563e8b8b31SMatthew G Knepley SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1573e8b8b31SMatthew G Knepley } 1583b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1593b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1603b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1613b224e63SBarry Smith if (ilink->fields) { 1623b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1633b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1643b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1653b224e63SBarry Smith if (j > 0) { 1663b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1673b224e63SBarry Smith } 1683b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1693b224e63SBarry Smith } 1703b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1713b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1723b224e63SBarry Smith } else { 1733b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1743b224e63SBarry Smith } 1753b224e63SBarry Smith ilink = ilink->next; 1763b224e63SBarry Smith } 177435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1783b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 179443836d0SMatthew G Knepley ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 1803b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 181443836d0SMatthew G Knepley if (jac->kspupper != jac->head->ksp) { 182443836d0SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr); 183443836d0SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 184443836d0SMatthew G Knepley ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 185443836d0SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 186443836d0SMatthew G Knepley } 187435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1883b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18912cae6f2SJed Brown if (jac->kspschur) { 1903b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 19112cae6f2SJed Brown } else { 19212cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 19312cae6f2SJed Brown } 1943b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1953b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1963b224e63SBarry Smith } else { 19765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1983b224e63SBarry Smith } 1993b224e63SBarry Smith PetscFunctionReturn(0); 2003b224e63SBarry Smith } 2013b224e63SBarry Smith 2023b224e63SBarry Smith #undef __FUNCT__ 2036c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2046c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2056c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2066c924f48SJed Brown { 2076c924f48SJed Brown PetscErrorCode ierr; 2086c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2095d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2105d4c12cdSJungho Lee PetscBool flg,flg_col; 2115d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2126c924f48SJed Brown 2136c924f48SJed Brown PetscFunctionBegin; 2146c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2155d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2166c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2178caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2188caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2198caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2206c924f48SJed Brown nfields = jac->bs; 22129499fbbSJungho Lee nfields_col = jac->bs; 2226c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2235d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2246c924f48SJed Brown if (!flg) break; 2255d4c12cdSJungho Lee else if (flg && !flg_col) { 2265d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2275d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2285d4c12cdSJungho Lee } 2295d4c12cdSJungho Lee else { 2305d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2315d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2325d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2335d4c12cdSJungho Lee } 2346c924f48SJed Brown } 2356c924f48SJed Brown if (i > 0) { 2366c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2376c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2386c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2396c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2406c924f48SJed Brown } 2416c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2425d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2436c924f48SJed Brown PetscFunctionReturn(0); 2446c924f48SJed Brown } 2456c924f48SJed Brown 2466c924f48SJed Brown #undef __FUNCT__ 24769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 24869a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2490971522cSBarry Smith { 2500971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2510971522cSBarry Smith PetscErrorCode ierr; 2525a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2537287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 2546c924f48SJed Brown PetscInt i; 2550971522cSBarry Smith 2560971522cSBarry Smith PetscFunctionBegin; 2577287d2a3SDmitry Karpeev /* 2587287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 2597287d2a3SDmitry Karpeev Should probably be rewritten. 2607287d2a3SDmitry Karpeev */ 2617287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 262d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 263d0af7cd3SBarry Smith if (pc->dm && !stokes) { 264bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 2650784a22cSJed Brown char **fieldNames; 2667b62db95SJungho Lee IS *fields; 267e7c4fc90SDmitry Karpeev DM *dms; 268bafc1b83SMatthew G Knepley DM subdm[128]; 269bafc1b83SMatthew G Knepley PetscBool flg; 270bafc1b83SMatthew G Knepley 27116621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 272bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 273bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE; ; i++) { 274bafc1b83SMatthew G Knepley PetscInt ifields[128]; 275bafc1b83SMatthew G Knepley IS compField; 276bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 277bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 278bafc1b83SMatthew G Knepley 2798caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 280bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 281bafc1b83SMatthew G Knepley if (!flg) break; 282bafc1b83SMatthew G Knepley if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 283bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 284bafc1b83SMatthew G Knepley if (nfields == 1) { 285bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 2864d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 2874d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 288bafc1b83SMatthew G Knepley } else { 2898caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 290bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 2914d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr); 2924d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 2937287d2a3SDmitry Karpeev } 294bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 295bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 296bafc1b83SMatthew G Knepley f = ifields[j]; 2977b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 2987b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 2997b62db95SJungho Lee } 300bafc1b83SMatthew G Knepley } 301bafc1b83SMatthew G Knepley if (i == 0) { 302bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 303bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 304bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 305bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 306bafc1b83SMatthew G Knepley } 307bafc1b83SMatthew G Knepley } else { 308bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 309bafc1b83SMatthew G Knepley for (j = 0; j < i; ++j) { 310bafc1b83SMatthew G Knepley dms[j] = subdm[j]; 311bafc1b83SMatthew G Knepley } 312bafc1b83SMatthew G Knepley } 3137b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3147b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 315e7c4fc90SDmitry Karpeev if (dms) { 3168b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 317bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3187287d2a3SDmitry Karpeev const char *prefix; 3197287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr); 3207287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix); CHKERRQ(ierr); 3217b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3227b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3237287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr); 324e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 3252fa5ba8aSJed Brown } 3267b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3278b8307b2SJed Brown } 32866ffff09SJed Brown } else { 329521d7252SBarry Smith if (jac->bs <= 0) { 330704ba839SBarry Smith if (pc->pmat) { 331521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 332704ba839SBarry Smith } else { 333704ba839SBarry Smith jac->bs = 1; 334704ba839SBarry Smith } 335521d7252SBarry Smith } 336d32f9abdSBarry Smith 3377287d2a3SDmitry Karpeev ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr); 3386ce1633cSBarry Smith if (stokes) { 3396ce1633cSBarry Smith IS zerodiags,rest; 3406ce1633cSBarry Smith PetscInt nmin,nmax; 3416ce1633cSBarry Smith 3426ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3436ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3446ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 3457287d2a3SDmitry Karpeev if (jac->reset) { 3467287d2a3SDmitry Karpeev jac->head->is = rest; 3477287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 3487287d2a3SDmitry Karpeev } 3497287d2a3SDmitry Karpeev else { 3506ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 3516ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 3527287d2a3SDmitry Karpeev } 353fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 354fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 3556ce1633cSBarry Smith } else { 3567287d2a3SDmitry Karpeev if (jac->reset) 3577287d2a3SDmitry Karpeev SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 3587287d2a3SDmitry Karpeev if (!fieldsplit_default) { 359d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 360d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 3616c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 3626c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 363d32f9abdSBarry Smith } 3647287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 365d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 366db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 3676c924f48SJed Brown char splitname[8]; 3688caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 3695d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 37079416396SBarry Smith } 3715d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 372521d7252SBarry Smith } 37366ffff09SJed Brown } 3746ce1633cSBarry Smith } 375edf189efSBarry Smith } else if (jac->nsplits == 1) { 376edf189efSBarry Smith if (ilink->is) { 377edf189efSBarry Smith IS is2; 378edf189efSBarry Smith PetscInt nmin,nmax; 379edf189efSBarry Smith 380edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 381edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 382db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 383fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 3847b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 385edf189efSBarry Smith } 386d0af7cd3SBarry Smith 387d0af7cd3SBarry Smith 3887b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 38969a612a9SBarry Smith PetscFunctionReturn(0); 39069a612a9SBarry Smith } 39169a612a9SBarry Smith 392514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 393514bf10dSMatthew G Knepley 39469a612a9SBarry Smith #undef __FUNCT__ 39569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 39669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 39769a612a9SBarry Smith { 39869a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 39969a612a9SBarry Smith PetscErrorCode ierr; 4005a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4012c9966d7SBarry Smith PetscInt i,nsplit; 40269a612a9SBarry Smith MatStructure flag = pc->flag; 4035f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 40469a612a9SBarry Smith 40569a612a9SBarry Smith PetscFunctionBegin; 40669a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 40797bbdb24SBarry Smith nsplit = jac->nsplits; 4085a9f2f41SSatish Balay ilink = jac->head; 40997bbdb24SBarry Smith 41097bbdb24SBarry Smith /* get the matrices for each split */ 411704ba839SBarry Smith if (!jac->issetup) { 4121b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 41397bbdb24SBarry Smith 414704ba839SBarry Smith jac->issetup = PETSC_TRUE; 415704ba839SBarry Smith 4165d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4172c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4182c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4192c9966d7SBarry Smith } 42051f519a2SBarry Smith bs = jac->bs; 42197bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4221b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4231b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4241b9fc7fcSBarry Smith if (jac->defaultsplit) { 425704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4265f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 427704ba839SBarry Smith } else if (!ilink->is) { 428ccb205f8SBarry Smith if (ilink->nfields > 1) { 4295f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4305a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4315f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4321b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4331b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4341b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4355f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 43697bbdb24SBarry Smith } 43797bbdb24SBarry Smith } 438d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4395f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 440ccb205f8SBarry Smith } else { 441704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4425f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 443ccb205f8SBarry Smith } 4443e197d65SBarry Smith } 445704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4465f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4475f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 448704ba839SBarry Smith ilink = ilink->next; 4491b9fc7fcSBarry Smith } 4501b9fc7fcSBarry Smith } 4511b9fc7fcSBarry Smith 452704ba839SBarry Smith ilink = jac->head; 45397bbdb24SBarry Smith if (!jac->pmat) { 454bdddcaaaSMatthew G Knepley Vec xtmp; 455bdddcaaaSMatthew G Knepley 456bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 457cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 458bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 459cf502942SBarry Smith for (i=0; i<nsplit; i++) { 460bdddcaaaSMatthew G Knepley MatNullSpace sp; 461bdddcaaaSMatthew G Knepley 462a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 463a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 464a3df900dSMatthew G Knepley if (jac->pmat[i]) { 465a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 466a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 467a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 468a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 469a3df900dSMatthew G Knepley } 470a3df900dSMatthew G Knepley } else { 4715f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 472a3df900dSMatthew G Knepley } 473bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 474bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 475443836d0SMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = PETSC_NULL; 476bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 477bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 478ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 479bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr); 480bafc1b83SMatthew G Knepley if (sp) { 481bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 482bafc1b83SMatthew G Knepley } 483ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 484ed1f0337SMatthew G Knepley if (sp) { 485ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 486ed1f0337SMatthew G Knepley } 487704ba839SBarry Smith ilink = ilink->next; 488cf502942SBarry Smith } 489bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 49097bbdb24SBarry Smith } else { 491cf502942SBarry Smith for (i=0; i<nsplit; i++) { 492a3df900dSMatthew G Knepley Mat pmat; 493a3df900dSMatthew G Knepley 494a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 495a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 496a3df900dSMatthew G Knepley if (!pmat) { 4975f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 498a3df900dSMatthew G Knepley } 499704ba839SBarry Smith ilink = ilink->next; 500cf502942SBarry Smith } 50197bbdb24SBarry Smith } 502519d70e2SJed Brown if (jac->realdiagonal) { 503519d70e2SJed Brown ilink = jac->head; 504519d70e2SJed Brown if (!jac->mat) { 505519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 506519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5075f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 508519d70e2SJed Brown ilink = ilink->next; 509519d70e2SJed Brown } 510519d70e2SJed Brown } else { 511519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5125f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 513519d70e2SJed Brown ilink = ilink->next; 514519d70e2SJed Brown } 515519d70e2SJed Brown } 516519d70e2SJed Brown } else { 517519d70e2SJed Brown jac->mat = jac->pmat; 518519d70e2SJed Brown } 51997bbdb24SBarry Smith 5206c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 52168dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 52268dd23aaSBarry Smith ilink = jac->head; 52368dd23aaSBarry Smith if (!jac->Afield) { 52468dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 52568dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5264aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 52768dd23aaSBarry Smith ilink = ilink->next; 52868dd23aaSBarry Smith } 52968dd23aaSBarry Smith } else { 53068dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5314aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 53268dd23aaSBarry Smith ilink = ilink->next; 53368dd23aaSBarry Smith } 53468dd23aaSBarry Smith } 53568dd23aaSBarry Smith } 53668dd23aaSBarry Smith 5373b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5383b224e63SBarry Smith IS ccis; 5394aa3045dSJed Brown PetscInt rstart,rend; 540093c86ffSJed Brown char lscname[256]; 541093c86ffSJed Brown PetscObject LSC_L; 542e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 54368dd23aaSBarry Smith 544e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 545e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 546e6cab6aaSJed Brown 5473b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5483b224e63SBarry Smith if (jac->schur) { 549*a3314f2cSMatthew G Knepley KSP kspA = jac->head->ksp, kspInner = PETSC_NULL, kspUpper = jac->kspupper; 550443836d0SMatthew G Knepley 551fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 5523b224e63SBarry Smith ilink = jac->head; 55349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5544aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 555fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5563b224e63SBarry Smith ilink = ilink->next; 55749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5584aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 559fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 560a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 561443836d0SMatthew G Knepley if (kspA != kspInner) { 562443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 563443836d0SMatthew G Knepley } 564443836d0SMatthew G Knepley if (kspUpper != kspA) { 565443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 566443836d0SMatthew G Knepley } 567084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 5683b224e63SBarry Smith } else { 5691cee3971SBarry Smith KSP ksp; 570bafc1b83SMatthew G Knepley const char *Dprefix; 571bafc1b83SMatthew G Knepley char schurprefix[256]; 572514bf10dSMatthew G Knepley char schurtestoption[256]; 573bdddcaaaSMatthew G Knepley MatNullSpace sp; 574514bf10dSMatthew G Knepley PetscBool flg; 5753b224e63SBarry Smith 576a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 5773b224e63SBarry Smith ilink = jac->head; 57849bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5794aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 580fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5813b224e63SBarry Smith ilink = ilink->next; 58249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5834aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 584fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 58520252d06SBarry Smith 58620252d06SBarry Smith /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */ 58720252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 58820252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 58920252d06SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr); 59020252d06SBarry Smith ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 59120252d06SBarry Smith /* Indent this deeper to emphasize the "inner" nature of this solver. */ 59220252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr); 59320252d06SBarry Smith ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr); 59420252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 59520252d06SBarry Smith 596bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 59720252d06SBarry Smith if (sp) { 59820252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 59920252d06SBarry Smith } 60020252d06SBarry Smith 60120252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 602514bf10dSMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 603514bf10dSMatthew G Knepley if (flg) { 604514bf10dSMatthew G Knepley DM dmInner; 605514bf10dSMatthew G Knepley 606514bf10dSMatthew G Knepley /* Set DM for new solver */ 607bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 608bafc1b83SMatthew G Knepley ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr); 6097287d2a3SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 610443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 611443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 612514bf10dSMatthew G Knepley } else { 613514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr); 614514bf10dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 615bafc1b83SMatthew G Knepley } 61620b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 6173b224e63SBarry Smith 618443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 619443836d0SMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 620443836d0SMatthew G Knepley if (flg) { 621443836d0SMatthew G Knepley DM dmInner; 622443836d0SMatthew G Knepley 623443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 624443836d0SMatthew G Knepley ierr = KSPCreate(((PetscObject) pc)->comm, &jac->kspupper);CHKERRQ(ierr); 625443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 626443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 627443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 628443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 629443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 630443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 631443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 632443836d0SMatthew G Knepley } else { 633443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 634443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 635443836d0SMatthew G Knepley } 636443836d0SMatthew G Knepley 6373b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur); CHKERRQ(ierr); 6389005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 63920252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1); CHKERRQ(ierr); 640084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 641084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 6427233a360SDmitry Karpeev PC pcschur; 6437233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 6447233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 645084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 646e69d4d44SBarry Smith } 6477287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr); 6487287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix); CHKERRQ(ierr); 6493b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 65020b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 65120b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 6523b224e63SBarry Smith } 653093c86ffSJed Brown 654093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 6558caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 656093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 657093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 658093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 6598caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 660093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 661093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 662093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 6633b224e63SBarry Smith } else { 66468bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 66597bbdb24SBarry Smith i = 0; 6665a9f2f41SSatish Balay ilink = jac->head; 6675a9f2f41SSatish Balay while (ilink) { 668519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 6693b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 670c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 67197bbdb24SBarry Smith i++; 6725a9f2f41SSatish Balay ilink = ilink->next; 6730971522cSBarry Smith } 6743b224e63SBarry Smith } 6753b224e63SBarry Smith 676c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 6770971522cSBarry Smith PetscFunctionReturn(0); 6780971522cSBarry Smith } 6790971522cSBarry Smith 6805a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 681ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 682ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 6835a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 684ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 685ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 68679416396SBarry Smith 6870971522cSBarry Smith #undef __FUNCT__ 6883b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 6893b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 6903b224e63SBarry Smith { 6913b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6923b224e63SBarry Smith PetscErrorCode ierr; 6933b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 694443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 6953b224e63SBarry Smith 6963b224e63SBarry Smith PetscFunctionBegin; 697c5d2311dSJed Brown switch (jac->schurfactorization) { 698c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 699a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 700c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 701c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 702c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 703443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 704c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 705c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 706c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 707c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 708c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 709c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 710c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 711c5d2311dSJed Brown break; 712c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 713a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 714c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 715c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 716443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 717c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 718c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 719c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 720c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 721c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 722c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 723c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 724c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 725c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 726c5d2311dSJed Brown break; 727c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 728a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 729c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 730c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 731c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 732c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 733c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 734c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 735c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 736c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 737443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 738c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 739c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 740c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 741c5d2311dSJed Brown break; 742c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 743a04f6461SBarry Smith /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1], an exact solve if applied exactly, needs one extra solve with A */ 7443b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7453b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 746443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 7473b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 7483b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 7493b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7503b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7513b224e63SBarry Smith 7523b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 7533b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7543b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7553b224e63SBarry Smith 756443836d0SMatthew G Knepley if (kspUpper == kspA) { 7573b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 7583b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 759443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 760443836d0SMatthew G Knepley } else { 761443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 762443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 763443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 764443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 765443836d0SMatthew G Knepley } 7663b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7673b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 768c5d2311dSJed Brown } 7693b224e63SBarry Smith PetscFunctionReturn(0); 7703b224e63SBarry Smith } 7713b224e63SBarry Smith 7723b224e63SBarry Smith #undef __FUNCT__ 7730971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 7740971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 7750971522cSBarry Smith { 7760971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7770971522cSBarry Smith PetscErrorCode ierr; 7785a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 779939b8a20SBarry Smith PetscInt cnt,bs; 7800971522cSBarry Smith 7810971522cSBarry Smith PetscFunctionBegin; 7824442daceSBarry Smith x->map->bs = jac->bs; 7834442daceSBarry Smith y->map->bs = jac->bs; 78451f519a2SBarry Smith CHKMEMQ; 78579416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 7861b9fc7fcSBarry Smith if (jac->defaultsplit) { 787939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 788939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 789939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 790939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 7910971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 7925a9f2f41SSatish Balay while (ilink) { 7935a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7945a9f2f41SSatish Balay ilink = ilink->next; 7950971522cSBarry Smith } 7960971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 7971b9fc7fcSBarry Smith } else { 798efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7995a9f2f41SSatish Balay while (ilink) { 8005a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8015a9f2f41SSatish Balay ilink = ilink->next; 8021b9fc7fcSBarry Smith } 8031b9fc7fcSBarry Smith } 80416913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 80579416396SBarry Smith if (!jac->w1) { 80679416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 80779416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 80879416396SBarry Smith } 809efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8105a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8113e197d65SBarry Smith cnt = 1; 8125a9f2f41SSatish Balay while (ilink->next) { 8135a9f2f41SSatish Balay ilink = ilink->next; 8143e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8153e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8163e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8173e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8183e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8193e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8203e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8213e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8223e197d65SBarry Smith } 82351f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 82411755939SBarry Smith cnt -= 2; 82551f519a2SBarry Smith while (ilink->previous) { 82651f519a2SBarry Smith ilink = ilink->previous; 82711755939SBarry Smith /* compute the residual only over the part of the vector needed */ 82811755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 82911755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 83011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 83111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 83211755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 83311755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 83411755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 83551f519a2SBarry Smith } 83611755939SBarry Smith } 83765e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 83851f519a2SBarry Smith CHKMEMQ; 8390971522cSBarry Smith PetscFunctionReturn(0); 8400971522cSBarry Smith } 8410971522cSBarry Smith 842421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 843ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 844ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 845421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 846ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 847ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 848421e10b8SBarry Smith 849421e10b8SBarry Smith #undef __FUNCT__ 8508c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 851421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 852421e10b8SBarry Smith { 853421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 854421e10b8SBarry Smith PetscErrorCode ierr; 855421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 856939b8a20SBarry Smith PetscInt bs; 857421e10b8SBarry Smith 858421e10b8SBarry Smith PetscFunctionBegin; 859421e10b8SBarry Smith CHKMEMQ; 860421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 861421e10b8SBarry Smith if (jac->defaultsplit) { 862939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 863939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 864939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 865939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 866421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 867421e10b8SBarry Smith while (ilink) { 868421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 869421e10b8SBarry Smith ilink = ilink->next; 870421e10b8SBarry Smith } 871421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 872421e10b8SBarry Smith } else { 873421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 874421e10b8SBarry Smith while (ilink) { 875421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 876421e10b8SBarry Smith ilink = ilink->next; 877421e10b8SBarry Smith } 878421e10b8SBarry Smith } 879421e10b8SBarry Smith } else { 880421e10b8SBarry Smith if (!jac->w1) { 881421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 882421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 883421e10b8SBarry Smith } 884421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 885421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 886421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 887421e10b8SBarry Smith while (ilink->next) { 888421e10b8SBarry Smith ilink = ilink->next; 8899989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 890421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 891421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 892421e10b8SBarry Smith } 893421e10b8SBarry Smith while (ilink->previous) { 894421e10b8SBarry Smith ilink = ilink->previous; 8959989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 896421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 897421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 898421e10b8SBarry Smith } 899421e10b8SBarry Smith } else { 900421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 901421e10b8SBarry Smith ilink = ilink->next; 902421e10b8SBarry Smith } 903421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 904421e10b8SBarry Smith while (ilink->previous) { 905421e10b8SBarry Smith ilink = ilink->previous; 9069989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 907421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 908421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 909421e10b8SBarry Smith } 910421e10b8SBarry Smith } 911421e10b8SBarry Smith } 912421e10b8SBarry Smith CHKMEMQ; 913421e10b8SBarry Smith PetscFunctionReturn(0); 914421e10b8SBarry Smith } 915421e10b8SBarry Smith 9160971522cSBarry Smith #undef __FUNCT__ 917574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 918574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9190971522cSBarry Smith { 9200971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9210971522cSBarry Smith PetscErrorCode ierr; 9225a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9230971522cSBarry Smith 9240971522cSBarry Smith PetscFunctionBegin; 9255a9f2f41SSatish Balay while (ilink) { 926574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 927fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 928fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 929443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 930fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 931fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 932b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9335a9f2f41SSatish Balay next = ilink->next; 9345a9f2f41SSatish Balay ilink = next; 9350971522cSBarry Smith } 93605b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 937574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 938574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 939574deadeSBarry Smith } else if (jac->mat) { 940574deadeSBarry Smith jac->mat = PETSC_NULL; 941574deadeSBarry Smith } 94297bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 94368dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 9446bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 9456bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 9466bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 9476bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 9486bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 9496bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 9506bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 95163ec74ffSBarry Smith jac->reset = PETSC_TRUE; 952574deadeSBarry Smith PetscFunctionReturn(0); 953574deadeSBarry Smith } 954574deadeSBarry Smith 955574deadeSBarry Smith #undef __FUNCT__ 956574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 957574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 958574deadeSBarry Smith { 959574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 960574deadeSBarry Smith PetscErrorCode ierr; 961574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 962574deadeSBarry Smith 963574deadeSBarry Smith PetscFunctionBegin; 964574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 965574deadeSBarry Smith while (ilink) { 9666bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 967574deadeSBarry Smith next = ilink->next; 968574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 969574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 9705d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 971574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 972574deadeSBarry Smith ilink = next; 973574deadeSBarry Smith } 974574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 975c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 9767b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 9777b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 9787b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 9797b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 9807b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 981ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 982c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 9830971522cSBarry Smith PetscFunctionReturn(0); 9840971522cSBarry Smith } 9850971522cSBarry Smith 9860971522cSBarry Smith #undef __FUNCT__ 9870971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 9880971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 9890971522cSBarry Smith { 9901b9fc7fcSBarry Smith PetscErrorCode ierr; 9916c924f48SJed Brown PetscInt bs; 992bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 9939dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9943b224e63SBarry Smith PCCompositeType ctype; 995e7c4fc90SDmitry Karpeev DM ddm; 996e7c4fc90SDmitry Karpeev char ddm_name[1025]; 9971b9fc7fcSBarry Smith 9980971522cSBarry Smith PetscFunctionBegin; 9991b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 1000e7c4fc90SDmitry Karpeev if (pc->dm) { 1001e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 1002e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 1003731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 1004e7c4fc90SDmitry Karpeev if (flg) { 100516621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 1006e7c4fc90SDmitry Karpeev if (!ddm) { 100716621825SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 1008e7c4fc90SDmitry Karpeev } 100916621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 1010e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 1011e7c4fc90SDmitry Karpeev } 1012e7c4fc90SDmitry Karpeev } 1013acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 101451f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 101551f519a2SBarry Smith if (flg) { 101651f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 101751f519a2SBarry Smith } 1018704ba839SBarry Smith 1019435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 1020c0adfefeSBarry Smith if (stokes) { 1021c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1022c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1023c0adfefeSBarry Smith } 1024c0adfefeSBarry Smith 10253b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10263b224e63SBarry Smith if (flg) { 10273b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10283b224e63SBarry Smith } 1029c30613efSMatthew Knepley /* Only setup fields once */ 1030c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1031d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1032d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10336c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10346c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1035d32f9abdSBarry Smith } 1036c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1037c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1038c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 1039c9c6ffaaSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_fact_type","Which off-diagonal parts of the block factorization to use","PCFieldSplitSetSchurFactType",PCFieldSplitSchurFactTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 1040084e4875SJed Brown 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); 1041c5d2311dSJed Brown } 10421b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10430971522cSBarry Smith PetscFunctionReturn(0); 10440971522cSBarry Smith } 10450971522cSBarry Smith 10460971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10470971522cSBarry Smith 10480971522cSBarry Smith EXTERN_C_BEGIN 10490971522cSBarry Smith #undef __FUNCT__ 10500971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 10515d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10520971522cSBarry Smith { 105397bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10540971522cSBarry Smith PetscErrorCode ierr; 10555a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 105669a612a9SBarry Smith char prefix[128]; 10575d4c12cdSJungho Lee PetscInt i; 10580971522cSBarry Smith 10590971522cSBarry Smith PetscFunctionBegin; 10606c924f48SJed Brown if (jac->splitdefined) { 10616c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10626c924f48SJed Brown PetscFunctionReturn(0); 10636c924f48SJed Brown } 106451f519a2SBarry Smith for (i=0; i<n; i++) { 1065e32f2f54SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 1066e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 106751f519a2SBarry Smith } 1068704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1069a04f6461SBarry Smith if (splitname) { 1070db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1071a04f6461SBarry Smith } else { 1072a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1073a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1074a04f6461SBarry Smith } 1075704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 10765a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 10775d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 10785d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 10795a9f2f41SSatish Balay ilink->nfields = n; 10805a9f2f41SSatish Balay ilink->next = PETSC_NULL; 10817adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 108220252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 10835a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10849005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 108569a612a9SBarry Smith 10868caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 10875a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 10880971522cSBarry Smith 10890971522cSBarry Smith if (!next) { 10905a9f2f41SSatish Balay jac->head = ilink; 109151f519a2SBarry Smith ilink->previous = PETSC_NULL; 10920971522cSBarry Smith } else { 10930971522cSBarry Smith while (next->next) { 10940971522cSBarry Smith next = next->next; 10950971522cSBarry Smith } 10965a9f2f41SSatish Balay next->next = ilink; 109751f519a2SBarry Smith ilink->previous = next; 10980971522cSBarry Smith } 10990971522cSBarry Smith jac->nsplits++; 11000971522cSBarry Smith PetscFunctionReturn(0); 11010971522cSBarry Smith } 11020971522cSBarry Smith EXTERN_C_END 11030971522cSBarry Smith 1104e69d4d44SBarry Smith EXTERN_C_BEGIN 1105e69d4d44SBarry Smith #undef __FUNCT__ 1106e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11077087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1108e69d4d44SBarry Smith { 1109e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1110e69d4d44SBarry Smith PetscErrorCode ierr; 1111e69d4d44SBarry Smith 1112e69d4d44SBarry Smith PetscFunctionBegin; 1113519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1114e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1115e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 111613e0d083SBarry Smith if (n) *n = jac->nsplits; 1117e69d4d44SBarry Smith PetscFunctionReturn(0); 1118e69d4d44SBarry Smith } 1119e69d4d44SBarry Smith EXTERN_C_END 11200971522cSBarry Smith 11210971522cSBarry Smith EXTERN_C_BEGIN 11220971522cSBarry Smith #undef __FUNCT__ 112369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11247087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11250971522cSBarry Smith { 11260971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11270971522cSBarry Smith PetscErrorCode ierr; 11280971522cSBarry Smith PetscInt cnt = 0; 11295a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11300971522cSBarry Smith 11310971522cSBarry Smith PetscFunctionBegin; 11325d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11335a9f2f41SSatish Balay while (ilink) { 11345a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11355a9f2f41SSatish Balay ilink = ilink->next; 11360971522cSBarry Smith } 11375d480477SMatthew G Knepley if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number of splits in linked list %D does not match number in object %D",cnt,jac->nsplits); 113813e0d083SBarry Smith if (n) *n = jac->nsplits; 11390971522cSBarry Smith PetscFunctionReturn(0); 11400971522cSBarry Smith } 11410971522cSBarry Smith EXTERN_C_END 11420971522cSBarry Smith 1143704ba839SBarry Smith EXTERN_C_BEGIN 1144704ba839SBarry Smith #undef __FUNCT__ 1145704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11467087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1147704ba839SBarry Smith { 1148704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1149704ba839SBarry Smith PetscErrorCode ierr; 1150704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1151704ba839SBarry Smith char prefix[128]; 1152704ba839SBarry Smith 1153704ba839SBarry Smith PetscFunctionBegin; 11546c924f48SJed Brown if (jac->splitdefined) { 11556c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11566c924f48SJed Brown PetscFunctionReturn(0); 11576c924f48SJed Brown } 115816913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1159a04f6461SBarry Smith if (splitname) { 1160db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1161a04f6461SBarry Smith } else { 1162b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1163b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1164a04f6461SBarry Smith } 11651ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1166b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1167b5787286SJed Brown ilink->is = is; 1168b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1169b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1170b5787286SJed Brown ilink->is_col = is; 1171704ba839SBarry Smith ilink->next = PETSC_NULL; 1172704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 117320252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1174704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11759005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1176704ba839SBarry Smith 11778caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1178704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1179704ba839SBarry Smith 1180704ba839SBarry Smith if (!next) { 1181704ba839SBarry Smith jac->head = ilink; 1182704ba839SBarry Smith ilink->previous = PETSC_NULL; 1183704ba839SBarry Smith } else { 1184704ba839SBarry Smith while (next->next) { 1185704ba839SBarry Smith next = next->next; 1186704ba839SBarry Smith } 1187704ba839SBarry Smith next->next = ilink; 1188704ba839SBarry Smith ilink->previous = next; 1189704ba839SBarry Smith } 1190704ba839SBarry Smith jac->nsplits++; 1191704ba839SBarry Smith 1192704ba839SBarry Smith PetscFunctionReturn(0); 1193704ba839SBarry Smith } 1194704ba839SBarry Smith EXTERN_C_END 1195704ba839SBarry Smith 11960971522cSBarry Smith #undef __FUNCT__ 11970971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 11980971522cSBarry Smith /*@ 11990971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12000971522cSBarry Smith 1201ad4df100SBarry Smith Logically Collective on PC 12020971522cSBarry Smith 12030971522cSBarry Smith Input Parameters: 12040971522cSBarry Smith + pc - the preconditioner context 1205a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 12060971522cSBarry Smith . n - the number of fields in this split 1207db4c96c1SJed Brown - fields - the fields in this split 12080971522cSBarry Smith 12090971522cSBarry Smith Level: intermediate 12100971522cSBarry Smith 1211d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1212d32f9abdSBarry Smith 12137287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1214d32f9abdSBarry 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 1215d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1216d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1217d32f9abdSBarry Smith 1218db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1219db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1220db4c96c1SJed Brown 12215d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12225d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12235d4c12cdSJungho Lee available when this routine is called. 12245d4c12cdSJungho Lee 1225d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12260971522cSBarry Smith 12270971522cSBarry Smith @*/ 12285d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12290971522cSBarry Smith { 12304ac538c5SBarry Smith PetscErrorCode ierr; 12310971522cSBarry Smith 12320971522cSBarry Smith PetscFunctionBegin; 12330700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1234db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1235db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1236db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12375d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12380971522cSBarry Smith PetscFunctionReturn(0); 12390971522cSBarry Smith } 12400971522cSBarry Smith 12410971522cSBarry Smith #undef __FUNCT__ 1242704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1243704ba839SBarry Smith /*@ 1244704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1245704ba839SBarry Smith 1246ad4df100SBarry Smith Logically Collective on PC 1247704ba839SBarry Smith 1248704ba839SBarry Smith Input Parameters: 1249704ba839SBarry Smith + pc - the preconditioner context 1250a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1251db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1252704ba839SBarry Smith 1253d32f9abdSBarry Smith 1254a6ffb8dbSJed Brown Notes: 1255a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1256a6ffb8dbSJed Brown 1257db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1258db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1259d32f9abdSBarry Smith 1260704ba839SBarry Smith Level: intermediate 1261704ba839SBarry Smith 1262704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1263704ba839SBarry Smith 1264704ba839SBarry Smith @*/ 12657087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1266704ba839SBarry Smith { 12674ac538c5SBarry Smith PetscErrorCode ierr; 1268704ba839SBarry Smith 1269704ba839SBarry Smith PetscFunctionBegin; 12700700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12717b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1272db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 12734ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1274704ba839SBarry Smith PetscFunctionReturn(0); 1275704ba839SBarry Smith } 1276704ba839SBarry Smith 1277704ba839SBarry Smith #undef __FUNCT__ 127857a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 127957a9adfeSMatthew G Knepley /*@ 128057a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 128157a9adfeSMatthew G Knepley 128257a9adfeSMatthew G Knepley Logically Collective on PC 128357a9adfeSMatthew G Knepley 128457a9adfeSMatthew G Knepley Input Parameters: 128557a9adfeSMatthew G Knepley + pc - the preconditioner context 128657a9adfeSMatthew G Knepley - splitname - name of this split 128757a9adfeSMatthew G Knepley 128857a9adfeSMatthew G Knepley Output Parameter: 128957a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 129057a9adfeSMatthew G Knepley 129157a9adfeSMatthew G Knepley Level: intermediate 129257a9adfeSMatthew G Knepley 129357a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 129457a9adfeSMatthew G Knepley 129557a9adfeSMatthew G Knepley @*/ 129657a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 129757a9adfeSMatthew G Knepley { 129857a9adfeSMatthew G Knepley PetscErrorCode ierr; 129957a9adfeSMatthew G Knepley 130057a9adfeSMatthew G Knepley PetscFunctionBegin; 130157a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 130257a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 130357a9adfeSMatthew G Knepley PetscValidPointer(is,3); 130457a9adfeSMatthew G Knepley { 130557a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 130657a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 130757a9adfeSMatthew G Knepley PetscBool found; 130857a9adfeSMatthew G Knepley 130957a9adfeSMatthew G Knepley *is = PETSC_NULL; 131057a9adfeSMatthew G Knepley while(ilink) { 131157a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 131257a9adfeSMatthew G Knepley if (found) { 131357a9adfeSMatthew G Knepley *is = ilink->is; 131457a9adfeSMatthew G Knepley break; 131557a9adfeSMatthew G Knepley } 131657a9adfeSMatthew G Knepley ilink = ilink->next; 131757a9adfeSMatthew G Knepley } 131857a9adfeSMatthew G Knepley } 131957a9adfeSMatthew G Knepley PetscFunctionReturn(0); 132057a9adfeSMatthew G Knepley } 132157a9adfeSMatthew G Knepley 132257a9adfeSMatthew G Knepley #undef __FUNCT__ 132351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 132451f519a2SBarry Smith /*@ 132551f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 132651f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 132751f519a2SBarry Smith 1328ad4df100SBarry Smith Logically Collective on PC 132951f519a2SBarry Smith 133051f519a2SBarry Smith Input Parameters: 133151f519a2SBarry Smith + pc - the preconditioner context 133251f519a2SBarry Smith - bs - the block size 133351f519a2SBarry Smith 133451f519a2SBarry Smith Level: intermediate 133551f519a2SBarry Smith 133651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 133751f519a2SBarry Smith 133851f519a2SBarry Smith @*/ 13397087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 134051f519a2SBarry Smith { 13414ac538c5SBarry Smith PetscErrorCode ierr; 134251f519a2SBarry Smith 134351f519a2SBarry Smith PetscFunctionBegin; 13440700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1345c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13464ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 134751f519a2SBarry Smith PetscFunctionReturn(0); 134851f519a2SBarry Smith } 134951f519a2SBarry Smith 135051f519a2SBarry Smith #undef __FUNCT__ 135169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 13520971522cSBarry Smith /*@C 135369a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 13540971522cSBarry Smith 135569a612a9SBarry Smith Collective on KSP 13560971522cSBarry Smith 13570971522cSBarry Smith Input Parameter: 13580971522cSBarry Smith . pc - the preconditioner context 13590971522cSBarry Smith 13600971522cSBarry Smith Output Parameters: 136113e0d083SBarry Smith + n - the number of splits 136269a612a9SBarry Smith - pc - the array of KSP contexts 13630971522cSBarry Smith 13640971522cSBarry Smith Note: 1365d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1366d32f9abdSBarry Smith (not the KSP just the array that contains them). 13670971522cSBarry Smith 136869a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 13690971522cSBarry Smith 1370196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 1371196cc216SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the 1372196cc216SBarry Smith KSP array must be. 1373196cc216SBarry Smith 1374196cc216SBarry Smith 13750971522cSBarry Smith Level: advanced 13760971522cSBarry Smith 13770971522cSBarry Smith .seealso: PCFIELDSPLIT 13780971522cSBarry Smith @*/ 13797087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 13800971522cSBarry Smith { 13814ac538c5SBarry Smith PetscErrorCode ierr; 13820971522cSBarry Smith 13830971522cSBarry Smith PetscFunctionBegin; 13840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 138513e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 13864ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 13870971522cSBarry Smith PetscFunctionReturn(0); 13880971522cSBarry Smith } 13890971522cSBarry Smith 1390e69d4d44SBarry Smith #undef __FUNCT__ 1391e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1392e69d4d44SBarry Smith /*@ 1393e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1394a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1395e69d4d44SBarry Smith 1396e69d4d44SBarry Smith Collective on PC 1397e69d4d44SBarry Smith 1398e69d4d44SBarry Smith Input Parameters: 1399e69d4d44SBarry Smith + pc - the preconditioner context 1400fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1401084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1402084e4875SJed Brown 1403e69d4d44SBarry Smith Options Database: 1404084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1405e69d4d44SBarry Smith 1406fd1303e9SJungho Lee Notes: 1407fd1303e9SJungho Lee $ If ptype is 1408fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1409fd1303e9SJungho Lee $ to this function). 1410fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1411fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1412fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1413fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1414fd1303e9SJungho Lee $ preconditioner 1415fd1303e9SJungho Lee 1416fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1417fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1418fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1419fd1303e9SJungho Lee 1420fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1421fd1303e9SJungho Lee the name since it sets a proceedure to use. 1422fd1303e9SJungho Lee 1423e69d4d44SBarry Smith Level: intermediate 1424e69d4d44SBarry Smith 1425fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1426e69d4d44SBarry Smith 1427e69d4d44SBarry Smith @*/ 14287087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1429e69d4d44SBarry Smith { 14304ac538c5SBarry Smith PetscErrorCode ierr; 1431e69d4d44SBarry Smith 1432e69d4d44SBarry Smith PetscFunctionBegin; 14330700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14344ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1435e69d4d44SBarry Smith PetscFunctionReturn(0); 1436e69d4d44SBarry Smith } 1437e69d4d44SBarry Smith 1438e69d4d44SBarry Smith EXTERN_C_BEGIN 1439e69d4d44SBarry Smith #undef __FUNCT__ 1440e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14417087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1442e69d4d44SBarry Smith { 1443e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1444084e4875SJed Brown PetscErrorCode ierr; 1445e69d4d44SBarry Smith 1446e69d4d44SBarry Smith PetscFunctionBegin; 1447084e4875SJed Brown jac->schurpre = ptype; 1448084e4875SJed Brown if (pre) { 14496bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1450084e4875SJed Brown jac->schur_user = pre; 1451084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1452084e4875SJed Brown } 1453e69d4d44SBarry Smith PetscFunctionReturn(0); 1454e69d4d44SBarry Smith } 1455e69d4d44SBarry Smith EXTERN_C_END 1456e69d4d44SBarry Smith 145730ad9308SMatthew Knepley #undef __FUNCT__ 1458c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1459ab1df9f5SJed Brown /*@ 1460c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1461ab1df9f5SJed Brown 1462ab1df9f5SJed Brown Collective on PC 1463ab1df9f5SJed Brown 1464ab1df9f5SJed Brown Input Parameters: 1465ab1df9f5SJed Brown + pc - the preconditioner context 1466c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1467ab1df9f5SJed Brown 1468ab1df9f5SJed Brown Options Database: 1469c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1470ab1df9f5SJed Brown 1471ab1df9f5SJed Brown 1472ab1df9f5SJed Brown Level: intermediate 1473ab1df9f5SJed Brown 1474ab1df9f5SJed Brown Notes: 1475ab1df9f5SJed Brown The FULL factorization is 1476ab1df9f5SJed Brown 1477ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1478ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1479ab1df9f5SJed Brown 14806be4592eSBarry Smith where S = D - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D, 14816be4592eSBarry Smith and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES). 1482ab1df9f5SJed Brown 14836be4592eSBarry Smith If applied exactly, FULL factorization is a direct solver. The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial 14846be4592eSBarry Smith of degree 2, so KSPGMRES converges in 2 iterations. If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner 14856be4592eSBarry Smith application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice. With DIAG, 14866be4592eSBarry Smith the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations. 1487ab1df9f5SJed Brown 14886be4592eSBarry Smith For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES. Note that a flexible method like KSPFGMRES 14896be4592eSBarry Smith or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used inside a split). 1490ab1df9f5SJed Brown 1491ab1df9f5SJed Brown References: 1492ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1493ab1df9f5SJed Brown 1494ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1495ab1df9f5SJed Brown 1496ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1497ab1df9f5SJed Brown @*/ 1498c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1499ab1df9f5SJed Brown { 1500ab1df9f5SJed Brown PetscErrorCode ierr; 1501ab1df9f5SJed Brown 1502ab1df9f5SJed Brown PetscFunctionBegin; 1503ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1504c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1505ab1df9f5SJed Brown PetscFunctionReturn(0); 1506ab1df9f5SJed Brown } 1507ab1df9f5SJed Brown 1508ab1df9f5SJed Brown #undef __FUNCT__ 1509c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1510c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1511ab1df9f5SJed Brown { 1512ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1513ab1df9f5SJed Brown 1514ab1df9f5SJed Brown PetscFunctionBegin; 1515ab1df9f5SJed Brown jac->schurfactorization = ftype; 1516ab1df9f5SJed Brown PetscFunctionReturn(0); 1517ab1df9f5SJed Brown } 1518ab1df9f5SJed Brown 1519ab1df9f5SJed Brown #undef __FUNCT__ 152030ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 152130ad9308SMatthew Knepley /*@C 15228c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 152330ad9308SMatthew Knepley 152430ad9308SMatthew Knepley Collective on KSP 152530ad9308SMatthew Knepley 152630ad9308SMatthew Knepley Input Parameter: 152730ad9308SMatthew Knepley . pc - the preconditioner context 152830ad9308SMatthew Knepley 152930ad9308SMatthew Knepley Output Parameters: 1530a04f6461SBarry Smith + A00 - the (0,0) block 1531a04f6461SBarry Smith . A01 - the (0,1) block 1532a04f6461SBarry Smith . A10 - the (1,0) block 1533a04f6461SBarry Smith - A11 - the (1,1) block 153430ad9308SMatthew Knepley 153530ad9308SMatthew Knepley Level: advanced 153630ad9308SMatthew Knepley 153730ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 153830ad9308SMatthew Knepley @*/ 1539a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 154030ad9308SMatthew Knepley { 154130ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 154230ad9308SMatthew Knepley 154330ad9308SMatthew Knepley PetscFunctionBegin; 15440700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1545c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1546a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1547a04f6461SBarry Smith if (A01) *A01 = jac->B; 1548a04f6461SBarry Smith if (A10) *A10 = jac->C; 1549a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 155030ad9308SMatthew Knepley PetscFunctionReturn(0); 155130ad9308SMatthew Knepley } 155230ad9308SMatthew Knepley 155379416396SBarry Smith EXTERN_C_BEGIN 155479416396SBarry Smith #undef __FUNCT__ 155579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 15567087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 155779416396SBarry Smith { 155879416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1559e69d4d44SBarry Smith PetscErrorCode ierr; 156079416396SBarry Smith 156179416396SBarry Smith PetscFunctionBegin; 156279416396SBarry Smith jac->type = type; 15633b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 15643b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 15653b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1566e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1567e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1568c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1569e69d4d44SBarry Smith 15703b224e63SBarry Smith } else { 15713b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 15723b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1573e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 15749e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1575c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 15763b224e63SBarry Smith } 157779416396SBarry Smith PetscFunctionReturn(0); 157879416396SBarry Smith } 157979416396SBarry Smith EXTERN_C_END 158079416396SBarry Smith 158151f519a2SBarry Smith EXTERN_C_BEGIN 158251f519a2SBarry Smith #undef __FUNCT__ 158351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 15847087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 158551f519a2SBarry Smith { 158651f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 158751f519a2SBarry Smith 158851f519a2SBarry Smith PetscFunctionBegin; 158965e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 159065e19b50SBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 159151f519a2SBarry Smith jac->bs = bs; 159251f519a2SBarry Smith PetscFunctionReturn(0); 159351f519a2SBarry Smith } 159451f519a2SBarry Smith EXTERN_C_END 159551f519a2SBarry Smith 159679416396SBarry Smith #undef __FUNCT__ 159779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1598bc08b0f1SBarry Smith /*@ 159979416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 160079416396SBarry Smith 160179416396SBarry Smith Collective on PC 160279416396SBarry Smith 160379416396SBarry Smith Input Parameter: 160479416396SBarry Smith . pc - the preconditioner context 160581540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 160679416396SBarry Smith 160779416396SBarry Smith Options Database Key: 1608a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 160979416396SBarry Smith 1610b02e2d75SMatthew G Knepley Level: Intermediate 161179416396SBarry Smith 161279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 161379416396SBarry Smith 161479416396SBarry Smith .seealso: PCCompositeSetType() 161579416396SBarry Smith 161679416396SBarry Smith @*/ 16177087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 161879416396SBarry Smith { 16194ac538c5SBarry Smith PetscErrorCode ierr; 162079416396SBarry Smith 162179416396SBarry Smith PetscFunctionBegin; 16220700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16234ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 162479416396SBarry Smith PetscFunctionReturn(0); 162579416396SBarry Smith } 162679416396SBarry Smith 1627b02e2d75SMatthew G Knepley #undef __FUNCT__ 1628b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1629b02e2d75SMatthew G Knepley /*@ 1630b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1631b02e2d75SMatthew G Knepley 1632b02e2d75SMatthew G Knepley Not collective 1633b02e2d75SMatthew G Knepley 1634b02e2d75SMatthew G Knepley Input Parameter: 1635b02e2d75SMatthew G Knepley . pc - the preconditioner context 1636b02e2d75SMatthew G Knepley 1637b02e2d75SMatthew G Knepley Output Parameter: 1638b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1639b02e2d75SMatthew G Knepley 1640b02e2d75SMatthew G Knepley Level: Intermediate 1641b02e2d75SMatthew G Knepley 1642b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1643b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1644b02e2d75SMatthew G Knepley @*/ 1645b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1646b02e2d75SMatthew G Knepley { 1647b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1648b02e2d75SMatthew G Knepley 1649b02e2d75SMatthew G Knepley PetscFunctionBegin; 1650b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1651b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1652b02e2d75SMatthew G Knepley *type = jac->type; 1653b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1654b02e2d75SMatthew G Knepley } 1655b02e2d75SMatthew G Knepley 16560971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 16570971522cSBarry Smith /*MC 1658a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1659a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 16600971522cSBarry Smith 1661edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1662edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 16630971522cSBarry Smith 1664a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 166569a612a9SBarry Smith and set the options directly on the resulting KSP object 16660971522cSBarry Smith 16670971522cSBarry Smith Level: intermediate 16680971522cSBarry Smith 166979416396SBarry Smith Options Database Keys: 167081540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 167181540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 167281540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 167381540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 16740f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 16750f188ba9SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag 1676435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 167779416396SBarry Smith 16785d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 16795d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 16805d4c12cdSJungho Lee 1681c8a0d604SMatthew G Knepley Notes: 1682c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1683d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1684d32f9abdSBarry Smith 1685d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1686d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1687d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1688d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1689d32f9abdSBarry Smith 1690c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1691c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1692c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1693c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1694c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1695a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1696c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1697c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1698c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1699a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1700c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1701c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 17025668aaf4SBarry Smith diag gives 1703c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1704c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 17055668aaf4SBarry Smith note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of 1706c8a0d604SMatthew G Knepley $ ( A00 0 ) 1707c8a0d604SMatthew G Knepley $ ( A10 S ) 1708c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1709c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1710c8a0d604SMatthew G Knepley $ ( 0 S ) 1711c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1712e69d4d44SBarry Smith 1713edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1714edf189efSBarry Smith is used automatically for a second block. 1715edf189efSBarry Smith 1716ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1717ff218e97SBarry Smith Generally it should be used with the AIJ format. 1718ff218e97SBarry Smith 1719ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1720ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1721ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 17220716a85fSBarry Smith 1723a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 17240971522cSBarry Smith 17257e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1726e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 17270971522cSBarry Smith M*/ 17280971522cSBarry Smith 17290971522cSBarry Smith EXTERN_C_BEGIN 17300971522cSBarry Smith #undef __FUNCT__ 17310971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 17327087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 17330971522cSBarry Smith { 17340971522cSBarry Smith PetscErrorCode ierr; 17350971522cSBarry Smith PC_FieldSplit *jac; 17360971522cSBarry Smith 17370971522cSBarry Smith PetscFunctionBegin; 173838f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 17390971522cSBarry Smith jac->bs = -1; 17400971522cSBarry Smith jac->nsplits = 0; 17413e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1742e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1743c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 174451f519a2SBarry Smith 17450971522cSBarry Smith pc->data = (void*)jac; 17460971522cSBarry Smith 17470971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1748421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 17490971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1750574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 17510971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 17520971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 17530971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 17540971522cSBarry Smith pc->ops->applyrichardson = 0; 17550971522cSBarry Smith 175669a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 175769a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 17580971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 17590971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1760704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1761704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 176279416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 176379416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 176451f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 176551f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 17660971522cSBarry Smith PetscFunctionReturn(0); 17670971522cSBarry Smith } 17680971522cSBarry Smith EXTERN_C_END 17690971522cSBarry Smith 17700971522cSBarry Smith 1771a541d17aSBarry Smith 1772