1dba47a55SKris Buschelman 2c6db04a5SJed Brown #include <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}; 6f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0}; 7c5d2311dSJed Brown 8c5d2311dSJed Brown typedef enum { 9c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG, 10c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER, 11c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER, 12c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL 13c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType; 14084e4875SJed Brown 150971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 160971522cSBarry Smith struct _PC_FieldSplitLink { 1769a612a9SBarry Smith KSP ksp; 180971522cSBarry Smith Vec x,y; 19db4c96c1SJed Brown char *splitname; 200971522cSBarry Smith PetscInt nfields; 215d4c12cdSJungho Lee PetscInt *fields,*fields_col; 221b9fc7fcSBarry Smith VecScatter sctx; 235d4c12cdSJungho Lee IS is,is_col; 2451f519a2SBarry Smith PC_FieldSplitLink next,previous; 250971522cSBarry Smith }; 260971522cSBarry Smith 270971522cSBarry Smith typedef struct { 2868dd23aaSBarry Smith PCCompositeType type; 29ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 30ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 31ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 3230ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 3330ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 3479416396SBarry Smith Vec *x,*y,w1,w2; 35519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 36519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3730ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 38ace3abfcSBarry Smith PetscBool issetup; 3930ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 4030ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 4130ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 42a04f6461SBarry Smith Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01 */ 43084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 44084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 45c5d2311dSJed Brown PCFieldSplitSchurFactorizationType schurfactorization; 4630ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 4797bbdb24SBarry Smith PC_FieldSplitLink head; 4863ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 49c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 500971522cSBarry Smith } PC_FieldSplit; 510971522cSBarry Smith 5216913363SBarry Smith /* 5316913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5416913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5516913363SBarry Smith PC you could change this. 5616913363SBarry Smith */ 57084e4875SJed Brown 58e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 59084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 60084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 61084e4875SJed Brown { 62084e4875SJed Brown switch (jac->schurpre) { 63084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 64084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 65084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 66084e4875SJed Brown default: 67084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 68084e4875SJed Brown } 69084e4875SJed Brown } 70084e4875SJed Brown 71084e4875SJed Brown 720971522cSBarry Smith #undef __FUNCT__ 730971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 740971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 750971522cSBarry Smith { 760971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 770971522cSBarry Smith PetscErrorCode ierr; 78ace3abfcSBarry Smith PetscBool iascii; 790971522cSBarry Smith PetscInt i,j; 805a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 810971522cSBarry Smith 820971522cSBarry Smith PetscFunctionBegin; 832692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 840971522cSBarry Smith if (iascii) { 8551f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 8669a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 870971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 880971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 891ab39975SBarry Smith if (ilink->fields) { 900971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 9179416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 925a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9379416396SBarry Smith if (j > 0) { 9479416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9579416396SBarry Smith } 965a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 970971522cSBarry Smith } 980971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1001ab39975SBarry Smith } else { 1011ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1021ab39975SBarry Smith } 1035a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1045a9f2f41SSatish Balay ilink = ilink->next; 1050971522cSBarry Smith } 1060971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1070971522cSBarry Smith } else { 10865e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1090971522cSBarry Smith } 1100971522cSBarry Smith PetscFunctionReturn(0); 1110971522cSBarry Smith } 1120971522cSBarry Smith 1130971522cSBarry Smith #undef __FUNCT__ 1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1163b224e63SBarry Smith { 1173b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1183b224e63SBarry Smith PetscErrorCode ierr; 119ace3abfcSBarry Smith PetscBool iascii; 1203b224e63SBarry Smith PetscInt i,j; 1213b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1223b224e63SBarry Smith KSP ksp; 1233b224e63SBarry Smith 1243b224e63SBarry Smith PetscFunctionBegin; 1252692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1263b224e63SBarry Smith if (iascii) { 127c5d2311dSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr); 1283b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1293b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1303b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1313b224e63SBarry Smith if (ilink->fields) { 1323b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1333b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1343b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1353b224e63SBarry Smith if (j > 0) { 1363b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1373b224e63SBarry Smith } 1383b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1393b224e63SBarry Smith } 1403b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1413b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1423b224e63SBarry Smith } else { 1433b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1443b224e63SBarry Smith } 1453b224e63SBarry Smith ilink = ilink->next; 1463b224e63SBarry Smith } 147435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1483b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14912cae6f2SJed Brown if (jac->schur) { 15012cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1513b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 15212cae6f2SJed Brown } else { 15312cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 15412cae6f2SJed Brown } 1553b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 156435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1573b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15812cae6f2SJed Brown if (jac->kspschur) { 1593b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 16012cae6f2SJed Brown } else { 16112cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 16212cae6f2SJed Brown } 1633b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1643b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1653b224e63SBarry Smith } else { 16665e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1673b224e63SBarry Smith } 1683b224e63SBarry Smith PetscFunctionReturn(0); 1693b224e63SBarry Smith } 1703b224e63SBarry Smith 1713b224e63SBarry Smith #undef __FUNCT__ 1726c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1736c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1746c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1756c924f48SJed Brown { 1766c924f48SJed Brown PetscErrorCode ierr; 1776c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1785d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 1795d4c12cdSJungho Lee PetscBool flg,flg_col; 1805d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 1816c924f48SJed Brown 1826c924f48SJed Brown PetscFunctionBegin; 1836c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 1845d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 1856c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 1866c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 1876c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 1885d4c12cdSJungho Lee ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 1896c924f48SJed Brown nfields = jac->bs; 190*29499fbbSJungho Lee nfields_col = jac->bs; 1916c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 1925d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 1936c924f48SJed Brown if (!flg) break; 1945d4c12cdSJungho Lee else if (flg && !flg_col) { 1955d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1965d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 1975d4c12cdSJungho Lee } 1985d4c12cdSJungho Lee else { 1995d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2005d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2015d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2025d4c12cdSJungho Lee } 2036c924f48SJed Brown } 2046c924f48SJed Brown if (i > 0) { 2056c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2066c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2076c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2086c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2096c924f48SJed Brown } 2106c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2115d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2126c924f48SJed Brown PetscFunctionReturn(0); 2136c924f48SJed Brown } 2146c924f48SJed Brown 2156c924f48SJed Brown #undef __FUNCT__ 21669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 21769a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2180971522cSBarry Smith { 2190971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2200971522cSBarry Smith PetscErrorCode ierr; 2215a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2226ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2236c924f48SJed Brown PetscInt i; 2240971522cSBarry Smith 2250971522cSBarry Smith PetscFunctionBegin; 226d32f9abdSBarry Smith if (!ilink) { 227d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 228d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2294d343eeaSMatthew G Knepley PetscBool dmcomposite; 2304d343eeaSMatthew G Knepley ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 2314d343eeaSMatthew G Knepley if (dmcomposite) { 2324d343eeaSMatthew G Knepley PetscInt nDM; 2335d4c12cdSJungho Lee IS *fields; 2348b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2358b8307b2SJed Brown ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 2365d4c12cdSJungho Lee ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 2375d4c12cdSJungho Lee for (i=0; i<nDM; i++) { 2385d4c12cdSJungho Lee char splitname[8]; 2395d4c12cdSJungho Lee ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2405d4c12cdSJungho Lee ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr); 2415d4c12cdSJungho Lee ierr = ISDestroy(&fields[i]);CHKERRQ(ierr); 2422fa5ba8aSJed Brown } 2435d4c12cdSJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 2448b8307b2SJed Brown } 24566ffff09SJed Brown } else { 246521d7252SBarry Smith if (jac->bs <= 0) { 247704ba839SBarry Smith if (pc->pmat) { 248521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 249704ba839SBarry Smith } else { 250704ba839SBarry Smith jac->bs = 1; 251704ba839SBarry Smith } 252521d7252SBarry Smith } 253d32f9abdSBarry Smith 254acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2556ce1633cSBarry Smith if (stokes) { 2566ce1633cSBarry Smith IS zerodiags,rest; 2576ce1633cSBarry Smith PetscInt nmin,nmax; 2586ce1633cSBarry Smith 2596ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2606ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2616ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2626ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2636ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 264fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 265fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2666ce1633cSBarry Smith } else { 267d32f9abdSBarry Smith if (!flg) { 268d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 269d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2706c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2716c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 272d32f9abdSBarry Smith } 2736c924f48SJed Brown if (flg || !jac->splitdefined) { 274d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 275db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2766c924f48SJed Brown char splitname[8]; 2776c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2785d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 27979416396SBarry Smith } 2805d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 281521d7252SBarry Smith } 28266ffff09SJed Brown } 2836ce1633cSBarry Smith } 284edf189efSBarry Smith } else if (jac->nsplits == 1) { 285edf189efSBarry Smith if (ilink->is) { 286edf189efSBarry Smith IS is2; 287edf189efSBarry Smith PetscInt nmin,nmax; 288edf189efSBarry Smith 289edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 290edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 291db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 292fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 2935d4c12cdSJungho Lee } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 29463ec74ffSBarry Smith } else if (jac->reset) { 29563ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 296d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 297d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 298d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 299d0af7cd3SBarry Smith if (pc->dm && !stokes) { 300d0af7cd3SBarry Smith PetscBool dmcomposite; 301d0af7cd3SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 302d0af7cd3SBarry Smith if (dmcomposite) { 303d0af7cd3SBarry Smith PetscInt nDM; 304d0af7cd3SBarry Smith IS *fields; 305d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 306d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 307d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 308d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 309d0af7cd3SBarry Smith ilink->is = fields[i]; 310d0af7cd3SBarry Smith ilink = ilink->next; 311edf189efSBarry Smith } 312d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 313d0af7cd3SBarry Smith } 314d0af7cd3SBarry Smith } else { 315d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 316d0af7cd3SBarry Smith if (stokes) { 317d0af7cd3SBarry Smith IS zerodiags,rest; 318d0af7cd3SBarry Smith PetscInt nmin,nmax; 319d0af7cd3SBarry Smith 320d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 321d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 322d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 32320b26d62SBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 32420b26d62SBarry Smith ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr); 325d0af7cd3SBarry Smith ilink->is = rest; 326d0af7cd3SBarry Smith ilink->next->is = zerodiags; 32720b26d62SBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 328d0af7cd3SBarry Smith } 329d0af7cd3SBarry Smith } 330d0af7cd3SBarry Smith 3315d4c12cdSJungho Lee if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 33269a612a9SBarry Smith PetscFunctionReturn(0); 33369a612a9SBarry Smith } 33469a612a9SBarry Smith 33569a612a9SBarry Smith #undef __FUNCT__ 33669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 33769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 33869a612a9SBarry Smith { 33969a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 34069a612a9SBarry Smith PetscErrorCode ierr; 3415a9f2f41SSatish Balay PC_FieldSplitLink ilink; 342cf502942SBarry Smith PetscInt i,nsplit,ccsize; 34369a612a9SBarry Smith MatStructure flag = pc->flag; 344ace3abfcSBarry Smith PetscBool sorted; 34569a612a9SBarry Smith 34669a612a9SBarry Smith PetscFunctionBegin; 34769a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 34897bbdb24SBarry Smith nsplit = jac->nsplits; 3495a9f2f41SSatish Balay ilink = jac->head; 35097bbdb24SBarry Smith 35197bbdb24SBarry Smith /* get the matrices for each split */ 352704ba839SBarry Smith if (!jac->issetup) { 3531b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 35497bbdb24SBarry Smith 355704ba839SBarry Smith jac->issetup = PETSC_TRUE; 356704ba839SBarry Smith 3575d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 35851f519a2SBarry Smith bs = jac->bs; 35997bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 360cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 3611b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3621b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3631b9fc7fcSBarry Smith if (jac->defaultsplit) { 364704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 365704ba839SBarry Smith } else if (!ilink->is) { 366ccb205f8SBarry Smith if (ilink->nfields > 1) { 3675a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 3685a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3691b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3701b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3711b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 37297bbdb24SBarry Smith } 37397bbdb24SBarry Smith } 374d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 375ccb205f8SBarry Smith } else { 376704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 377ccb205f8SBarry Smith } 3783e197d65SBarry Smith } 379704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 380e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 381704ba839SBarry Smith ilink = ilink->next; 3821b9fc7fcSBarry Smith } 3831b9fc7fcSBarry Smith } 3841b9fc7fcSBarry Smith 385704ba839SBarry Smith ilink = jac->head; 38697bbdb24SBarry Smith if (!jac->pmat) { 387cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 388cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3894aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 390704ba839SBarry Smith ilink = ilink->next; 391cf502942SBarry Smith } 39297bbdb24SBarry Smith } else { 393cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3944aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 395704ba839SBarry Smith ilink = ilink->next; 396cf502942SBarry Smith } 39797bbdb24SBarry Smith } 398519d70e2SJed Brown if (jac->realdiagonal) { 399519d70e2SJed Brown ilink = jac->head; 400519d70e2SJed Brown if (!jac->mat) { 401519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 402519d70e2SJed Brown for (i=0; i<nsplit; i++) { 403519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 404519d70e2SJed Brown ilink = ilink->next; 405519d70e2SJed Brown } 406519d70e2SJed Brown } else { 407519d70e2SJed Brown for (i=0; i<nsplit; i++) { 408966e74d7SJed Brown if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 409519d70e2SJed Brown ilink = ilink->next; 410519d70e2SJed Brown } 411519d70e2SJed Brown } 412519d70e2SJed Brown } else { 413519d70e2SJed Brown jac->mat = jac->pmat; 414519d70e2SJed Brown } 41597bbdb24SBarry Smith 4166c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 41768dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 41868dd23aaSBarry Smith ilink = jac->head; 41968dd23aaSBarry Smith if (!jac->Afield) { 42068dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 42168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4224aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 42368dd23aaSBarry Smith ilink = ilink->next; 42468dd23aaSBarry Smith } 42568dd23aaSBarry Smith } else { 42668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4274aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 42868dd23aaSBarry Smith ilink = ilink->next; 42968dd23aaSBarry Smith } 43068dd23aaSBarry Smith } 43168dd23aaSBarry Smith } 43268dd23aaSBarry Smith 4333b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 4343b224e63SBarry Smith IS ccis; 4354aa3045dSJed Brown PetscInt rstart,rend; 436e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 43768dd23aaSBarry Smith 438e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 439e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 440e6cab6aaSJed Brown 4413b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4423b224e63SBarry Smith if (jac->schur) { 4433b224e63SBarry Smith ilink = jac->head; 4444aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4454aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 446fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4473b224e63SBarry Smith ilink = ilink->next; 4484aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4494aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 450fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 451519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 452084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4533b224e63SBarry Smith 4543b224e63SBarry Smith } else { 4551cee3971SBarry Smith KSP ksp; 4566c924f48SJed Brown char schurprefix[256]; 4573b224e63SBarry Smith 458a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4593b224e63SBarry Smith ilink = jac->head; 4604aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4614aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 462fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4633b224e63SBarry Smith ilink = ilink->next; 4644aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4654aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 466fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4677d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 468519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 469a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4701cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 471e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 472a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 473a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 47420b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 47520b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4763b224e63SBarry Smith 4773b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4789005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4791cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 480084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 481084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 482084e4875SJed Brown PC pc; 483cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 484084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 485084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 486e69d4d44SBarry Smith } 4876c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 4886c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 4893b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 49020b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 49120b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 4923b224e63SBarry Smith 4933b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 4943b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 4953b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 4963b224e63SBarry Smith ilink = jac->head; 4973b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 4983b224e63SBarry Smith ilink = ilink->next; 4993b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 5003b224e63SBarry Smith } 5013b224e63SBarry Smith } else { 50297bbdb24SBarry Smith /* set up the individual PCs */ 50397bbdb24SBarry Smith i = 0; 5045a9f2f41SSatish Balay ilink = jac->head; 5055a9f2f41SSatish Balay while (ilink) { 506519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 5073b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 508c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 5095a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 51097bbdb24SBarry Smith i++; 5115a9f2f41SSatish Balay ilink = ilink->next; 5120971522cSBarry Smith } 51397bbdb24SBarry Smith 51497bbdb24SBarry Smith /* create work vectors for each split */ 5151b9fc7fcSBarry Smith if (!jac->x) { 51697bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5175a9f2f41SSatish Balay ilink = jac->head; 51897bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 519906ed7ccSBarry Smith Vec *vl,*vr; 520906ed7ccSBarry Smith 521906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 522906ed7ccSBarry Smith ilink->x = *vr; 523906ed7ccSBarry Smith ilink->y = *vl; 524906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 525906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5265a9f2f41SSatish Balay jac->x[i] = ilink->x; 5275a9f2f41SSatish Balay jac->y[i] = ilink->y; 5285a9f2f41SSatish Balay ilink = ilink->next; 52997bbdb24SBarry Smith } 5303b224e63SBarry Smith } 5313b224e63SBarry Smith } 5323b224e63SBarry Smith 5333b224e63SBarry Smith 534d0f46423SBarry Smith if (!jac->head->sctx) { 5353b224e63SBarry Smith Vec xtmp; 5363b224e63SBarry Smith 53779416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5381b9fc7fcSBarry Smith 5395a9f2f41SSatish Balay ilink = jac->head; 5401b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5411b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 542704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5435a9f2f41SSatish Balay ilink = ilink->next; 5441b9fc7fcSBarry Smith } 5456bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5461b9fc7fcSBarry Smith } 547c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 5480971522cSBarry Smith PetscFunctionReturn(0); 5490971522cSBarry Smith } 5500971522cSBarry Smith 5515a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 552ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 553ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5545a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 555ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 556ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 55779416396SBarry Smith 5580971522cSBarry Smith #undef __FUNCT__ 5593b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5603b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5613b224e63SBarry Smith { 5623b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5633b224e63SBarry Smith PetscErrorCode ierr; 5643b224e63SBarry Smith KSP ksp; 5653b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5663b224e63SBarry Smith 5673b224e63SBarry Smith PetscFunctionBegin; 5683b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5693b224e63SBarry Smith 570c5d2311dSJed Brown switch (jac->schurfactorization) { 571c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 572a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 573c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 574c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 575c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 576c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 577c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 578c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 580c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 581c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 582c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 583c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 584c5d2311dSJed Brown break; 585c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 586a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 587c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 588c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 589c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 590c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 591c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 592c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 593c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 594c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 595c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 596c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 598c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 599c5d2311dSJed Brown break; 600c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 601a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 602c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 603c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 604c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 605c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 606c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 607c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 608c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 609c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 610c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 611c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 612c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 613c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 614c5d2311dSJed Brown break; 615c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 616a04f6461SBarry 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 */ 6173b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6183b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6193b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6203b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6213b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6223b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6233b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6243b224e63SBarry Smith 6253b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6263b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6273b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6283b224e63SBarry Smith 6293b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6303b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6313b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6323b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6333b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 634c5d2311dSJed Brown } 6353b224e63SBarry Smith PetscFunctionReturn(0); 6363b224e63SBarry Smith } 6373b224e63SBarry Smith 6383b224e63SBarry Smith #undef __FUNCT__ 6390971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6400971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6410971522cSBarry Smith { 6420971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6430971522cSBarry Smith PetscErrorCode ierr; 6445a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 645af93645bSJed Brown PetscInt cnt; 6460971522cSBarry Smith 6470971522cSBarry Smith PetscFunctionBegin; 64851f519a2SBarry Smith CHKMEMQ; 64951f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 65051f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 65151f519a2SBarry Smith 65279416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6531b9fc7fcSBarry Smith if (jac->defaultsplit) { 6540971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6555a9f2f41SSatish Balay while (ilink) { 6565a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6575a9f2f41SSatish Balay ilink = ilink->next; 6580971522cSBarry Smith } 6590971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6601b9fc7fcSBarry Smith } else { 661efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6625a9f2f41SSatish Balay while (ilink) { 6635a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6645a9f2f41SSatish Balay ilink = ilink->next; 6651b9fc7fcSBarry Smith } 6661b9fc7fcSBarry Smith } 66716913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 66879416396SBarry Smith if (!jac->w1) { 66979416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 67079416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 67179416396SBarry Smith } 672efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6735a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6743e197d65SBarry Smith cnt = 1; 6755a9f2f41SSatish Balay while (ilink->next) { 6765a9f2f41SSatish Balay ilink = ilink->next; 6773e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6783e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6793e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6803e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6813e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6823e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6833e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6843e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6853e197d65SBarry Smith } 68651f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 68711755939SBarry Smith cnt -= 2; 68851f519a2SBarry Smith while (ilink->previous) { 68951f519a2SBarry Smith ilink = ilink->previous; 69011755939SBarry Smith /* compute the residual only over the part of the vector needed */ 69111755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 69211755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 69311755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 69411755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 69511755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 69611755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 69711755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 69851f519a2SBarry Smith } 69911755939SBarry Smith } 70065e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 70151f519a2SBarry Smith CHKMEMQ; 7020971522cSBarry Smith PetscFunctionReturn(0); 7030971522cSBarry Smith } 7040971522cSBarry Smith 705421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 706ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 707ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 708421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 709ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 710ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 711421e10b8SBarry Smith 712421e10b8SBarry Smith #undef __FUNCT__ 7138c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 714421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 715421e10b8SBarry Smith { 716421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 717421e10b8SBarry Smith PetscErrorCode ierr; 718421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 719421e10b8SBarry Smith 720421e10b8SBarry Smith PetscFunctionBegin; 721421e10b8SBarry Smith CHKMEMQ; 722421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 723421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 724421e10b8SBarry Smith 725421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 726421e10b8SBarry Smith if (jac->defaultsplit) { 727421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 728421e10b8SBarry Smith while (ilink) { 729421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 730421e10b8SBarry Smith ilink = ilink->next; 731421e10b8SBarry Smith } 732421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 733421e10b8SBarry Smith } else { 734421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 735421e10b8SBarry Smith while (ilink) { 736421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 737421e10b8SBarry Smith ilink = ilink->next; 738421e10b8SBarry Smith } 739421e10b8SBarry Smith } 740421e10b8SBarry Smith } else { 741421e10b8SBarry Smith if (!jac->w1) { 742421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 743421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 744421e10b8SBarry Smith } 745421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 746421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 747421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 748421e10b8SBarry Smith while (ilink->next) { 749421e10b8SBarry Smith ilink = ilink->next; 7509989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 751421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 752421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 753421e10b8SBarry Smith } 754421e10b8SBarry Smith while (ilink->previous) { 755421e10b8SBarry Smith ilink = ilink->previous; 7569989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 757421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 758421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 759421e10b8SBarry Smith } 760421e10b8SBarry Smith } else { 761421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 762421e10b8SBarry Smith ilink = ilink->next; 763421e10b8SBarry Smith } 764421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 765421e10b8SBarry Smith while (ilink->previous) { 766421e10b8SBarry Smith ilink = ilink->previous; 7679989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 768421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 769421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 770421e10b8SBarry Smith } 771421e10b8SBarry Smith } 772421e10b8SBarry Smith } 773421e10b8SBarry Smith CHKMEMQ; 774421e10b8SBarry Smith PetscFunctionReturn(0); 775421e10b8SBarry Smith } 776421e10b8SBarry Smith 7770971522cSBarry Smith #undef __FUNCT__ 778574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 779574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 7800971522cSBarry Smith { 7810971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7820971522cSBarry Smith PetscErrorCode ierr; 7835a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 7840971522cSBarry Smith 7850971522cSBarry Smith PetscFunctionBegin; 7865a9f2f41SSatish Balay while (ilink) { 787574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 788fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 789fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 790fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 791fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 7925a9f2f41SSatish Balay next = ilink->next; 7935a9f2f41SSatish Balay ilink = next; 7940971522cSBarry Smith } 79505b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 796574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 797574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 798574deadeSBarry Smith } else if (jac->mat) { 799574deadeSBarry Smith jac->mat = PETSC_NULL; 800574deadeSBarry Smith } 80197bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 80268dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8036bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8046bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8056bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8066bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8076bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8086bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8096bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 81063ec74ffSBarry Smith jac->reset = PETSC_TRUE; 811574deadeSBarry Smith PetscFunctionReturn(0); 812574deadeSBarry Smith } 813574deadeSBarry Smith 814574deadeSBarry Smith #undef __FUNCT__ 815574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 816574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 817574deadeSBarry Smith { 818574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 819574deadeSBarry Smith PetscErrorCode ierr; 820574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 821574deadeSBarry Smith 822574deadeSBarry Smith PetscFunctionBegin; 823574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 824574deadeSBarry Smith while (ilink) { 8256bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 826574deadeSBarry Smith next = ilink->next; 827574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 828574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 8295d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 830574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 831574deadeSBarry Smith ilink = next; 832574deadeSBarry Smith } 833574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 834c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8350971522cSBarry Smith PetscFunctionReturn(0); 8360971522cSBarry Smith } 8370971522cSBarry Smith 8380971522cSBarry Smith #undef __FUNCT__ 8390971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8400971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8410971522cSBarry Smith { 8421b9fc7fcSBarry Smith PetscErrorCode ierr; 8436c924f48SJed Brown PetscInt bs; 844bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8459dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8463b224e63SBarry Smith PCCompositeType ctype; 8471b9fc7fcSBarry Smith 8480971522cSBarry Smith PetscFunctionBegin; 8491b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 850acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 85151f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 85251f519a2SBarry Smith if (flg) { 85351f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 85451f519a2SBarry Smith } 855704ba839SBarry Smith 856435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 857c0adfefeSBarry Smith if (stokes) { 858c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 859c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 860c0adfefeSBarry Smith } 861c0adfefeSBarry Smith 8623b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 8633b224e63SBarry Smith if (flg) { 8643b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 8653b224e63SBarry Smith } 866d32f9abdSBarry Smith 867c30613efSMatthew Knepley /* Only setup fields once */ 868c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 869d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 870d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 8716c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 8726c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 873d32f9abdSBarry Smith } 874c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 875c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 876084e4875SJed 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); 877c5d2311dSJed Brown } 8781b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8790971522cSBarry Smith PetscFunctionReturn(0); 8800971522cSBarry Smith } 8810971522cSBarry Smith 8820971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 8830971522cSBarry Smith 8840971522cSBarry Smith EXTERN_C_BEGIN 8850971522cSBarry Smith #undef __FUNCT__ 8860971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 8875d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 8880971522cSBarry Smith { 88997bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8900971522cSBarry Smith PetscErrorCode ierr; 8915a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 89269a612a9SBarry Smith char prefix[128]; 8935d4c12cdSJungho Lee PetscInt i; 8940971522cSBarry Smith 8950971522cSBarry Smith PetscFunctionBegin; 8966c924f48SJed Brown if (jac->splitdefined) { 8976c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 8986c924f48SJed Brown PetscFunctionReturn(0); 8996c924f48SJed Brown } 90051f519a2SBarry Smith for (i=0; i<n; i++) { 901e32f2f54SBarry 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); 902e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 90351f519a2SBarry Smith } 904704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 905a04f6461SBarry Smith if (splitname) { 906db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 907a04f6461SBarry Smith } else { 908a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 909a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 910a04f6461SBarry Smith } 911704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 9125a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 9135d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 9145d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 9155a9f2f41SSatish Balay ilink->nfields = n; 9165a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9177adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9181cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9195a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9209005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 92169a612a9SBarry Smith 922a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9235a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9240971522cSBarry Smith 9250971522cSBarry Smith if (!next) { 9265a9f2f41SSatish Balay jac->head = ilink; 92751f519a2SBarry Smith ilink->previous = PETSC_NULL; 9280971522cSBarry Smith } else { 9290971522cSBarry Smith while (next->next) { 9300971522cSBarry Smith next = next->next; 9310971522cSBarry Smith } 9325a9f2f41SSatish Balay next->next = ilink; 93351f519a2SBarry Smith ilink->previous = next; 9340971522cSBarry Smith } 9350971522cSBarry Smith jac->nsplits++; 9360971522cSBarry Smith PetscFunctionReturn(0); 9370971522cSBarry Smith } 9380971522cSBarry Smith EXTERN_C_END 9390971522cSBarry Smith 940e69d4d44SBarry Smith EXTERN_C_BEGIN 941e69d4d44SBarry Smith #undef __FUNCT__ 942e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9437087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 944e69d4d44SBarry Smith { 945e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 946e69d4d44SBarry Smith PetscErrorCode ierr; 947e69d4d44SBarry Smith 948e69d4d44SBarry Smith PetscFunctionBegin; 949519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 950e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 951e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 95213e0d083SBarry Smith if (n) *n = jac->nsplits; 953e69d4d44SBarry Smith PetscFunctionReturn(0); 954e69d4d44SBarry Smith } 955e69d4d44SBarry Smith EXTERN_C_END 9560971522cSBarry Smith 9570971522cSBarry Smith EXTERN_C_BEGIN 9580971522cSBarry Smith #undef __FUNCT__ 95969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 9607087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 9610971522cSBarry Smith { 9620971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9630971522cSBarry Smith PetscErrorCode ierr; 9640971522cSBarry Smith PetscInt cnt = 0; 9655a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 9660971522cSBarry Smith 9670971522cSBarry Smith PetscFunctionBegin; 9685d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 9695a9f2f41SSatish Balay while (ilink) { 9705a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 9715a9f2f41SSatish Balay ilink = ilink->next; 9720971522cSBarry Smith } 9735d480477SMatthew 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); 97413e0d083SBarry Smith if (n) *n = jac->nsplits; 9750971522cSBarry Smith PetscFunctionReturn(0); 9760971522cSBarry Smith } 9770971522cSBarry Smith EXTERN_C_END 9780971522cSBarry Smith 979704ba839SBarry Smith EXTERN_C_BEGIN 980704ba839SBarry Smith #undef __FUNCT__ 981704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 9827087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 983704ba839SBarry Smith { 984704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 985704ba839SBarry Smith PetscErrorCode ierr; 986704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 987704ba839SBarry Smith char prefix[128]; 988704ba839SBarry Smith 989704ba839SBarry Smith PetscFunctionBegin; 9906c924f48SJed Brown if (jac->splitdefined) { 9916c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9926c924f48SJed Brown PetscFunctionReturn(0); 9936c924f48SJed Brown } 99416913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 995a04f6461SBarry Smith if (splitname) { 996db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 997a04f6461SBarry Smith } else { 998a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 9995d4c12cdSJungho Lee ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1000a04f6461SBarry Smith } 10011ab39975SBarry Smith ilink->is = is; 10021ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1003704ba839SBarry Smith ilink->next = PETSC_NULL; 1004704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10051cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1006704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10079005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1008704ba839SBarry Smith 1009a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1010704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1011704ba839SBarry Smith 1012704ba839SBarry Smith if (!next) { 1013704ba839SBarry Smith jac->head = ilink; 1014704ba839SBarry Smith ilink->previous = PETSC_NULL; 1015704ba839SBarry Smith } else { 1016704ba839SBarry Smith while (next->next) { 1017704ba839SBarry Smith next = next->next; 1018704ba839SBarry Smith } 1019704ba839SBarry Smith next->next = ilink; 1020704ba839SBarry Smith ilink->previous = next; 1021704ba839SBarry Smith } 1022704ba839SBarry Smith jac->nsplits++; 1023704ba839SBarry Smith 1024704ba839SBarry Smith PetscFunctionReturn(0); 1025704ba839SBarry Smith } 1026704ba839SBarry Smith EXTERN_C_END 1027704ba839SBarry Smith 10280971522cSBarry Smith #undef __FUNCT__ 10290971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10300971522cSBarry Smith /*@ 10310971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10320971522cSBarry Smith 1033ad4df100SBarry Smith Logically Collective on PC 10340971522cSBarry Smith 10350971522cSBarry Smith Input Parameters: 10360971522cSBarry Smith + pc - the preconditioner context 1037a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10380971522cSBarry Smith . n - the number of fields in this split 1039db4c96c1SJed Brown - fields - the fields in this split 10400971522cSBarry Smith 10410971522cSBarry Smith Level: intermediate 10420971522cSBarry Smith 1043d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1044d32f9abdSBarry Smith 1045d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1046d32f9abdSBarry 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 1047d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1048d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1049d32f9abdSBarry Smith 1050db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1051db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1052db4c96c1SJed Brown 10535d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 10545d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 10555d4c12cdSJungho Lee available when this routine is called. 10565d4c12cdSJungho Lee 1057d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 10580971522cSBarry Smith 10590971522cSBarry Smith @*/ 10605d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10610971522cSBarry Smith { 10624ac538c5SBarry Smith PetscErrorCode ierr; 10630971522cSBarry Smith 10640971522cSBarry Smith PetscFunctionBegin; 10650700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1066db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1067db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1068db4c96c1SJed Brown PetscValidIntPointer(fields,3); 10695d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 10700971522cSBarry Smith PetscFunctionReturn(0); 10710971522cSBarry Smith } 10720971522cSBarry Smith 10730971522cSBarry Smith #undef __FUNCT__ 1074704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1075704ba839SBarry Smith /*@ 1076704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1077704ba839SBarry Smith 1078ad4df100SBarry Smith Logically Collective on PC 1079704ba839SBarry Smith 1080704ba839SBarry Smith Input Parameters: 1081704ba839SBarry Smith + pc - the preconditioner context 1082a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1083db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1084704ba839SBarry Smith 1085d32f9abdSBarry Smith 1086a6ffb8dbSJed Brown Notes: 1087a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1088a6ffb8dbSJed Brown 1089db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1090db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1091d32f9abdSBarry Smith 1092704ba839SBarry Smith Level: intermediate 1093704ba839SBarry Smith 1094704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1095704ba839SBarry Smith 1096704ba839SBarry Smith @*/ 10977087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1098704ba839SBarry Smith { 10994ac538c5SBarry Smith PetscErrorCode ierr; 1100704ba839SBarry Smith 1101704ba839SBarry Smith PetscFunctionBegin; 11020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11035d4c12cdSJungho Lee PetscValidCharPointer(splitname,2); 1104db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11054ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1106704ba839SBarry Smith PetscFunctionReturn(0); 1107704ba839SBarry Smith } 1108704ba839SBarry Smith 1109704ba839SBarry Smith #undef __FUNCT__ 111057a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 111157a9adfeSMatthew G Knepley /*@ 111257a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 111357a9adfeSMatthew G Knepley 111457a9adfeSMatthew G Knepley Logically Collective on PC 111557a9adfeSMatthew G Knepley 111657a9adfeSMatthew G Knepley Input Parameters: 111757a9adfeSMatthew G Knepley + pc - the preconditioner context 111857a9adfeSMatthew G Knepley - splitname - name of this split 111957a9adfeSMatthew G Knepley 112057a9adfeSMatthew G Knepley Output Parameter: 112157a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 112257a9adfeSMatthew G Knepley 112357a9adfeSMatthew G Knepley Level: intermediate 112457a9adfeSMatthew G Knepley 112557a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 112657a9adfeSMatthew G Knepley 112757a9adfeSMatthew G Knepley @*/ 112857a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 112957a9adfeSMatthew G Knepley { 113057a9adfeSMatthew G Knepley PetscErrorCode ierr; 113157a9adfeSMatthew G Knepley 113257a9adfeSMatthew G Knepley PetscFunctionBegin; 113357a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 113457a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 113557a9adfeSMatthew G Knepley PetscValidPointer(is,3); 113657a9adfeSMatthew G Knepley { 113757a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 113857a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 113957a9adfeSMatthew G Knepley PetscBool found; 114057a9adfeSMatthew G Knepley 114157a9adfeSMatthew G Knepley *is = PETSC_NULL; 114257a9adfeSMatthew G Knepley while(ilink) { 114357a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 114457a9adfeSMatthew G Knepley if (found) { 114557a9adfeSMatthew G Knepley *is = ilink->is; 114657a9adfeSMatthew G Knepley break; 114757a9adfeSMatthew G Knepley } 114857a9adfeSMatthew G Knepley ilink = ilink->next; 114957a9adfeSMatthew G Knepley } 115057a9adfeSMatthew G Knepley } 115157a9adfeSMatthew G Knepley PetscFunctionReturn(0); 115257a9adfeSMatthew G Knepley } 115357a9adfeSMatthew G Knepley 115457a9adfeSMatthew G Knepley #undef __FUNCT__ 115551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 115651f519a2SBarry Smith /*@ 115751f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 115851f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 115951f519a2SBarry Smith 1160ad4df100SBarry Smith Logically Collective on PC 116151f519a2SBarry Smith 116251f519a2SBarry Smith Input Parameters: 116351f519a2SBarry Smith + pc - the preconditioner context 116451f519a2SBarry Smith - bs - the block size 116551f519a2SBarry Smith 116651f519a2SBarry Smith Level: intermediate 116751f519a2SBarry Smith 116851f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 116951f519a2SBarry Smith 117051f519a2SBarry Smith @*/ 11717087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 117251f519a2SBarry Smith { 11734ac538c5SBarry Smith PetscErrorCode ierr; 117451f519a2SBarry Smith 117551f519a2SBarry Smith PetscFunctionBegin; 11760700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1177c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 11784ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 117951f519a2SBarry Smith PetscFunctionReturn(0); 118051f519a2SBarry Smith } 118151f519a2SBarry Smith 118251f519a2SBarry Smith #undef __FUNCT__ 118369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 11840971522cSBarry Smith /*@C 118569a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 11860971522cSBarry Smith 118769a612a9SBarry Smith Collective on KSP 11880971522cSBarry Smith 11890971522cSBarry Smith Input Parameter: 11900971522cSBarry Smith . pc - the preconditioner context 11910971522cSBarry Smith 11920971522cSBarry Smith Output Parameters: 119313e0d083SBarry Smith + n - the number of splits 119469a612a9SBarry Smith - pc - the array of KSP contexts 11950971522cSBarry Smith 11960971522cSBarry Smith Note: 1197d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1198d32f9abdSBarry Smith (not the KSP just the array that contains them). 11990971522cSBarry Smith 120069a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12010971522cSBarry Smith 12020971522cSBarry Smith Level: advanced 12030971522cSBarry Smith 12040971522cSBarry Smith .seealso: PCFIELDSPLIT 12050971522cSBarry Smith @*/ 12067087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12070971522cSBarry Smith { 12084ac538c5SBarry Smith PetscErrorCode ierr; 12090971522cSBarry Smith 12100971522cSBarry Smith PetscFunctionBegin; 12110700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 121213e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 12134ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 12140971522cSBarry Smith PetscFunctionReturn(0); 12150971522cSBarry Smith } 12160971522cSBarry Smith 1217e69d4d44SBarry Smith #undef __FUNCT__ 1218e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1219e69d4d44SBarry Smith /*@ 1220e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1221a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1222e69d4d44SBarry Smith 1223e69d4d44SBarry Smith Collective on PC 1224e69d4d44SBarry Smith 1225e69d4d44SBarry Smith Input Parameters: 1226e69d4d44SBarry Smith + pc - the preconditioner context 1227fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1228084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1229084e4875SJed Brown 1230e69d4d44SBarry Smith Options Database: 1231084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1232e69d4d44SBarry Smith 1233fd1303e9SJungho Lee Notes: 1234fd1303e9SJungho Lee $ If ptype is 1235fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1236fd1303e9SJungho Lee $ to this function). 1237fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1238fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1239fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1240fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1241fd1303e9SJungho Lee $ preconditioner 1242fd1303e9SJungho Lee 1243fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1244fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1245fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1246fd1303e9SJungho Lee 1247fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1248fd1303e9SJungho Lee the name since it sets a proceedure to use. 1249fd1303e9SJungho Lee 1250e69d4d44SBarry Smith Level: intermediate 1251e69d4d44SBarry Smith 1252fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1253e69d4d44SBarry Smith 1254e69d4d44SBarry Smith @*/ 12557087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1256e69d4d44SBarry Smith { 12574ac538c5SBarry Smith PetscErrorCode ierr; 1258e69d4d44SBarry Smith 1259e69d4d44SBarry Smith PetscFunctionBegin; 12600700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12614ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1262e69d4d44SBarry Smith PetscFunctionReturn(0); 1263e69d4d44SBarry Smith } 1264e69d4d44SBarry Smith 1265e69d4d44SBarry Smith EXTERN_C_BEGIN 1266e69d4d44SBarry Smith #undef __FUNCT__ 1267e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 12687087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1269e69d4d44SBarry Smith { 1270e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1271084e4875SJed Brown PetscErrorCode ierr; 1272e69d4d44SBarry Smith 1273e69d4d44SBarry Smith PetscFunctionBegin; 1274084e4875SJed Brown jac->schurpre = ptype; 1275084e4875SJed Brown if (pre) { 12766bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1277084e4875SJed Brown jac->schur_user = pre; 1278084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1279084e4875SJed Brown } 1280e69d4d44SBarry Smith PetscFunctionReturn(0); 1281e69d4d44SBarry Smith } 1282e69d4d44SBarry Smith EXTERN_C_END 1283e69d4d44SBarry Smith 128430ad9308SMatthew Knepley #undef __FUNCT__ 128530ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 128630ad9308SMatthew Knepley /*@C 12878c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 128830ad9308SMatthew Knepley 128930ad9308SMatthew Knepley Collective on KSP 129030ad9308SMatthew Knepley 129130ad9308SMatthew Knepley Input Parameter: 129230ad9308SMatthew Knepley . pc - the preconditioner context 129330ad9308SMatthew Knepley 129430ad9308SMatthew Knepley Output Parameters: 1295a04f6461SBarry Smith + A00 - the (0,0) block 1296a04f6461SBarry Smith . A01 - the (0,1) block 1297a04f6461SBarry Smith . A10 - the (1,0) block 1298a04f6461SBarry Smith - A11 - the (1,1) block 129930ad9308SMatthew Knepley 130030ad9308SMatthew Knepley Level: advanced 130130ad9308SMatthew Knepley 130230ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 130330ad9308SMatthew Knepley @*/ 1304a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 130530ad9308SMatthew Knepley { 130630ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 130730ad9308SMatthew Knepley 130830ad9308SMatthew Knepley PetscFunctionBegin; 13090700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1310c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1311a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1312a04f6461SBarry Smith if (A01) *A01 = jac->B; 1313a04f6461SBarry Smith if (A10) *A10 = jac->C; 1314a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 131530ad9308SMatthew Knepley PetscFunctionReturn(0); 131630ad9308SMatthew Knepley } 131730ad9308SMatthew Knepley 131879416396SBarry Smith EXTERN_C_BEGIN 131979416396SBarry Smith #undef __FUNCT__ 132079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 13217087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 132279416396SBarry Smith { 132379416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1324e69d4d44SBarry Smith PetscErrorCode ierr; 132579416396SBarry Smith 132679416396SBarry Smith PetscFunctionBegin; 132779416396SBarry Smith jac->type = type; 13283b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 13293b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 13303b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1331e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1332e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1333e69d4d44SBarry Smith 13343b224e63SBarry Smith } else { 13353b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 13363b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1337e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13389e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 13393b224e63SBarry Smith } 134079416396SBarry Smith PetscFunctionReturn(0); 134179416396SBarry Smith } 134279416396SBarry Smith EXTERN_C_END 134379416396SBarry Smith 134451f519a2SBarry Smith EXTERN_C_BEGIN 134551f519a2SBarry Smith #undef __FUNCT__ 134651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 13477087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 134851f519a2SBarry Smith { 134951f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 135051f519a2SBarry Smith 135151f519a2SBarry Smith PetscFunctionBegin; 135265e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 135365e19b50SBarry 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); 135451f519a2SBarry Smith jac->bs = bs; 135551f519a2SBarry Smith PetscFunctionReturn(0); 135651f519a2SBarry Smith } 135751f519a2SBarry Smith EXTERN_C_END 135851f519a2SBarry Smith 135979416396SBarry Smith #undef __FUNCT__ 136079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1361bc08b0f1SBarry Smith /*@ 136279416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 136379416396SBarry Smith 136479416396SBarry Smith Collective on PC 136579416396SBarry Smith 136679416396SBarry Smith Input Parameter: 136779416396SBarry Smith . pc - the preconditioner context 136881540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 136979416396SBarry Smith 137079416396SBarry Smith Options Database Key: 1371a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 137279416396SBarry Smith 137379416396SBarry Smith Level: Developer 137479416396SBarry Smith 137579416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 137679416396SBarry Smith 137779416396SBarry Smith .seealso: PCCompositeSetType() 137879416396SBarry Smith 137979416396SBarry Smith @*/ 13807087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 138179416396SBarry Smith { 13824ac538c5SBarry Smith PetscErrorCode ierr; 138379416396SBarry Smith 138479416396SBarry Smith PetscFunctionBegin; 13850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13864ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 138779416396SBarry Smith PetscFunctionReturn(0); 138879416396SBarry Smith } 138979416396SBarry Smith 13900971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 13910971522cSBarry Smith /*MC 1392a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1393a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 13940971522cSBarry Smith 1395edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1396edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 13970971522cSBarry Smith 1398a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 139969a612a9SBarry Smith and set the options directly on the resulting KSP object 14000971522cSBarry Smith 14010971522cSBarry Smith Level: intermediate 14020971522cSBarry Smith 140379416396SBarry Smith Options Database Keys: 140481540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 140581540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 140681540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 140781540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 140881540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1409e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1410435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 141179416396SBarry Smith 14125d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 14135d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 14145d4c12cdSJungho Lee 14155d4c12cdSJungho Lee Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1416d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1417d32f9abdSBarry Smith 1418d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1419d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1420d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1421d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1422d32f9abdSBarry Smith 14235d4c12cdSJungho Lee For the Schur complement preconditioner if J = ( A00 A01 ) 14245d4c12cdSJungho Lee ( A10 A11 ) 14255d4c12cdSJungho Lee the preconditioner is 14265d4c12cdSJungho Lee (I -B ksp(A00)) ( inv(A00) 0 (I 0 ) 14275d4c12cdSJungho Lee (0 I ) ( 0 ksp(S) ) (-A10 ksp(A00) I ) 1428a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 14295d4c12cdSJungho Lee ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given in providing the SECOND split or 1 if not give). 14305d4c12cdSJungho Lee For PCFieldSplitGetKSP() when field number is 14315d4c12cdSJungho Lee 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1432a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 14335d4c12cdSJungho Lee option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc 1434e69d4d44SBarry Smith 1435edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1436edf189efSBarry Smith is used automatically for a second block. 1437edf189efSBarry Smith 1438ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1439ff218e97SBarry Smith Generally it should be used with the AIJ format. 1440ff218e97SBarry Smith 1441ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1442ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1443ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 14440716a85fSBarry Smith 1445a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 14460971522cSBarry Smith 14477e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1448e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 14490971522cSBarry Smith M*/ 14500971522cSBarry Smith 14510971522cSBarry Smith EXTERN_C_BEGIN 14520971522cSBarry Smith #undef __FUNCT__ 14530971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 14547087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 14550971522cSBarry Smith { 14560971522cSBarry Smith PetscErrorCode ierr; 14570971522cSBarry Smith PC_FieldSplit *jac; 14580971522cSBarry Smith 14590971522cSBarry Smith PetscFunctionBegin; 146038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 14610971522cSBarry Smith jac->bs = -1; 14620971522cSBarry Smith jac->nsplits = 0; 14633e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1464e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1465c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 146651f519a2SBarry Smith 14670971522cSBarry Smith pc->data = (void*)jac; 14680971522cSBarry Smith 14690971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1470421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 14710971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1472574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 14730971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 14740971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 14750971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 14760971522cSBarry Smith pc->ops->applyrichardson = 0; 14770971522cSBarry Smith 147869a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 147969a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14800971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 14810971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1482704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1483704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 148479416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 148579416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 148651f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 148751f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 14880971522cSBarry Smith PetscFunctionReturn(0); 14890971522cSBarry Smith } 14900971522cSBarry Smith EXTERN_C_END 14910971522cSBarry Smith 14920971522cSBarry Smith 1493a541d17aSBarry Smith 1494