1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 33c48a1e8SJed Brown #include <petscdmcomposite.h> /*I "petscdmcomposite.h" I*/ 40971522cSBarry Smith 5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0}; 7c5d2311dSJed Brown 80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 90971522cSBarry Smith struct _PC_FieldSplitLink { 1069a612a9SBarry Smith KSP ksp; 110971522cSBarry Smith Vec x,y; 12db4c96c1SJed Brown char *splitname; 130971522cSBarry Smith PetscInt nfields; 145d4c12cdSJungho Lee PetscInt *fields,*fields_col; 151b9fc7fcSBarry Smith VecScatter sctx; 165d4c12cdSJungho Lee IS is,is_col; 1751f519a2SBarry Smith PC_FieldSplitLink next,previous; 180971522cSBarry Smith }; 190971522cSBarry Smith 200971522cSBarry Smith typedef struct { 2168dd23aaSBarry Smith PCCompositeType type; 22ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 23ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 24ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 2530ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2630ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2779416396SBarry Smith Vec *x,*y,w1,w2; 28519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 29519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3030ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 31ace3abfcSBarry Smith PetscBool issetup; 3230ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3330ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3430ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 35a04f6461SBarry Smith Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01 */ 36084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 37084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 38c9c6ffaaSJed Brown PCFieldSplitSchurFactType schurfactorization; 3930ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 4097bbdb24SBarry Smith PC_FieldSplitLink head; 4163ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 42c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 430971522cSBarry Smith } PC_FieldSplit; 440971522cSBarry Smith 4516913363SBarry Smith /* 4616913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4716913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4816913363SBarry Smith PC you could change this. 4916913363SBarry Smith */ 50084e4875SJed Brown 51e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 52084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 53084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 54084e4875SJed Brown { 55084e4875SJed Brown switch (jac->schurpre) { 56084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 57084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 58084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 59084e4875SJed Brown default: 60084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 61084e4875SJed Brown } 62084e4875SJed Brown } 63084e4875SJed Brown 64084e4875SJed Brown 650971522cSBarry Smith #undef __FUNCT__ 660971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 670971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 680971522cSBarry Smith { 690971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 700971522cSBarry Smith PetscErrorCode ierr; 71ace3abfcSBarry Smith PetscBool iascii; 720971522cSBarry Smith PetscInt i,j; 735a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 740971522cSBarry Smith 750971522cSBarry Smith PetscFunctionBegin; 762692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 770971522cSBarry Smith if (iascii) { 7851f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 7969a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 800971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 810971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 821ab39975SBarry Smith if (ilink->fields) { 830971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 8479416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 855a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 8679416396SBarry Smith if (j > 0) { 8779416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 8879416396SBarry Smith } 895a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 900971522cSBarry Smith } 910971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9279416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 931ab39975SBarry Smith } else { 941ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 951ab39975SBarry Smith } 965a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 975a9f2f41SSatish Balay ilink = ilink->next; 980971522cSBarry Smith } 990971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1000971522cSBarry Smith } else { 10165e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1020971522cSBarry Smith } 1030971522cSBarry Smith PetscFunctionReturn(0); 1040971522cSBarry Smith } 1050971522cSBarry Smith 1060971522cSBarry Smith #undef __FUNCT__ 1073b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1083b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1093b224e63SBarry Smith { 1103b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1113b224e63SBarry Smith PetscErrorCode ierr; 112ace3abfcSBarry Smith PetscBool iascii; 1133b224e63SBarry Smith PetscInt i,j; 1143b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1153b224e63SBarry Smith KSP ksp; 1163b224e63SBarry Smith 1173b224e63SBarry Smith PetscFunctionBegin; 1182692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1193b224e63SBarry Smith if (iascii) { 120c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1213b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1223b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1233b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1243b224e63SBarry Smith if (ilink->fields) { 1253b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1263b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1273b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1283b224e63SBarry Smith if (j > 0) { 1293b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1303b224e63SBarry Smith } 1313b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1323b224e63SBarry Smith } 1333b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1343b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1353b224e63SBarry Smith } else { 1363b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1373b224e63SBarry Smith } 1383b224e63SBarry Smith ilink = ilink->next; 1393b224e63SBarry Smith } 140435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1413b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14212cae6f2SJed Brown if (jac->schur) { 14312cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1443b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 14512cae6f2SJed Brown } else { 14612cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 14712cae6f2SJed Brown } 1483b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 149435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1503b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15112cae6f2SJed Brown if (jac->kspschur) { 1523b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 15312cae6f2SJed Brown } else { 15412cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 15512cae6f2SJed Brown } 1563b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1573b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1583b224e63SBarry Smith } else { 15965e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1603b224e63SBarry Smith } 1613b224e63SBarry Smith PetscFunctionReturn(0); 1623b224e63SBarry Smith } 1633b224e63SBarry Smith 1643b224e63SBarry Smith #undef __FUNCT__ 1656c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1666c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1676c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1686c924f48SJed Brown { 1696c924f48SJed Brown PetscErrorCode ierr; 1706c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1715d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 1725d4c12cdSJungho Lee PetscBool flg,flg_col; 1735d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 1746c924f48SJed Brown 1756c924f48SJed Brown PetscFunctionBegin; 1766c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 1775d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 1786c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 1796c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 1806c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 1815d4c12cdSJungho Lee ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 1826c924f48SJed Brown nfields = jac->bs; 18329499fbbSJungho Lee nfields_col = jac->bs; 1846c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 1855d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 1866c924f48SJed Brown if (!flg) break; 1875d4c12cdSJungho Lee else if (flg && !flg_col) { 1885d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1895d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 1905d4c12cdSJungho Lee } 1915d4c12cdSJungho Lee else { 1925d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1935d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 1945d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 1955d4c12cdSJungho Lee } 1966c924f48SJed Brown } 1976c924f48SJed Brown if (i > 0) { 1986c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 1996c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2006c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2016c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2026c924f48SJed Brown } 2036c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2045d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2056c924f48SJed Brown PetscFunctionReturn(0); 2066c924f48SJed Brown } 2076c924f48SJed Brown 2086c924f48SJed Brown #undef __FUNCT__ 20969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 21069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2110971522cSBarry Smith { 2120971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2130971522cSBarry Smith PetscErrorCode ierr; 2145a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2156ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2166c924f48SJed Brown PetscInt i; 2170971522cSBarry Smith 2180971522cSBarry Smith PetscFunctionBegin; 219d32f9abdSBarry Smith if (!ilink) { 220d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 221d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2227b62db95SJungho Lee PetscInt numFields, f; 2230784a22cSJed Brown char **fieldNames; 2247b62db95SJungho Lee IS *fields; 225*e7c4fc90SDmitry Karpeev DM *dms; 226*e7c4fc90SDmitry Karpeev ierr = DMCreateDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms); CHKERRQ(ierr); 2277b62db95SJungho Lee for(f = 0; f < numFields; ++f) { 2287b62db95SJungho Lee ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 2297b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 2307b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 2317b62db95SJungho Lee } 2327b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 2337b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 234*e7c4fc90SDmitry Karpeev if(dms) { 2358b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2367b62db95SJungho Lee for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) { 2377b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 2387b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 239*e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 2402fa5ba8aSJed Brown } 2417b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 2428b8307b2SJed Brown } 24366ffff09SJed Brown } else { 244521d7252SBarry Smith if (jac->bs <= 0) { 245704ba839SBarry Smith if (pc->pmat) { 246521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 247704ba839SBarry Smith } else { 248704ba839SBarry Smith jac->bs = 1; 249704ba839SBarry Smith } 250521d7252SBarry Smith } 251d32f9abdSBarry Smith 252acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2536ce1633cSBarry Smith if (stokes) { 2546ce1633cSBarry Smith IS zerodiags,rest; 2556ce1633cSBarry Smith PetscInt nmin,nmax; 2566ce1633cSBarry Smith 2576ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2586ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2596ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2606ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2616ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 262fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 263fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2646ce1633cSBarry Smith } else { 265d32f9abdSBarry Smith if (!flg) { 266d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 267d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2686c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2696c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 270d32f9abdSBarry Smith } 2716c924f48SJed Brown if (flg || !jac->splitdefined) { 272d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 273db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2746c924f48SJed Brown char splitname[8]; 2756c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2765d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 27779416396SBarry Smith } 2785d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 279521d7252SBarry Smith } 28066ffff09SJed Brown } 2816ce1633cSBarry Smith } 282edf189efSBarry Smith } else if (jac->nsplits == 1) { 283edf189efSBarry Smith if (ilink->is) { 284edf189efSBarry Smith IS is2; 285edf189efSBarry Smith PetscInt nmin,nmax; 286edf189efSBarry Smith 287edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 288edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 289db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 290fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 2917b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 29263ec74ffSBarry Smith } else if (jac->reset) { 29363ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 294d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 295d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 296d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 297d0af7cd3SBarry Smith if (pc->dm && !stokes) { 298d0af7cd3SBarry Smith PetscBool dmcomposite; 299d0af7cd3SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 300d0af7cd3SBarry Smith if (dmcomposite) { 301d0af7cd3SBarry Smith PetscInt nDM; 302d0af7cd3SBarry Smith IS *fields; 3037b62db95SJungho Lee DM *dms; 304d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 305d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 306d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 3077b62db95SJungho Lee ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr); 3087b62db95SJungho Lee ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr); 309d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 3107b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr); 3117b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr); 312d0af7cd3SBarry Smith ilink->is = fields[i]; 313d0af7cd3SBarry Smith ilink = ilink->next; 314edf189efSBarry Smith } 315d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 3167b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 317d0af7cd3SBarry Smith } 318d0af7cd3SBarry Smith } else { 319d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 320d0af7cd3SBarry Smith if (stokes) { 321d0af7cd3SBarry Smith IS zerodiags,rest; 322d0af7cd3SBarry Smith PetscInt nmin,nmax; 323d0af7cd3SBarry Smith 324d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 325d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 326d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 32720b26d62SBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 32820b26d62SBarry Smith ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr); 329d0af7cd3SBarry Smith ilink->is = rest; 330d0af7cd3SBarry Smith ilink->next->is = zerodiags; 33120b26d62SBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 332d0af7cd3SBarry Smith } 333d0af7cd3SBarry Smith } 334d0af7cd3SBarry Smith 3357b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 33669a612a9SBarry Smith PetscFunctionReturn(0); 33769a612a9SBarry Smith } 33869a612a9SBarry Smith 33969a612a9SBarry Smith #undef __FUNCT__ 34069a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 34169a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 34269a612a9SBarry Smith { 34369a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 34469a612a9SBarry Smith PetscErrorCode ierr; 3455a9f2f41SSatish Balay PC_FieldSplitLink ilink; 346cf502942SBarry Smith PetscInt i,nsplit,ccsize; 34769a612a9SBarry Smith MatStructure flag = pc->flag; 3485f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 34969a612a9SBarry Smith 35069a612a9SBarry Smith PetscFunctionBegin; 35169a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 35297bbdb24SBarry Smith nsplit = jac->nsplits; 3535a9f2f41SSatish Balay ilink = jac->head; 35497bbdb24SBarry Smith 35597bbdb24SBarry Smith /* get the matrices for each split */ 356704ba839SBarry Smith if (!jac->issetup) { 3571b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 35897bbdb24SBarry Smith 359704ba839SBarry Smith jac->issetup = PETSC_TRUE; 360704ba839SBarry Smith 3615d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 36251f519a2SBarry Smith bs = jac->bs; 36397bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 364cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 3651b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3661b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3671b9fc7fcSBarry Smith if (jac->defaultsplit) { 368704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 3695f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 370704ba839SBarry Smith } else if (!ilink->is) { 371ccb205f8SBarry Smith if (ilink->nfields > 1) { 3725f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 3735a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3745f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 3751b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3761b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3771b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 3785f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 37997bbdb24SBarry Smith } 38097bbdb24SBarry Smith } 381d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 3825f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 383ccb205f8SBarry Smith } else { 384704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 3855f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 386ccb205f8SBarry Smith } 3873e197d65SBarry Smith } 388704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 3895f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 3905f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 391704ba839SBarry Smith ilink = ilink->next; 3921b9fc7fcSBarry Smith } 3931b9fc7fcSBarry Smith } 3941b9fc7fcSBarry Smith 395704ba839SBarry Smith ilink = jac->head; 39697bbdb24SBarry Smith if (!jac->pmat) { 397cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 398cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3995f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 400704ba839SBarry Smith ilink = ilink->next; 401cf502942SBarry Smith } 40297bbdb24SBarry Smith } else { 403cf502942SBarry Smith for (i=0; i<nsplit; i++) { 4045f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 405704ba839SBarry Smith ilink = ilink->next; 406cf502942SBarry Smith } 40797bbdb24SBarry Smith } 408519d70e2SJed Brown if (jac->realdiagonal) { 409519d70e2SJed Brown ilink = jac->head; 410519d70e2SJed Brown if (!jac->mat) { 411519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 412519d70e2SJed Brown for (i=0; i<nsplit; i++) { 4135f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 414519d70e2SJed Brown ilink = ilink->next; 415519d70e2SJed Brown } 416519d70e2SJed Brown } else { 417519d70e2SJed Brown for (i=0; i<nsplit; i++) { 4185f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 419519d70e2SJed Brown ilink = ilink->next; 420519d70e2SJed Brown } 421519d70e2SJed Brown } 422519d70e2SJed Brown } else { 423519d70e2SJed Brown jac->mat = jac->pmat; 424519d70e2SJed Brown } 42597bbdb24SBarry Smith 4266c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 42768dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 42868dd23aaSBarry Smith ilink = jac->head; 42968dd23aaSBarry Smith if (!jac->Afield) { 43068dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 43168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4324aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 43368dd23aaSBarry Smith ilink = ilink->next; 43468dd23aaSBarry Smith } 43568dd23aaSBarry Smith } else { 43668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4374aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 43868dd23aaSBarry Smith ilink = ilink->next; 43968dd23aaSBarry Smith } 44068dd23aaSBarry Smith } 44168dd23aaSBarry Smith } 44268dd23aaSBarry Smith 4433b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 4443b224e63SBarry Smith IS ccis; 4454aa3045dSJed Brown PetscInt rstart,rend; 446093c86ffSJed Brown char lscname[256]; 447093c86ffSJed Brown PetscObject LSC_L; 448e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 44968dd23aaSBarry Smith 450e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 451e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 452e6cab6aaSJed Brown 4533b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4543b224e63SBarry Smith if (jac->schur) { 4553b224e63SBarry Smith ilink = jac->head; 45649bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4574aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 458fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4593b224e63SBarry Smith ilink = ilink->next; 46049bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4614aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 462fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 463519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 464084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4653b224e63SBarry Smith 4663b224e63SBarry Smith } else { 4671cee3971SBarry Smith KSP ksp; 4686c924f48SJed Brown char schurprefix[256]; 4693b224e63SBarry Smith 470a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4713b224e63SBarry Smith ilink = jac->head; 47249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4734aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 474fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4753b224e63SBarry Smith ilink = ilink->next; 47649bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4774aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 478fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4797d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 480519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 481a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4821cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 483e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 484a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 485a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 48620b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 48720b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4887b62db95SJungho Lee ierr = MatSetUp(jac->schur);CHKERRQ(ierr); 4893b224e63SBarry Smith 4903b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4919005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4921cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 493084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 494084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 495084e4875SJed Brown PC pc; 496cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 497084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 498084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 499e69d4d44SBarry Smith } 5006c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 5016c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 5023b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 50320b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 50420b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 5053b224e63SBarry Smith 5063b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 5073b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 5083b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 5093b224e63SBarry Smith ilink = jac->head; 5103b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 5113b224e63SBarry Smith ilink = ilink->next; 5123b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 5133b224e63SBarry Smith } 514093c86ffSJed Brown 515093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 516093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 517093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 518093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 519093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 520093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 521093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 522093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 523093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 5243b224e63SBarry Smith } else { 52597bbdb24SBarry Smith /* set up the individual PCs */ 52697bbdb24SBarry Smith i = 0; 5275a9f2f41SSatish Balay ilink = jac->head; 5285a9f2f41SSatish Balay while (ilink) { 529519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 5303b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 531c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 5325a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 53397bbdb24SBarry Smith i++; 5345a9f2f41SSatish Balay ilink = ilink->next; 5350971522cSBarry Smith } 53697bbdb24SBarry Smith 53797bbdb24SBarry Smith /* create work vectors for each split */ 5381b9fc7fcSBarry Smith if (!jac->x) { 53997bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5405a9f2f41SSatish Balay ilink = jac->head; 54197bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 542906ed7ccSBarry Smith Vec *vl,*vr; 543906ed7ccSBarry Smith 544906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 545906ed7ccSBarry Smith ilink->x = *vr; 546906ed7ccSBarry Smith ilink->y = *vl; 547906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 548906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5495a9f2f41SSatish Balay jac->x[i] = ilink->x; 5505a9f2f41SSatish Balay jac->y[i] = ilink->y; 5515a9f2f41SSatish Balay ilink = ilink->next; 55297bbdb24SBarry Smith } 5533b224e63SBarry Smith } 5543b224e63SBarry Smith } 5553b224e63SBarry Smith 5563b224e63SBarry Smith 557d0f46423SBarry Smith if (!jac->head->sctx) { 5583b224e63SBarry Smith Vec xtmp; 5593b224e63SBarry Smith 56079416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5611b9fc7fcSBarry Smith 5625a9f2f41SSatish Balay ilink = jac->head; 5631b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5641b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 565704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5665a9f2f41SSatish Balay ilink = ilink->next; 5671b9fc7fcSBarry Smith } 5686bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5691b9fc7fcSBarry Smith } 570c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 5710971522cSBarry Smith PetscFunctionReturn(0); 5720971522cSBarry Smith } 5730971522cSBarry Smith 5745a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 575ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 576ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5775a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 578ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 579ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 58079416396SBarry Smith 5810971522cSBarry Smith #undef __FUNCT__ 5823b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5833b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5843b224e63SBarry Smith { 5853b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5863b224e63SBarry Smith PetscErrorCode ierr; 5873b224e63SBarry Smith KSP ksp; 5883b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5893b224e63SBarry Smith 5903b224e63SBarry Smith PetscFunctionBegin; 5913b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5923b224e63SBarry Smith 593c5d2311dSJed Brown switch (jac->schurfactorization) { 594c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 595a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 596c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 598c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 599c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 600c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 601c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 602c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 603c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 604c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 605c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 606c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 607c5d2311dSJed Brown break; 608c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 609a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 610c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 611c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 612c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 613c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 614c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 615c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 616c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 617c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 618c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 619c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 620c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 621c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 622c5d2311dSJed Brown break; 623c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 624a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 625c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 626c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 627c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 628c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 629c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 630c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 631c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 632c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 633c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 634c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 635c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 636c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 637c5d2311dSJed Brown break; 638c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 639a04f6461SBarry 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 */ 6403b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6413b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6423b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6433b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6443b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6453b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6463b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6473b224e63SBarry Smith 6483b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6493b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6503b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6513b224e63SBarry Smith 6523b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6533b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6543b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6553b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6563b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 657c5d2311dSJed Brown } 6583b224e63SBarry Smith PetscFunctionReturn(0); 6593b224e63SBarry Smith } 6603b224e63SBarry Smith 6613b224e63SBarry Smith #undef __FUNCT__ 6620971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6630971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6640971522cSBarry Smith { 6650971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6660971522cSBarry Smith PetscErrorCode ierr; 6675a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 668af93645bSJed Brown PetscInt cnt; 6690971522cSBarry Smith 6700971522cSBarry Smith PetscFunctionBegin; 67151f519a2SBarry Smith CHKMEMQ; 67251f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 67351f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 67451f519a2SBarry Smith 67579416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6761b9fc7fcSBarry Smith if (jac->defaultsplit) { 6770971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6785a9f2f41SSatish Balay while (ilink) { 6795a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6805a9f2f41SSatish Balay ilink = ilink->next; 6810971522cSBarry Smith } 6820971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6831b9fc7fcSBarry Smith } else { 684efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6855a9f2f41SSatish Balay while (ilink) { 6865a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6875a9f2f41SSatish Balay ilink = ilink->next; 6881b9fc7fcSBarry Smith } 6891b9fc7fcSBarry Smith } 69016913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 69179416396SBarry Smith if (!jac->w1) { 69279416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 69379416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 69479416396SBarry Smith } 695efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6965a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6973e197d65SBarry Smith cnt = 1; 6985a9f2f41SSatish Balay while (ilink->next) { 6995a9f2f41SSatish Balay ilink = ilink->next; 7003e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 7013e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 7023e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 7033e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7043e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7053e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7063e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7073e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7083e197d65SBarry Smith } 70951f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 71011755939SBarry Smith cnt -= 2; 71151f519a2SBarry Smith while (ilink->previous) { 71251f519a2SBarry Smith ilink = ilink->previous; 71311755939SBarry Smith /* compute the residual only over the part of the vector needed */ 71411755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 71511755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 71611755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 71711755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 71811755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 71911755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 72011755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 72151f519a2SBarry Smith } 72211755939SBarry Smith } 72365e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 72451f519a2SBarry Smith CHKMEMQ; 7250971522cSBarry Smith PetscFunctionReturn(0); 7260971522cSBarry Smith } 7270971522cSBarry Smith 728421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 729ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 730ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 731421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 732ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 733ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 734421e10b8SBarry Smith 735421e10b8SBarry Smith #undef __FUNCT__ 7368c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 737421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 738421e10b8SBarry Smith { 739421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 740421e10b8SBarry Smith PetscErrorCode ierr; 741421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 742421e10b8SBarry Smith 743421e10b8SBarry Smith PetscFunctionBegin; 744421e10b8SBarry Smith CHKMEMQ; 745421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 746421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 747421e10b8SBarry Smith 748421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 749421e10b8SBarry Smith if (jac->defaultsplit) { 750421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 751421e10b8SBarry Smith while (ilink) { 752421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 753421e10b8SBarry Smith ilink = ilink->next; 754421e10b8SBarry Smith } 755421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 756421e10b8SBarry Smith } else { 757421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 758421e10b8SBarry Smith while (ilink) { 759421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 760421e10b8SBarry Smith ilink = ilink->next; 761421e10b8SBarry Smith } 762421e10b8SBarry Smith } 763421e10b8SBarry Smith } else { 764421e10b8SBarry Smith if (!jac->w1) { 765421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 766421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 767421e10b8SBarry Smith } 768421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 769421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 770421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 771421e10b8SBarry Smith while (ilink->next) { 772421e10b8SBarry Smith ilink = ilink->next; 7739989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 774421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 775421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 776421e10b8SBarry Smith } 777421e10b8SBarry Smith while (ilink->previous) { 778421e10b8SBarry Smith ilink = ilink->previous; 7799989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 780421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 781421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 782421e10b8SBarry Smith } 783421e10b8SBarry Smith } else { 784421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 785421e10b8SBarry Smith ilink = ilink->next; 786421e10b8SBarry Smith } 787421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 788421e10b8SBarry Smith while (ilink->previous) { 789421e10b8SBarry Smith ilink = ilink->previous; 7909989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 791421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 792421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 793421e10b8SBarry Smith } 794421e10b8SBarry Smith } 795421e10b8SBarry Smith } 796421e10b8SBarry Smith CHKMEMQ; 797421e10b8SBarry Smith PetscFunctionReturn(0); 798421e10b8SBarry Smith } 799421e10b8SBarry Smith 8000971522cSBarry Smith #undef __FUNCT__ 801574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 802574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 8030971522cSBarry Smith { 8040971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8050971522cSBarry Smith PetscErrorCode ierr; 8065a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 8070971522cSBarry Smith 8080971522cSBarry Smith PetscFunctionBegin; 8095a9f2f41SSatish Balay while (ilink) { 810574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 811fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 812fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 813fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 814fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 815b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 8165a9f2f41SSatish Balay next = ilink->next; 8175a9f2f41SSatish Balay ilink = next; 8180971522cSBarry Smith } 81905b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 820574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 821574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 822574deadeSBarry Smith } else if (jac->mat) { 823574deadeSBarry Smith jac->mat = PETSC_NULL; 824574deadeSBarry Smith } 82597bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 82668dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8276bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8286bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8296bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8306bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8316bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8326bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8336bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 83463ec74ffSBarry Smith jac->reset = PETSC_TRUE; 835574deadeSBarry Smith PetscFunctionReturn(0); 836574deadeSBarry Smith } 837574deadeSBarry Smith 838574deadeSBarry Smith #undef __FUNCT__ 839574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 840574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 841574deadeSBarry Smith { 842574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 843574deadeSBarry Smith PetscErrorCode ierr; 844574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 845574deadeSBarry Smith 846574deadeSBarry Smith PetscFunctionBegin; 847574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 848574deadeSBarry Smith while (ilink) { 8496bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 850574deadeSBarry Smith next = ilink->next; 851574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 852574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 8535d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 854574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 855574deadeSBarry Smith ilink = next; 856574deadeSBarry Smith } 857574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 858c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8597b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 8607b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 8617b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 8627b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 8637b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 864ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 865c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 8660971522cSBarry Smith PetscFunctionReturn(0); 8670971522cSBarry Smith } 8680971522cSBarry Smith 8690971522cSBarry Smith #undef __FUNCT__ 8700971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8710971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8720971522cSBarry Smith { 8731b9fc7fcSBarry Smith PetscErrorCode ierr; 8746c924f48SJed Brown PetscInt bs; 875bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8769dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8773b224e63SBarry Smith PCCompositeType ctype; 878*e7c4fc90SDmitry Karpeev DM ddm; 879*e7c4fc90SDmitry Karpeev char ddm_name[1025]; 8801b9fc7fcSBarry Smith 8810971522cSBarry Smith PetscFunctionBegin; 8821b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 883*e7c4fc90SDmitry Karpeev if(pc->dm) { 884*e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 885*e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 886*e7c4fc90SDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition_dm", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 887*e7c4fc90SDmitry Karpeev if(flg) { 888*e7c4fc90SDmitry Karpeev ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 889*e7c4fc90SDmitry Karpeev if(!ddm) { 890*e7c4fc90SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name); 891*e7c4fc90SDmitry Karpeev } 892*e7c4fc90SDmitry Karpeev ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr); 893*e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 894*e7c4fc90SDmitry Karpeev } 895*e7c4fc90SDmitry Karpeev } 896acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 89751f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 89851f519a2SBarry Smith if (flg) { 89951f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 90051f519a2SBarry Smith } 901704ba839SBarry Smith 902435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 903c0adfefeSBarry Smith if (stokes) { 904c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 905c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 906c0adfefeSBarry Smith } 907c0adfefeSBarry Smith 9083b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 9093b224e63SBarry Smith if (flg) { 9103b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 9113b224e63SBarry Smith } 912c30613efSMatthew Knepley /* Only setup fields once */ 913c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 914d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 915d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 9166c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 9176c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 918d32f9abdSBarry Smith } 919c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 920c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 921c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 922c9c6ffaaSJed 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); 923084e4875SJed 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); 924c5d2311dSJed Brown } 9251b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 9260971522cSBarry Smith PetscFunctionReturn(0); 9270971522cSBarry Smith } 9280971522cSBarry Smith 9290971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 9300971522cSBarry Smith 9310971522cSBarry Smith EXTERN_C_BEGIN 9320971522cSBarry Smith #undef __FUNCT__ 9330971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 9345d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 9350971522cSBarry Smith { 93697bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9370971522cSBarry Smith PetscErrorCode ierr; 9385a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 93969a612a9SBarry Smith char prefix[128]; 9405d4c12cdSJungho Lee PetscInt i; 9410971522cSBarry Smith 9420971522cSBarry Smith PetscFunctionBegin; 9436c924f48SJed Brown if (jac->splitdefined) { 9446c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9456c924f48SJed Brown PetscFunctionReturn(0); 9466c924f48SJed Brown } 94751f519a2SBarry Smith for (i=0; i<n; i++) { 948e32f2f54SBarry 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); 949e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 95051f519a2SBarry Smith } 951704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 952a04f6461SBarry Smith if (splitname) { 953db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 954a04f6461SBarry Smith } else { 955a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 956a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 957a04f6461SBarry Smith } 958704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 9595a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 9605d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 9615d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 9625a9f2f41SSatish Balay ilink->nfields = n; 9635a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9647adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9651cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9665a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9679005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 96869a612a9SBarry Smith 969a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9705a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9710971522cSBarry Smith 9720971522cSBarry Smith if (!next) { 9735a9f2f41SSatish Balay jac->head = ilink; 97451f519a2SBarry Smith ilink->previous = PETSC_NULL; 9750971522cSBarry Smith } else { 9760971522cSBarry Smith while (next->next) { 9770971522cSBarry Smith next = next->next; 9780971522cSBarry Smith } 9795a9f2f41SSatish Balay next->next = ilink; 98051f519a2SBarry Smith ilink->previous = next; 9810971522cSBarry Smith } 9820971522cSBarry Smith jac->nsplits++; 9830971522cSBarry Smith PetscFunctionReturn(0); 9840971522cSBarry Smith } 9850971522cSBarry Smith EXTERN_C_END 9860971522cSBarry Smith 987e69d4d44SBarry Smith EXTERN_C_BEGIN 988e69d4d44SBarry Smith #undef __FUNCT__ 989e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9907087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 991e69d4d44SBarry Smith { 992e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 993e69d4d44SBarry Smith PetscErrorCode ierr; 994e69d4d44SBarry Smith 995e69d4d44SBarry Smith PetscFunctionBegin; 996519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 997e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 998e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 99913e0d083SBarry Smith if (n) *n = jac->nsplits; 1000e69d4d44SBarry Smith PetscFunctionReturn(0); 1001e69d4d44SBarry Smith } 1002e69d4d44SBarry Smith EXTERN_C_END 10030971522cSBarry Smith 10040971522cSBarry Smith EXTERN_C_BEGIN 10050971522cSBarry Smith #undef __FUNCT__ 100669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 10077087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 10080971522cSBarry Smith { 10090971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10100971522cSBarry Smith PetscErrorCode ierr; 10110971522cSBarry Smith PetscInt cnt = 0; 10125a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 10130971522cSBarry Smith 10140971522cSBarry Smith PetscFunctionBegin; 10155d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 10165a9f2f41SSatish Balay while (ilink) { 10175a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 10185a9f2f41SSatish Balay ilink = ilink->next; 10190971522cSBarry Smith } 10205d480477SMatthew 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); 102113e0d083SBarry Smith if (n) *n = jac->nsplits; 10220971522cSBarry Smith PetscFunctionReturn(0); 10230971522cSBarry Smith } 10240971522cSBarry Smith EXTERN_C_END 10250971522cSBarry Smith 1026704ba839SBarry Smith EXTERN_C_BEGIN 1027704ba839SBarry Smith #undef __FUNCT__ 1028704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 10297087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1030704ba839SBarry Smith { 1031704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1032704ba839SBarry Smith PetscErrorCode ierr; 1033704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1034704ba839SBarry Smith char prefix[128]; 1035704ba839SBarry Smith 1036704ba839SBarry Smith PetscFunctionBegin; 10376c924f48SJed Brown if (jac->splitdefined) { 10386c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10396c924f48SJed Brown PetscFunctionReturn(0); 10406c924f48SJed Brown } 104116913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1042a04f6461SBarry Smith if (splitname) { 1043db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1044a04f6461SBarry Smith } else { 1045b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1046b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1047a04f6461SBarry Smith } 10481ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1049b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1050b5787286SJed Brown ilink->is = is; 1051b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1052b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1053b5787286SJed Brown ilink->is_col = is; 1054704ba839SBarry Smith ilink->next = PETSC_NULL; 1055704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10561cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1057704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10589005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1059704ba839SBarry Smith 1060a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1061704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1062704ba839SBarry Smith 1063704ba839SBarry Smith if (!next) { 1064704ba839SBarry Smith jac->head = ilink; 1065704ba839SBarry Smith ilink->previous = PETSC_NULL; 1066704ba839SBarry Smith } else { 1067704ba839SBarry Smith while (next->next) { 1068704ba839SBarry Smith next = next->next; 1069704ba839SBarry Smith } 1070704ba839SBarry Smith next->next = ilink; 1071704ba839SBarry Smith ilink->previous = next; 1072704ba839SBarry Smith } 1073704ba839SBarry Smith jac->nsplits++; 1074704ba839SBarry Smith 1075704ba839SBarry Smith PetscFunctionReturn(0); 1076704ba839SBarry Smith } 1077704ba839SBarry Smith EXTERN_C_END 1078704ba839SBarry Smith 10790971522cSBarry Smith #undef __FUNCT__ 10800971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10810971522cSBarry Smith /*@ 10820971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10830971522cSBarry Smith 1084ad4df100SBarry Smith Logically Collective on PC 10850971522cSBarry Smith 10860971522cSBarry Smith Input Parameters: 10870971522cSBarry Smith + pc - the preconditioner context 1088a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10890971522cSBarry Smith . n - the number of fields in this split 1090db4c96c1SJed Brown - fields - the fields in this split 10910971522cSBarry Smith 10920971522cSBarry Smith Level: intermediate 10930971522cSBarry Smith 1094d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1095d32f9abdSBarry Smith 1096d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1097d32f9abdSBarry 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 1098d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1099d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1100d32f9abdSBarry Smith 1101db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1102db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1103db4c96c1SJed Brown 11045d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 11055d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 11065d4c12cdSJungho Lee available when this routine is called. 11075d4c12cdSJungho Lee 1108d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 11090971522cSBarry Smith 11100971522cSBarry Smith @*/ 11115d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11120971522cSBarry Smith { 11134ac538c5SBarry Smith PetscErrorCode ierr; 11140971522cSBarry Smith 11150971522cSBarry Smith PetscFunctionBegin; 11160700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1117db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1118db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1119db4c96c1SJed Brown PetscValidIntPointer(fields,3); 11205d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 11210971522cSBarry Smith PetscFunctionReturn(0); 11220971522cSBarry Smith } 11230971522cSBarry Smith 11240971522cSBarry Smith #undef __FUNCT__ 1125704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1126704ba839SBarry Smith /*@ 1127704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1128704ba839SBarry Smith 1129ad4df100SBarry Smith Logically Collective on PC 1130704ba839SBarry Smith 1131704ba839SBarry Smith Input Parameters: 1132704ba839SBarry Smith + pc - the preconditioner context 1133a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1134db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1135704ba839SBarry Smith 1136d32f9abdSBarry Smith 1137a6ffb8dbSJed Brown Notes: 1138a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1139a6ffb8dbSJed Brown 1140db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1141db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1142d32f9abdSBarry Smith 1143704ba839SBarry Smith Level: intermediate 1144704ba839SBarry Smith 1145704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1146704ba839SBarry Smith 1147704ba839SBarry Smith @*/ 11487087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1149704ba839SBarry Smith { 11504ac538c5SBarry Smith PetscErrorCode ierr; 1151704ba839SBarry Smith 1152704ba839SBarry Smith PetscFunctionBegin; 11530700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11547b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1155db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11564ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1157704ba839SBarry Smith PetscFunctionReturn(0); 1158704ba839SBarry Smith } 1159704ba839SBarry Smith 1160704ba839SBarry Smith #undef __FUNCT__ 116157a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 116257a9adfeSMatthew G Knepley /*@ 116357a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 116457a9adfeSMatthew G Knepley 116557a9adfeSMatthew G Knepley Logically Collective on PC 116657a9adfeSMatthew G Knepley 116757a9adfeSMatthew G Knepley Input Parameters: 116857a9adfeSMatthew G Knepley + pc - the preconditioner context 116957a9adfeSMatthew G Knepley - splitname - name of this split 117057a9adfeSMatthew G Knepley 117157a9adfeSMatthew G Knepley Output Parameter: 117257a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 117357a9adfeSMatthew G Knepley 117457a9adfeSMatthew G Knepley Level: intermediate 117557a9adfeSMatthew G Knepley 117657a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 117757a9adfeSMatthew G Knepley 117857a9adfeSMatthew G Knepley @*/ 117957a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 118057a9adfeSMatthew G Knepley { 118157a9adfeSMatthew G Knepley PetscErrorCode ierr; 118257a9adfeSMatthew G Knepley 118357a9adfeSMatthew G Knepley PetscFunctionBegin; 118457a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 118557a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 118657a9adfeSMatthew G Knepley PetscValidPointer(is,3); 118757a9adfeSMatthew G Knepley { 118857a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 118957a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 119057a9adfeSMatthew G Knepley PetscBool found; 119157a9adfeSMatthew G Knepley 119257a9adfeSMatthew G Knepley *is = PETSC_NULL; 119357a9adfeSMatthew G Knepley while(ilink) { 119457a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 119557a9adfeSMatthew G Knepley if (found) { 119657a9adfeSMatthew G Knepley *is = ilink->is; 119757a9adfeSMatthew G Knepley break; 119857a9adfeSMatthew G Knepley } 119957a9adfeSMatthew G Knepley ilink = ilink->next; 120057a9adfeSMatthew G Knepley } 120157a9adfeSMatthew G Knepley } 120257a9adfeSMatthew G Knepley PetscFunctionReturn(0); 120357a9adfeSMatthew G Knepley } 120457a9adfeSMatthew G Knepley 120557a9adfeSMatthew G Knepley #undef __FUNCT__ 120651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 120751f519a2SBarry Smith /*@ 120851f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 120951f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 121051f519a2SBarry Smith 1211ad4df100SBarry Smith Logically Collective on PC 121251f519a2SBarry Smith 121351f519a2SBarry Smith Input Parameters: 121451f519a2SBarry Smith + pc - the preconditioner context 121551f519a2SBarry Smith - bs - the block size 121651f519a2SBarry Smith 121751f519a2SBarry Smith Level: intermediate 121851f519a2SBarry Smith 121951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 122051f519a2SBarry Smith 122151f519a2SBarry Smith @*/ 12227087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 122351f519a2SBarry Smith { 12244ac538c5SBarry Smith PetscErrorCode ierr; 122551f519a2SBarry Smith 122651f519a2SBarry Smith PetscFunctionBegin; 12270700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1228c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 12294ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 123051f519a2SBarry Smith PetscFunctionReturn(0); 123151f519a2SBarry Smith } 123251f519a2SBarry Smith 123351f519a2SBarry Smith #undef __FUNCT__ 123469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 12350971522cSBarry Smith /*@C 123669a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 12370971522cSBarry Smith 123869a612a9SBarry Smith Collective on KSP 12390971522cSBarry Smith 12400971522cSBarry Smith Input Parameter: 12410971522cSBarry Smith . pc - the preconditioner context 12420971522cSBarry Smith 12430971522cSBarry Smith Output Parameters: 124413e0d083SBarry Smith + n - the number of splits 124569a612a9SBarry Smith - pc - the array of KSP contexts 12460971522cSBarry Smith 12470971522cSBarry Smith Note: 1248d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1249d32f9abdSBarry Smith (not the KSP just the array that contains them). 12500971522cSBarry Smith 125169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12520971522cSBarry Smith 12530971522cSBarry Smith Level: advanced 12540971522cSBarry Smith 12550971522cSBarry Smith .seealso: PCFIELDSPLIT 12560971522cSBarry Smith @*/ 12577087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12580971522cSBarry Smith { 12594ac538c5SBarry Smith PetscErrorCode ierr; 12600971522cSBarry Smith 12610971522cSBarry Smith PetscFunctionBegin; 12620700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 126313e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 12644ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 12650971522cSBarry Smith PetscFunctionReturn(0); 12660971522cSBarry Smith } 12670971522cSBarry Smith 1268e69d4d44SBarry Smith #undef __FUNCT__ 1269e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1270e69d4d44SBarry Smith /*@ 1271e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1272a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1273e69d4d44SBarry Smith 1274e69d4d44SBarry Smith Collective on PC 1275e69d4d44SBarry Smith 1276e69d4d44SBarry Smith Input Parameters: 1277e69d4d44SBarry Smith + pc - the preconditioner context 1278fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1279084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1280084e4875SJed Brown 1281e69d4d44SBarry Smith Options Database: 1282084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1283e69d4d44SBarry Smith 1284fd1303e9SJungho Lee Notes: 1285fd1303e9SJungho Lee $ If ptype is 1286fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1287fd1303e9SJungho Lee $ to this function). 1288fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1289fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1290fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1291fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1292fd1303e9SJungho Lee $ preconditioner 1293fd1303e9SJungho Lee 1294fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1295fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1296fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1297fd1303e9SJungho Lee 1298fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1299fd1303e9SJungho Lee the name since it sets a proceedure to use. 1300fd1303e9SJungho Lee 1301e69d4d44SBarry Smith Level: intermediate 1302e69d4d44SBarry Smith 1303fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1304e69d4d44SBarry Smith 1305e69d4d44SBarry Smith @*/ 13067087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1307e69d4d44SBarry Smith { 13084ac538c5SBarry Smith PetscErrorCode ierr; 1309e69d4d44SBarry Smith 1310e69d4d44SBarry Smith PetscFunctionBegin; 13110700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13124ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1313e69d4d44SBarry Smith PetscFunctionReturn(0); 1314e69d4d44SBarry Smith } 1315e69d4d44SBarry Smith 1316e69d4d44SBarry Smith EXTERN_C_BEGIN 1317e69d4d44SBarry Smith #undef __FUNCT__ 1318e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 13197087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1320e69d4d44SBarry Smith { 1321e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1322084e4875SJed Brown PetscErrorCode ierr; 1323e69d4d44SBarry Smith 1324e69d4d44SBarry Smith PetscFunctionBegin; 1325084e4875SJed Brown jac->schurpre = ptype; 1326084e4875SJed Brown if (pre) { 13276bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1328084e4875SJed Brown jac->schur_user = pre; 1329084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1330084e4875SJed Brown } 1331e69d4d44SBarry Smith PetscFunctionReturn(0); 1332e69d4d44SBarry Smith } 1333e69d4d44SBarry Smith EXTERN_C_END 1334e69d4d44SBarry Smith 133530ad9308SMatthew Knepley #undef __FUNCT__ 1336c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1337ab1df9f5SJed Brown /*@ 1338c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1339ab1df9f5SJed Brown 1340ab1df9f5SJed Brown Collective on PC 1341ab1df9f5SJed Brown 1342ab1df9f5SJed Brown Input Parameters: 1343ab1df9f5SJed Brown + pc - the preconditioner context 1344c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1345ab1df9f5SJed Brown 1346ab1df9f5SJed Brown Options Database: 1347c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1348ab1df9f5SJed Brown 1349ab1df9f5SJed Brown 1350ab1df9f5SJed Brown Level: intermediate 1351ab1df9f5SJed Brown 1352ab1df9f5SJed Brown Notes: 1353ab1df9f5SJed Brown The FULL factorization is 1354ab1df9f5SJed Brown 1355ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1356ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1357ab1df9f5SJed Brown 13586be4592eSBarry 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, 13596be4592eSBarry 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). 1360ab1df9f5SJed Brown 13616be4592eSBarry 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 13626be4592eSBarry 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 13636be4592eSBarry 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, 13646be4592eSBarry 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. 1365ab1df9f5SJed Brown 13666be4592eSBarry 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 13676be4592eSBarry 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). 1368ab1df9f5SJed Brown 1369ab1df9f5SJed Brown References: 1370ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1371ab1df9f5SJed Brown 1372ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1373ab1df9f5SJed Brown 1374ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1375ab1df9f5SJed Brown @*/ 1376c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1377ab1df9f5SJed Brown { 1378ab1df9f5SJed Brown PetscErrorCode ierr; 1379ab1df9f5SJed Brown 1380ab1df9f5SJed Brown PetscFunctionBegin; 1381ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1382c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1383ab1df9f5SJed Brown PetscFunctionReturn(0); 1384ab1df9f5SJed Brown } 1385ab1df9f5SJed Brown 1386ab1df9f5SJed Brown #undef __FUNCT__ 1387c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1388c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1389ab1df9f5SJed Brown { 1390ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1391ab1df9f5SJed Brown 1392ab1df9f5SJed Brown PetscFunctionBegin; 1393ab1df9f5SJed Brown jac->schurfactorization = ftype; 1394ab1df9f5SJed Brown PetscFunctionReturn(0); 1395ab1df9f5SJed Brown } 1396ab1df9f5SJed Brown 1397ab1df9f5SJed Brown #undef __FUNCT__ 139830ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 139930ad9308SMatthew Knepley /*@C 14008c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 140130ad9308SMatthew Knepley 140230ad9308SMatthew Knepley Collective on KSP 140330ad9308SMatthew Knepley 140430ad9308SMatthew Knepley Input Parameter: 140530ad9308SMatthew Knepley . pc - the preconditioner context 140630ad9308SMatthew Knepley 140730ad9308SMatthew Knepley Output Parameters: 1408a04f6461SBarry Smith + A00 - the (0,0) block 1409a04f6461SBarry Smith . A01 - the (0,1) block 1410a04f6461SBarry Smith . A10 - the (1,0) block 1411a04f6461SBarry Smith - A11 - the (1,1) block 141230ad9308SMatthew Knepley 141330ad9308SMatthew Knepley Level: advanced 141430ad9308SMatthew Knepley 141530ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 141630ad9308SMatthew Knepley @*/ 1417a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 141830ad9308SMatthew Knepley { 141930ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 142030ad9308SMatthew Knepley 142130ad9308SMatthew Knepley PetscFunctionBegin; 14220700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1423c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1424a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1425a04f6461SBarry Smith if (A01) *A01 = jac->B; 1426a04f6461SBarry Smith if (A10) *A10 = jac->C; 1427a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 142830ad9308SMatthew Knepley PetscFunctionReturn(0); 142930ad9308SMatthew Knepley } 143030ad9308SMatthew Knepley 143179416396SBarry Smith EXTERN_C_BEGIN 143279416396SBarry Smith #undef __FUNCT__ 143379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 14347087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 143579416396SBarry Smith { 143679416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1437e69d4d44SBarry Smith PetscErrorCode ierr; 143879416396SBarry Smith 143979416396SBarry Smith PetscFunctionBegin; 144079416396SBarry Smith jac->type = type; 14413b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 14423b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 14433b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1444e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1445e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1446c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1447e69d4d44SBarry Smith 14483b224e63SBarry Smith } else { 14493b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 14503b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1451e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14529e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1453c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 14543b224e63SBarry Smith } 145579416396SBarry Smith PetscFunctionReturn(0); 145679416396SBarry Smith } 145779416396SBarry Smith EXTERN_C_END 145879416396SBarry Smith 145951f519a2SBarry Smith EXTERN_C_BEGIN 146051f519a2SBarry Smith #undef __FUNCT__ 146151f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 14627087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 146351f519a2SBarry Smith { 146451f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 146551f519a2SBarry Smith 146651f519a2SBarry Smith PetscFunctionBegin; 146765e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 146865e19b50SBarry 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); 146951f519a2SBarry Smith jac->bs = bs; 147051f519a2SBarry Smith PetscFunctionReturn(0); 147151f519a2SBarry Smith } 147251f519a2SBarry Smith EXTERN_C_END 147351f519a2SBarry Smith 147479416396SBarry Smith #undef __FUNCT__ 147579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1476bc08b0f1SBarry Smith /*@ 147779416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 147879416396SBarry Smith 147979416396SBarry Smith Collective on PC 148079416396SBarry Smith 148179416396SBarry Smith Input Parameter: 148279416396SBarry Smith . pc - the preconditioner context 148381540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 148479416396SBarry Smith 148579416396SBarry Smith Options Database Key: 1486a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 148779416396SBarry Smith 148879416396SBarry Smith Level: Developer 148979416396SBarry Smith 149079416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 149179416396SBarry Smith 149279416396SBarry Smith .seealso: PCCompositeSetType() 149379416396SBarry Smith 149479416396SBarry Smith @*/ 14957087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 149679416396SBarry Smith { 14974ac538c5SBarry Smith PetscErrorCode ierr; 149879416396SBarry Smith 149979416396SBarry Smith PetscFunctionBegin; 15000700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15014ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 150279416396SBarry Smith PetscFunctionReturn(0); 150379416396SBarry Smith } 150479416396SBarry Smith 15050971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 15060971522cSBarry Smith /*MC 1507a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1508a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 15090971522cSBarry Smith 1510edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1511edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 15120971522cSBarry Smith 1513a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 151469a612a9SBarry Smith and set the options directly on the resulting KSP object 15150971522cSBarry Smith 15160971522cSBarry Smith Level: intermediate 15170971522cSBarry Smith 151879416396SBarry Smith Options Database Keys: 151981540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 152081540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 152181540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 152281540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 152381540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1524e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1525435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 152679416396SBarry Smith 15275d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 15285d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 15295d4c12cdSJungho Lee 1530c8a0d604SMatthew G Knepley Notes: 1531c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1532d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1533d32f9abdSBarry Smith 1534d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1535d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1536d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1537d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1538d32f9abdSBarry Smith 1539c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1540c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1541c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1542c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1543c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1544a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1545c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1546c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1547c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1548a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1549c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1550c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 15515668aaf4SBarry Smith diag gives 1552c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1553c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 15545668aaf4SBarry 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 1555c8a0d604SMatthew G Knepley $ ( A00 0 ) 1556c8a0d604SMatthew G Knepley $ ( A10 S ) 1557c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1558c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1559c8a0d604SMatthew G Knepley $ ( 0 S ) 1560c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1561e69d4d44SBarry Smith 1562edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1563edf189efSBarry Smith is used automatically for a second block. 1564edf189efSBarry Smith 1565ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1566ff218e97SBarry Smith Generally it should be used with the AIJ format. 1567ff218e97SBarry Smith 1568ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1569ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1570ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 15710716a85fSBarry Smith 1572a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 15730971522cSBarry Smith 15747e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1575e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 15760971522cSBarry Smith M*/ 15770971522cSBarry Smith 15780971522cSBarry Smith EXTERN_C_BEGIN 15790971522cSBarry Smith #undef __FUNCT__ 15800971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 15817087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 15820971522cSBarry Smith { 15830971522cSBarry Smith PetscErrorCode ierr; 15840971522cSBarry Smith PC_FieldSplit *jac; 15850971522cSBarry Smith 15860971522cSBarry Smith PetscFunctionBegin; 158738f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 15880971522cSBarry Smith jac->bs = -1; 15890971522cSBarry Smith jac->nsplits = 0; 15903e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1591e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1592c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 159351f519a2SBarry Smith 15940971522cSBarry Smith pc->data = (void*)jac; 15950971522cSBarry Smith 15960971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1597421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 15980971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1599574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 16000971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 16010971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 16020971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 16030971522cSBarry Smith pc->ops->applyrichardson = 0; 16040971522cSBarry Smith 160569a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 160669a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 16070971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 16080971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1609704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1610704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 161179416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 161279416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 161351f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 161451f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 16150971522cSBarry Smith PetscFunctionReturn(0); 16160971522cSBarry Smith } 16170971522cSBarry Smith EXTERN_C_END 16180971522cSBarry Smith 16190971522cSBarry Smith 1620a541d17aSBarry Smith 1621