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; 19029499fbbSJungho 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; 344*5f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 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); 365*5f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 366704ba839SBarry Smith } else if (!ilink->is) { 367ccb205f8SBarry Smith if (ilink->nfields > 1) { 368*5f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 3695a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 370*5f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 3711b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3721b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3731b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 374*5f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 37597bbdb24SBarry Smith } 37697bbdb24SBarry Smith } 377d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 378*5f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 379ccb205f8SBarry Smith } else { 380704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 381*5f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 382ccb205f8SBarry Smith } 3833e197d65SBarry Smith } 384704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 385*5f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 386*5f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 387704ba839SBarry Smith ilink = ilink->next; 3881b9fc7fcSBarry Smith } 3891b9fc7fcSBarry Smith } 3901b9fc7fcSBarry Smith 391704ba839SBarry Smith ilink = jac->head; 39297bbdb24SBarry Smith if (!jac->pmat) { 393cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 394cf502942SBarry Smith for (i=0; i<nsplit; i++) { 395*5f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 396704ba839SBarry Smith ilink = ilink->next; 397cf502942SBarry Smith } 39897bbdb24SBarry Smith } else { 399cf502942SBarry Smith for (i=0; i<nsplit; i++) { 400*5f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 401704ba839SBarry Smith ilink = ilink->next; 402cf502942SBarry Smith } 40397bbdb24SBarry Smith } 404519d70e2SJed Brown if (jac->realdiagonal) { 405519d70e2SJed Brown ilink = jac->head; 406519d70e2SJed Brown if (!jac->mat) { 407519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 408519d70e2SJed Brown for (i=0; i<nsplit; i++) { 409*5f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 410519d70e2SJed Brown ilink = ilink->next; 411519d70e2SJed Brown } 412519d70e2SJed Brown } else { 413519d70e2SJed Brown for (i=0; i<nsplit; i++) { 414*5f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 415519d70e2SJed Brown ilink = ilink->next; 416519d70e2SJed Brown } 417519d70e2SJed Brown } 418519d70e2SJed Brown } else { 419519d70e2SJed Brown jac->mat = jac->pmat; 420519d70e2SJed Brown } 42197bbdb24SBarry Smith 4226c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 42368dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 42468dd23aaSBarry Smith ilink = jac->head; 42568dd23aaSBarry Smith if (!jac->Afield) { 42668dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 42768dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4284aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 42968dd23aaSBarry Smith ilink = ilink->next; 43068dd23aaSBarry Smith } 43168dd23aaSBarry Smith } else { 43268dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4334aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 43468dd23aaSBarry Smith ilink = ilink->next; 43568dd23aaSBarry Smith } 43668dd23aaSBarry Smith } 43768dd23aaSBarry Smith } 43868dd23aaSBarry Smith 4393b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 4403b224e63SBarry Smith IS ccis; 4414aa3045dSJed Brown PetscInt rstart,rend; 442e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 44368dd23aaSBarry Smith 444e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 445e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 446e6cab6aaSJed Brown 4473b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4483b224e63SBarry Smith if (jac->schur) { 4493b224e63SBarry Smith ilink = jac->head; 4504aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4514aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 452fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4533b224e63SBarry Smith ilink = ilink->next; 4544aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4554aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 456fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 457519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 458084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4593b224e63SBarry Smith 4603b224e63SBarry Smith } else { 4611cee3971SBarry Smith KSP ksp; 4626c924f48SJed Brown char schurprefix[256]; 4633b224e63SBarry Smith 464a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4653b224e63SBarry Smith ilink = jac->head; 4664aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4674aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 468fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4693b224e63SBarry Smith ilink = ilink->next; 4704aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4714aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 472fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4737d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 474519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 475a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4761cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 477e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 478a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 479a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 48020b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 48120b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4823b224e63SBarry Smith 4833b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4849005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4851cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 486084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 487084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 488084e4875SJed Brown PC pc; 489cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 490084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 491084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 492e69d4d44SBarry Smith } 4936c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 4946c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 4953b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 49620b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 49720b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 4983b224e63SBarry Smith 4993b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 5003b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 5013b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 5023b224e63SBarry Smith ilink = jac->head; 5033b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 5043b224e63SBarry Smith ilink = ilink->next; 5053b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 5063b224e63SBarry Smith } 5073b224e63SBarry Smith } else { 50897bbdb24SBarry Smith /* set up the individual PCs */ 50997bbdb24SBarry Smith i = 0; 5105a9f2f41SSatish Balay ilink = jac->head; 5115a9f2f41SSatish Balay while (ilink) { 512519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 5133b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 514c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 5155a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 51697bbdb24SBarry Smith i++; 5175a9f2f41SSatish Balay ilink = ilink->next; 5180971522cSBarry Smith } 51997bbdb24SBarry Smith 52097bbdb24SBarry Smith /* create work vectors for each split */ 5211b9fc7fcSBarry Smith if (!jac->x) { 52297bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5235a9f2f41SSatish Balay ilink = jac->head; 52497bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 525906ed7ccSBarry Smith Vec *vl,*vr; 526906ed7ccSBarry Smith 527906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 528906ed7ccSBarry Smith ilink->x = *vr; 529906ed7ccSBarry Smith ilink->y = *vl; 530906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 531906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5325a9f2f41SSatish Balay jac->x[i] = ilink->x; 5335a9f2f41SSatish Balay jac->y[i] = ilink->y; 5345a9f2f41SSatish Balay ilink = ilink->next; 53597bbdb24SBarry Smith } 5363b224e63SBarry Smith } 5373b224e63SBarry Smith } 5383b224e63SBarry Smith 5393b224e63SBarry Smith 540d0f46423SBarry Smith if (!jac->head->sctx) { 5413b224e63SBarry Smith Vec xtmp; 5423b224e63SBarry Smith 54379416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5441b9fc7fcSBarry Smith 5455a9f2f41SSatish Balay ilink = jac->head; 5461b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5471b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 548704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5495a9f2f41SSatish Balay ilink = ilink->next; 5501b9fc7fcSBarry Smith } 5516bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5521b9fc7fcSBarry Smith } 553c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 5540971522cSBarry Smith PetscFunctionReturn(0); 5550971522cSBarry Smith } 5560971522cSBarry Smith 5575a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 558ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 559ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5605a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 561ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 562ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 56379416396SBarry Smith 5640971522cSBarry Smith #undef __FUNCT__ 5653b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5663b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5673b224e63SBarry Smith { 5683b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5693b224e63SBarry Smith PetscErrorCode ierr; 5703b224e63SBarry Smith KSP ksp; 5713b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5723b224e63SBarry Smith 5733b224e63SBarry Smith PetscFunctionBegin; 5743b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5753b224e63SBarry Smith 576c5d2311dSJed Brown switch (jac->schurfactorization) { 577c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 578a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 579c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 580c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 581c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 582c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 583c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 584c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 585c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 586c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 587c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 588c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 589c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 590c5d2311dSJed Brown break; 591c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 592a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 593c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 594c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 595c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 596c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 598c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 599c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 600c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 601c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 602c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 603c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 604c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 605c5d2311dSJed Brown break; 606c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 607a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 608c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 609c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 610c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 611c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 612c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 613c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 614c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 615c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 616c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 617c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 618c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 619c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 620c5d2311dSJed Brown break; 621c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 622a04f6461SBarry 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 */ 6233b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6243b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6253b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6263b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6273b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6283b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6293b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6303b224e63SBarry Smith 6313b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6323b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6333b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6343b224e63SBarry Smith 6353b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6363b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6373b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6383b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6393b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 640c5d2311dSJed Brown } 6413b224e63SBarry Smith PetscFunctionReturn(0); 6423b224e63SBarry Smith } 6433b224e63SBarry Smith 6443b224e63SBarry Smith #undef __FUNCT__ 6450971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6460971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6470971522cSBarry Smith { 6480971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6490971522cSBarry Smith PetscErrorCode ierr; 6505a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 651af93645bSJed Brown PetscInt cnt; 6520971522cSBarry Smith 6530971522cSBarry Smith PetscFunctionBegin; 65451f519a2SBarry Smith CHKMEMQ; 65551f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 65651f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 65751f519a2SBarry Smith 65879416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6591b9fc7fcSBarry Smith if (jac->defaultsplit) { 6600971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6615a9f2f41SSatish Balay while (ilink) { 6625a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6635a9f2f41SSatish Balay ilink = ilink->next; 6640971522cSBarry Smith } 6650971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6661b9fc7fcSBarry Smith } else { 667efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6685a9f2f41SSatish Balay while (ilink) { 6695a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6705a9f2f41SSatish Balay ilink = ilink->next; 6711b9fc7fcSBarry Smith } 6721b9fc7fcSBarry Smith } 67316913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 67479416396SBarry Smith if (!jac->w1) { 67579416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 67679416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 67779416396SBarry Smith } 678efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6795a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6803e197d65SBarry Smith cnt = 1; 6815a9f2f41SSatish Balay while (ilink->next) { 6825a9f2f41SSatish Balay ilink = ilink->next; 6833e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6843e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6853e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6863e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6873e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6883e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6893e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6903e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6913e197d65SBarry Smith } 69251f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 69311755939SBarry Smith cnt -= 2; 69451f519a2SBarry Smith while (ilink->previous) { 69551f519a2SBarry Smith ilink = ilink->previous; 69611755939SBarry Smith /* compute the residual only over the part of the vector needed */ 69711755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 69811755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 69911755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 70011755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 70111755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 70211755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 70311755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 70451f519a2SBarry Smith } 70511755939SBarry Smith } 70665e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 70751f519a2SBarry Smith CHKMEMQ; 7080971522cSBarry Smith PetscFunctionReturn(0); 7090971522cSBarry Smith } 7100971522cSBarry Smith 711421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 712ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 713ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 714421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 715ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 716ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 717421e10b8SBarry Smith 718421e10b8SBarry Smith #undef __FUNCT__ 7198c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 720421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 721421e10b8SBarry Smith { 722421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 723421e10b8SBarry Smith PetscErrorCode ierr; 724421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 725421e10b8SBarry Smith 726421e10b8SBarry Smith PetscFunctionBegin; 727421e10b8SBarry Smith CHKMEMQ; 728421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 729421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 730421e10b8SBarry Smith 731421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 732421e10b8SBarry Smith if (jac->defaultsplit) { 733421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 734421e10b8SBarry Smith while (ilink) { 735421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 736421e10b8SBarry Smith ilink = ilink->next; 737421e10b8SBarry Smith } 738421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 739421e10b8SBarry Smith } else { 740421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 741421e10b8SBarry Smith while (ilink) { 742421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 743421e10b8SBarry Smith ilink = ilink->next; 744421e10b8SBarry Smith } 745421e10b8SBarry Smith } 746421e10b8SBarry Smith } else { 747421e10b8SBarry Smith if (!jac->w1) { 748421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 749421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 750421e10b8SBarry Smith } 751421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 752421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 753421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 754421e10b8SBarry Smith while (ilink->next) { 755421e10b8SBarry Smith ilink = ilink->next; 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 while (ilink->previous) { 761421e10b8SBarry Smith ilink = ilink->previous; 7629989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 763421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 764421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 765421e10b8SBarry Smith } 766421e10b8SBarry Smith } else { 767421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 768421e10b8SBarry Smith ilink = ilink->next; 769421e10b8SBarry Smith } 770421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 771421e10b8SBarry Smith while (ilink->previous) { 772421e10b8SBarry Smith ilink = ilink->previous; 7739989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 774421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 775421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 776421e10b8SBarry Smith } 777421e10b8SBarry Smith } 778421e10b8SBarry Smith } 779421e10b8SBarry Smith CHKMEMQ; 780421e10b8SBarry Smith PetscFunctionReturn(0); 781421e10b8SBarry Smith } 782421e10b8SBarry Smith 7830971522cSBarry Smith #undef __FUNCT__ 784574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 785574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 7860971522cSBarry Smith { 7870971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7880971522cSBarry Smith PetscErrorCode ierr; 7895a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 7900971522cSBarry Smith 7910971522cSBarry Smith PetscFunctionBegin; 7925a9f2f41SSatish Balay while (ilink) { 793574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 794fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 795fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 796fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 797fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 798*5f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); } 7995a9f2f41SSatish Balay next = ilink->next; 8005a9f2f41SSatish Balay ilink = next; 8010971522cSBarry Smith } 80205b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 803574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 804574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 805574deadeSBarry Smith } else if (jac->mat) { 806574deadeSBarry Smith jac->mat = PETSC_NULL; 807574deadeSBarry Smith } 80897bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 80968dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8106bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8116bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8126bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8136bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8146bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8156bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8166bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 81763ec74ffSBarry Smith jac->reset = PETSC_TRUE; 818574deadeSBarry Smith PetscFunctionReturn(0); 819574deadeSBarry Smith } 820574deadeSBarry Smith 821574deadeSBarry Smith #undef __FUNCT__ 822574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 823574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 824574deadeSBarry Smith { 825574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 826574deadeSBarry Smith PetscErrorCode ierr; 827574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 828574deadeSBarry Smith 829574deadeSBarry Smith PetscFunctionBegin; 830574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 831574deadeSBarry Smith while (ilink) { 8326bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 833574deadeSBarry Smith next = ilink->next; 834574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 835574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 8365d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 837574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 838574deadeSBarry Smith ilink = next; 839574deadeSBarry Smith } 840574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 841c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8420971522cSBarry Smith PetscFunctionReturn(0); 8430971522cSBarry Smith } 8440971522cSBarry Smith 8450971522cSBarry Smith #undef __FUNCT__ 8460971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8470971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8480971522cSBarry Smith { 8491b9fc7fcSBarry Smith PetscErrorCode ierr; 8506c924f48SJed Brown PetscInt bs; 851bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8529dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8533b224e63SBarry Smith PCCompositeType ctype; 8541b9fc7fcSBarry Smith 8550971522cSBarry Smith PetscFunctionBegin; 8561b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 857acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 85851f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 85951f519a2SBarry Smith if (flg) { 86051f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 86151f519a2SBarry Smith } 862704ba839SBarry Smith 863435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 864c0adfefeSBarry Smith if (stokes) { 865c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 866c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 867c0adfefeSBarry Smith } 868c0adfefeSBarry Smith 8693b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 8703b224e63SBarry Smith if (flg) { 8713b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 8723b224e63SBarry Smith } 873d32f9abdSBarry Smith 874c30613efSMatthew Knepley /* Only setup fields once */ 875c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 876d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 877d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 8786c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 8796c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 880d32f9abdSBarry Smith } 881c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 882c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 883084e4875SJed 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); 884c5d2311dSJed Brown } 8851b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8860971522cSBarry Smith PetscFunctionReturn(0); 8870971522cSBarry Smith } 8880971522cSBarry Smith 8890971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 8900971522cSBarry Smith 8910971522cSBarry Smith EXTERN_C_BEGIN 8920971522cSBarry Smith #undef __FUNCT__ 8930971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 8945d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 8950971522cSBarry Smith { 89697bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8970971522cSBarry Smith PetscErrorCode ierr; 8985a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 89969a612a9SBarry Smith char prefix[128]; 9005d4c12cdSJungho Lee PetscInt i; 9010971522cSBarry Smith 9020971522cSBarry Smith PetscFunctionBegin; 9036c924f48SJed Brown if (jac->splitdefined) { 9046c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9056c924f48SJed Brown PetscFunctionReturn(0); 9066c924f48SJed Brown } 90751f519a2SBarry Smith for (i=0; i<n; i++) { 908e32f2f54SBarry 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); 909e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 91051f519a2SBarry Smith } 911704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 912a04f6461SBarry Smith if (splitname) { 913db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 914a04f6461SBarry Smith } else { 915a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 916a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 917a04f6461SBarry Smith } 918704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 9195a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 9205d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 9215d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 9225a9f2f41SSatish Balay ilink->nfields = n; 9235a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9247adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9251cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9265a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9279005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 92869a612a9SBarry Smith 929a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9305a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9310971522cSBarry Smith 9320971522cSBarry Smith if (!next) { 9335a9f2f41SSatish Balay jac->head = ilink; 93451f519a2SBarry Smith ilink->previous = PETSC_NULL; 9350971522cSBarry Smith } else { 9360971522cSBarry Smith while (next->next) { 9370971522cSBarry Smith next = next->next; 9380971522cSBarry Smith } 9395a9f2f41SSatish Balay next->next = ilink; 94051f519a2SBarry Smith ilink->previous = next; 9410971522cSBarry Smith } 9420971522cSBarry Smith jac->nsplits++; 9430971522cSBarry Smith PetscFunctionReturn(0); 9440971522cSBarry Smith } 9450971522cSBarry Smith EXTERN_C_END 9460971522cSBarry Smith 947e69d4d44SBarry Smith EXTERN_C_BEGIN 948e69d4d44SBarry Smith #undef __FUNCT__ 949e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9507087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 951e69d4d44SBarry Smith { 952e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 953e69d4d44SBarry Smith PetscErrorCode ierr; 954e69d4d44SBarry Smith 955e69d4d44SBarry Smith PetscFunctionBegin; 956519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 957e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 958e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 95913e0d083SBarry Smith if (n) *n = jac->nsplits; 960e69d4d44SBarry Smith PetscFunctionReturn(0); 961e69d4d44SBarry Smith } 962e69d4d44SBarry Smith EXTERN_C_END 9630971522cSBarry Smith 9640971522cSBarry Smith EXTERN_C_BEGIN 9650971522cSBarry Smith #undef __FUNCT__ 96669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 9677087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 9680971522cSBarry Smith { 9690971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9700971522cSBarry Smith PetscErrorCode ierr; 9710971522cSBarry Smith PetscInt cnt = 0; 9725a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 9730971522cSBarry Smith 9740971522cSBarry Smith PetscFunctionBegin; 9755d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 9765a9f2f41SSatish Balay while (ilink) { 9775a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 9785a9f2f41SSatish Balay ilink = ilink->next; 9790971522cSBarry Smith } 9805d480477SMatthew 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); 98113e0d083SBarry Smith if (n) *n = jac->nsplits; 9820971522cSBarry Smith PetscFunctionReturn(0); 9830971522cSBarry Smith } 9840971522cSBarry Smith EXTERN_C_END 9850971522cSBarry Smith 986704ba839SBarry Smith EXTERN_C_BEGIN 987704ba839SBarry Smith #undef __FUNCT__ 988704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 9897087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 990704ba839SBarry Smith { 991704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 992704ba839SBarry Smith PetscErrorCode ierr; 993704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 994704ba839SBarry Smith char prefix[128]; 995704ba839SBarry Smith 996704ba839SBarry Smith PetscFunctionBegin; 9976c924f48SJed Brown if (jac->splitdefined) { 9986c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9996c924f48SJed Brown PetscFunctionReturn(0); 10006c924f48SJed Brown } 100116913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1002a04f6461SBarry Smith if (splitname) { 1003db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1004a04f6461SBarry Smith } else { 1005a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 10065d4c12cdSJungho Lee ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1007a04f6461SBarry Smith } 10081ab39975SBarry Smith ilink->is = is; 1009*5f4ab4e1SJungho Lee ilink->is_col = is; 10101ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1011704ba839SBarry Smith ilink->next = PETSC_NULL; 1012704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10131cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1014704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10159005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1016704ba839SBarry Smith 1017a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1018704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1019704ba839SBarry Smith 1020704ba839SBarry Smith if (!next) { 1021704ba839SBarry Smith jac->head = ilink; 1022704ba839SBarry Smith ilink->previous = PETSC_NULL; 1023704ba839SBarry Smith } else { 1024704ba839SBarry Smith while (next->next) { 1025704ba839SBarry Smith next = next->next; 1026704ba839SBarry Smith } 1027704ba839SBarry Smith next->next = ilink; 1028704ba839SBarry Smith ilink->previous = next; 1029704ba839SBarry Smith } 1030704ba839SBarry Smith jac->nsplits++; 1031704ba839SBarry Smith 1032704ba839SBarry Smith PetscFunctionReturn(0); 1033704ba839SBarry Smith } 1034704ba839SBarry Smith EXTERN_C_END 1035704ba839SBarry Smith 10360971522cSBarry Smith #undef __FUNCT__ 10370971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10380971522cSBarry Smith /*@ 10390971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10400971522cSBarry Smith 1041ad4df100SBarry Smith Logically Collective on PC 10420971522cSBarry Smith 10430971522cSBarry Smith Input Parameters: 10440971522cSBarry Smith + pc - the preconditioner context 1045a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10460971522cSBarry Smith . n - the number of fields in this split 1047db4c96c1SJed Brown - fields - the fields in this split 10480971522cSBarry Smith 10490971522cSBarry Smith Level: intermediate 10500971522cSBarry Smith 1051d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1052d32f9abdSBarry Smith 1053d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1054d32f9abdSBarry 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 1055d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1056d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1057d32f9abdSBarry Smith 1058db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1059db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1060db4c96c1SJed Brown 10615d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 10625d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 10635d4c12cdSJungho Lee available when this routine is called. 10645d4c12cdSJungho Lee 1065d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 10660971522cSBarry Smith 10670971522cSBarry Smith @*/ 10685d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10690971522cSBarry Smith { 10704ac538c5SBarry Smith PetscErrorCode ierr; 10710971522cSBarry Smith 10720971522cSBarry Smith PetscFunctionBegin; 10730700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1074db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1075db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1076db4c96c1SJed Brown PetscValidIntPointer(fields,3); 10775d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 10780971522cSBarry Smith PetscFunctionReturn(0); 10790971522cSBarry Smith } 10800971522cSBarry Smith 10810971522cSBarry Smith #undef __FUNCT__ 1082704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1083704ba839SBarry Smith /*@ 1084704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1085704ba839SBarry Smith 1086ad4df100SBarry Smith Logically Collective on PC 1087704ba839SBarry Smith 1088704ba839SBarry Smith Input Parameters: 1089704ba839SBarry Smith + pc - the preconditioner context 1090a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1091db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1092704ba839SBarry Smith 1093d32f9abdSBarry Smith 1094a6ffb8dbSJed Brown Notes: 1095a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1096a6ffb8dbSJed Brown 1097db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1098db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1099d32f9abdSBarry Smith 1100704ba839SBarry Smith Level: intermediate 1101704ba839SBarry Smith 1102704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1103704ba839SBarry Smith 1104704ba839SBarry Smith @*/ 11057087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1106704ba839SBarry Smith { 11074ac538c5SBarry Smith PetscErrorCode ierr; 1108704ba839SBarry Smith 1109704ba839SBarry Smith PetscFunctionBegin; 11100700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11115d4c12cdSJungho Lee PetscValidCharPointer(splitname,2); 1112db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11134ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1114704ba839SBarry Smith PetscFunctionReturn(0); 1115704ba839SBarry Smith } 1116704ba839SBarry Smith 1117704ba839SBarry Smith #undef __FUNCT__ 111857a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 111957a9adfeSMatthew G Knepley /*@ 112057a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 112157a9adfeSMatthew G Knepley 112257a9adfeSMatthew G Knepley Logically Collective on PC 112357a9adfeSMatthew G Knepley 112457a9adfeSMatthew G Knepley Input Parameters: 112557a9adfeSMatthew G Knepley + pc - the preconditioner context 112657a9adfeSMatthew G Knepley - splitname - name of this split 112757a9adfeSMatthew G Knepley 112857a9adfeSMatthew G Knepley Output Parameter: 112957a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 113057a9adfeSMatthew G Knepley 113157a9adfeSMatthew G Knepley Level: intermediate 113257a9adfeSMatthew G Knepley 113357a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 113457a9adfeSMatthew G Knepley 113557a9adfeSMatthew G Knepley @*/ 113657a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 113757a9adfeSMatthew G Knepley { 113857a9adfeSMatthew G Knepley PetscErrorCode ierr; 113957a9adfeSMatthew G Knepley 114057a9adfeSMatthew G Knepley PetscFunctionBegin; 114157a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 114257a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 114357a9adfeSMatthew G Knepley PetscValidPointer(is,3); 114457a9adfeSMatthew G Knepley { 114557a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 114657a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 114757a9adfeSMatthew G Knepley PetscBool found; 114857a9adfeSMatthew G Knepley 114957a9adfeSMatthew G Knepley *is = PETSC_NULL; 115057a9adfeSMatthew G Knepley while(ilink) { 115157a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 115257a9adfeSMatthew G Knepley if (found) { 115357a9adfeSMatthew G Knepley *is = ilink->is; 115457a9adfeSMatthew G Knepley break; 115557a9adfeSMatthew G Knepley } 115657a9adfeSMatthew G Knepley ilink = ilink->next; 115757a9adfeSMatthew G Knepley } 115857a9adfeSMatthew G Knepley } 115957a9adfeSMatthew G Knepley PetscFunctionReturn(0); 116057a9adfeSMatthew G Knepley } 116157a9adfeSMatthew G Knepley 116257a9adfeSMatthew G Knepley #undef __FUNCT__ 116351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 116451f519a2SBarry Smith /*@ 116551f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 116651f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 116751f519a2SBarry Smith 1168ad4df100SBarry Smith Logically Collective on PC 116951f519a2SBarry Smith 117051f519a2SBarry Smith Input Parameters: 117151f519a2SBarry Smith + pc - the preconditioner context 117251f519a2SBarry Smith - bs - the block size 117351f519a2SBarry Smith 117451f519a2SBarry Smith Level: intermediate 117551f519a2SBarry Smith 117651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 117751f519a2SBarry Smith 117851f519a2SBarry Smith @*/ 11797087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 118051f519a2SBarry Smith { 11814ac538c5SBarry Smith PetscErrorCode ierr; 118251f519a2SBarry Smith 118351f519a2SBarry Smith PetscFunctionBegin; 11840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1185c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 11864ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 118751f519a2SBarry Smith PetscFunctionReturn(0); 118851f519a2SBarry Smith } 118951f519a2SBarry Smith 119051f519a2SBarry Smith #undef __FUNCT__ 119169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 11920971522cSBarry Smith /*@C 119369a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 11940971522cSBarry Smith 119569a612a9SBarry Smith Collective on KSP 11960971522cSBarry Smith 11970971522cSBarry Smith Input Parameter: 11980971522cSBarry Smith . pc - the preconditioner context 11990971522cSBarry Smith 12000971522cSBarry Smith Output Parameters: 120113e0d083SBarry Smith + n - the number of splits 120269a612a9SBarry Smith - pc - the array of KSP contexts 12030971522cSBarry Smith 12040971522cSBarry Smith Note: 1205d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1206d32f9abdSBarry Smith (not the KSP just the array that contains them). 12070971522cSBarry Smith 120869a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12090971522cSBarry Smith 12100971522cSBarry Smith Level: advanced 12110971522cSBarry Smith 12120971522cSBarry Smith .seealso: PCFIELDSPLIT 12130971522cSBarry Smith @*/ 12147087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12150971522cSBarry Smith { 12164ac538c5SBarry Smith PetscErrorCode ierr; 12170971522cSBarry Smith 12180971522cSBarry Smith PetscFunctionBegin; 12190700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 122013e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 12214ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 12220971522cSBarry Smith PetscFunctionReturn(0); 12230971522cSBarry Smith } 12240971522cSBarry Smith 1225e69d4d44SBarry Smith #undef __FUNCT__ 1226e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1227e69d4d44SBarry Smith /*@ 1228e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1229a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1230e69d4d44SBarry Smith 1231e69d4d44SBarry Smith Collective on PC 1232e69d4d44SBarry Smith 1233e69d4d44SBarry Smith Input Parameters: 1234e69d4d44SBarry Smith + pc - the preconditioner context 1235fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1236084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1237084e4875SJed Brown 1238e69d4d44SBarry Smith Options Database: 1239084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1240e69d4d44SBarry Smith 1241fd1303e9SJungho Lee Notes: 1242fd1303e9SJungho Lee $ If ptype is 1243fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1244fd1303e9SJungho Lee $ to this function). 1245fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1246fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1247fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1248fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1249fd1303e9SJungho Lee $ preconditioner 1250fd1303e9SJungho Lee 1251fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1252fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1253fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1254fd1303e9SJungho Lee 1255fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1256fd1303e9SJungho Lee the name since it sets a proceedure to use. 1257fd1303e9SJungho Lee 1258e69d4d44SBarry Smith Level: intermediate 1259e69d4d44SBarry Smith 1260fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1261e69d4d44SBarry Smith 1262e69d4d44SBarry Smith @*/ 12637087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1264e69d4d44SBarry Smith { 12654ac538c5SBarry Smith PetscErrorCode ierr; 1266e69d4d44SBarry Smith 1267e69d4d44SBarry Smith PetscFunctionBegin; 12680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12694ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1270e69d4d44SBarry Smith PetscFunctionReturn(0); 1271e69d4d44SBarry Smith } 1272e69d4d44SBarry Smith 1273e69d4d44SBarry Smith EXTERN_C_BEGIN 1274e69d4d44SBarry Smith #undef __FUNCT__ 1275e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 12767087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1277e69d4d44SBarry Smith { 1278e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1279084e4875SJed Brown PetscErrorCode ierr; 1280e69d4d44SBarry Smith 1281e69d4d44SBarry Smith PetscFunctionBegin; 1282084e4875SJed Brown jac->schurpre = ptype; 1283084e4875SJed Brown if (pre) { 12846bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1285084e4875SJed Brown jac->schur_user = pre; 1286084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1287084e4875SJed Brown } 1288e69d4d44SBarry Smith PetscFunctionReturn(0); 1289e69d4d44SBarry Smith } 1290e69d4d44SBarry Smith EXTERN_C_END 1291e69d4d44SBarry Smith 129230ad9308SMatthew Knepley #undef __FUNCT__ 129330ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 129430ad9308SMatthew Knepley /*@C 12958c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 129630ad9308SMatthew Knepley 129730ad9308SMatthew Knepley Collective on KSP 129830ad9308SMatthew Knepley 129930ad9308SMatthew Knepley Input Parameter: 130030ad9308SMatthew Knepley . pc - the preconditioner context 130130ad9308SMatthew Knepley 130230ad9308SMatthew Knepley Output Parameters: 1303a04f6461SBarry Smith + A00 - the (0,0) block 1304a04f6461SBarry Smith . A01 - the (0,1) block 1305a04f6461SBarry Smith . A10 - the (1,0) block 1306a04f6461SBarry Smith - A11 - the (1,1) block 130730ad9308SMatthew Knepley 130830ad9308SMatthew Knepley Level: advanced 130930ad9308SMatthew Knepley 131030ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 131130ad9308SMatthew Knepley @*/ 1312a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 131330ad9308SMatthew Knepley { 131430ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 131530ad9308SMatthew Knepley 131630ad9308SMatthew Knepley PetscFunctionBegin; 13170700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1318c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1319a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1320a04f6461SBarry Smith if (A01) *A01 = jac->B; 1321a04f6461SBarry Smith if (A10) *A10 = jac->C; 1322a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 132330ad9308SMatthew Knepley PetscFunctionReturn(0); 132430ad9308SMatthew Knepley } 132530ad9308SMatthew Knepley 132679416396SBarry Smith EXTERN_C_BEGIN 132779416396SBarry Smith #undef __FUNCT__ 132879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 13297087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 133079416396SBarry Smith { 133179416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1332e69d4d44SBarry Smith PetscErrorCode ierr; 133379416396SBarry Smith 133479416396SBarry Smith PetscFunctionBegin; 133579416396SBarry Smith jac->type = type; 13363b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 13373b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 13383b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1339e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1340e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1341e69d4d44SBarry Smith 13423b224e63SBarry Smith } else { 13433b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 13443b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1345e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13469e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 13473b224e63SBarry Smith } 134879416396SBarry Smith PetscFunctionReturn(0); 134979416396SBarry Smith } 135079416396SBarry Smith EXTERN_C_END 135179416396SBarry Smith 135251f519a2SBarry Smith EXTERN_C_BEGIN 135351f519a2SBarry Smith #undef __FUNCT__ 135451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 13557087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 135651f519a2SBarry Smith { 135751f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 135851f519a2SBarry Smith 135951f519a2SBarry Smith PetscFunctionBegin; 136065e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 136165e19b50SBarry 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); 136251f519a2SBarry Smith jac->bs = bs; 136351f519a2SBarry Smith PetscFunctionReturn(0); 136451f519a2SBarry Smith } 136551f519a2SBarry Smith EXTERN_C_END 136651f519a2SBarry Smith 136779416396SBarry Smith #undef __FUNCT__ 136879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1369bc08b0f1SBarry Smith /*@ 137079416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 137179416396SBarry Smith 137279416396SBarry Smith Collective on PC 137379416396SBarry Smith 137479416396SBarry Smith Input Parameter: 137579416396SBarry Smith . pc - the preconditioner context 137681540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 137779416396SBarry Smith 137879416396SBarry Smith Options Database Key: 1379a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 138079416396SBarry Smith 138179416396SBarry Smith Level: Developer 138279416396SBarry Smith 138379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 138479416396SBarry Smith 138579416396SBarry Smith .seealso: PCCompositeSetType() 138679416396SBarry Smith 138779416396SBarry Smith @*/ 13887087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 138979416396SBarry Smith { 13904ac538c5SBarry Smith PetscErrorCode ierr; 139179416396SBarry Smith 139279416396SBarry Smith PetscFunctionBegin; 13930700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13944ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 139579416396SBarry Smith PetscFunctionReturn(0); 139679416396SBarry Smith } 139779416396SBarry Smith 13980971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 13990971522cSBarry Smith /*MC 1400a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1401a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 14020971522cSBarry Smith 1403edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1404edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 14050971522cSBarry Smith 1406a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 140769a612a9SBarry Smith and set the options directly on the resulting KSP object 14080971522cSBarry Smith 14090971522cSBarry Smith Level: intermediate 14100971522cSBarry Smith 141179416396SBarry Smith Options Database Keys: 141281540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 141381540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 141481540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 141581540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 141681540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1417e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1418435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 141979416396SBarry Smith 14205d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 14215d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 14225d4c12cdSJungho Lee 1423c8a0d604SMatthew G Knepley Notes: 1424c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1425d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1426d32f9abdSBarry Smith 1427d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1428d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1429d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1430d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1431d32f9abdSBarry Smith 1432c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1433c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1434c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1435c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1436c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1437a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1438c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1439c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1440c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1441a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1442c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1443c8a0d604SMatthew G Knepley factorization type is set using -pc_fieldsplit_schur_factorization_type <diag, lower, upper, full>. The full is shown above, 1444c8a0d604SMatthew G Knepley but diag gives 1445c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1446c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 1447c8a0d604SMatthew G Knepley so that the preconditioner is positive definite. The lower factorization is the inverse of 1448c8a0d604SMatthew G Knepley $ ( A00 0 ) 1449c8a0d604SMatthew G Knepley $ ( A10 S ) 1450c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1451c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1452c8a0d604SMatthew G Knepley $ ( 0 S ) 1453c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1454e69d4d44SBarry Smith 1455edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1456edf189efSBarry Smith is used automatically for a second block. 1457edf189efSBarry Smith 1458ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1459ff218e97SBarry Smith Generally it should be used with the AIJ format. 1460ff218e97SBarry Smith 1461ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1462ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1463ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 14640716a85fSBarry Smith 1465a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 14660971522cSBarry Smith 14677e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1468e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 14690971522cSBarry Smith M*/ 14700971522cSBarry Smith 14710971522cSBarry Smith EXTERN_C_BEGIN 14720971522cSBarry Smith #undef __FUNCT__ 14730971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 14747087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 14750971522cSBarry Smith { 14760971522cSBarry Smith PetscErrorCode ierr; 14770971522cSBarry Smith PC_FieldSplit *jac; 14780971522cSBarry Smith 14790971522cSBarry Smith PetscFunctionBegin; 148038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 14810971522cSBarry Smith jac->bs = -1; 14820971522cSBarry Smith jac->nsplits = 0; 14833e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1484e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1485c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 148651f519a2SBarry Smith 14870971522cSBarry Smith pc->data = (void*)jac; 14880971522cSBarry Smith 14890971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1490421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 14910971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1492574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 14930971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 14940971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 14950971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 14960971522cSBarry Smith pc->ops->applyrichardson = 0; 14970971522cSBarry Smith 149869a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 149969a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 15000971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 15010971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1502704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1503704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 150479416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 150579416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 150651f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 150751f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 15080971522cSBarry Smith PetscFunctionReturn(0); 15090971522cSBarry Smith } 15100971522cSBarry Smith EXTERN_C_END 15110971522cSBarry Smith 15120971522cSBarry Smith 1513a541d17aSBarry Smith 1514