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; 3445f4ab4e1SJungho 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); 3655f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 366704ba839SBarry Smith } else if (!ilink->is) { 367ccb205f8SBarry Smith if (ilink->nfields > 1) { 3685f4ab4e1SJungho 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); 3705f4ab4e1SJungho 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]; 3745f4ab4e1SJungho 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); 3785f4ab4e1SJungho 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); 3815f4ab4e1SJungho 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); 3855f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 3865f4ab4e1SJungho 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++) { 3955f4ab4e1SJungho 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++) { 4005f4ab4e1SJungho 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++) { 4095f4ab4e1SJungho 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++) { 4145f4ab4e1SJungho 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; 442*093c86ffSJed Brown char lscname[256]; 443*093c86ffSJed Brown PetscObject LSC_L; 444e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 44568dd23aaSBarry Smith 446e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 447e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 448e6cab6aaSJed Brown 4493b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4503b224e63SBarry Smith if (jac->schur) { 4513b224e63SBarry Smith ilink = jac->head; 45249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4534aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 454fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4553b224e63SBarry Smith ilink = ilink->next; 45649bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4574aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 458fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 459519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 460084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4613b224e63SBarry Smith 4623b224e63SBarry Smith } else { 4631cee3971SBarry Smith KSP ksp; 4646c924f48SJed Brown char schurprefix[256]; 4653b224e63SBarry Smith 466a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4673b224e63SBarry Smith ilink = jac->head; 46849bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4694aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 470fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4713b224e63SBarry Smith ilink = ilink->next; 47249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 4734aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 474fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4757d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 476519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 477a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4781cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 479e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 480a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 481a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 48220b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 48320b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4843b224e63SBarry Smith 4853b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4869005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4871cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 488084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 489084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 490084e4875SJed Brown PC pc; 491cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 492084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 493084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 494e69d4d44SBarry Smith } 4956c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 4966c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 4973b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 49820b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 49920b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 5003b224e63SBarry Smith 5013b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 5023b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 5033b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 5043b224e63SBarry Smith ilink = jac->head; 5053b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 5063b224e63SBarry Smith ilink = ilink->next; 5073b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 5083b224e63SBarry Smith } 509*093c86ffSJed Brown 510*093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 511*093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 512*093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 513*093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 514*093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 515*093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 516*093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 517*093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 518*093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 5193b224e63SBarry Smith } else { 52097bbdb24SBarry Smith /* set up the individual PCs */ 52197bbdb24SBarry Smith i = 0; 5225a9f2f41SSatish Balay ilink = jac->head; 5235a9f2f41SSatish Balay while (ilink) { 524519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 5253b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 526c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 5275a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 52897bbdb24SBarry Smith i++; 5295a9f2f41SSatish Balay ilink = ilink->next; 5300971522cSBarry Smith } 53197bbdb24SBarry Smith 53297bbdb24SBarry Smith /* create work vectors for each split */ 5331b9fc7fcSBarry Smith if (!jac->x) { 53497bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5355a9f2f41SSatish Balay ilink = jac->head; 53697bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 537906ed7ccSBarry Smith Vec *vl,*vr; 538906ed7ccSBarry Smith 539906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 540906ed7ccSBarry Smith ilink->x = *vr; 541906ed7ccSBarry Smith ilink->y = *vl; 542906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 543906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5445a9f2f41SSatish Balay jac->x[i] = ilink->x; 5455a9f2f41SSatish Balay jac->y[i] = ilink->y; 5465a9f2f41SSatish Balay ilink = ilink->next; 54797bbdb24SBarry Smith } 5483b224e63SBarry Smith } 5493b224e63SBarry Smith } 5503b224e63SBarry Smith 5513b224e63SBarry Smith 552d0f46423SBarry Smith if (!jac->head->sctx) { 5533b224e63SBarry Smith Vec xtmp; 5543b224e63SBarry Smith 55579416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5561b9fc7fcSBarry Smith 5575a9f2f41SSatish Balay ilink = jac->head; 5581b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5591b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 560704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5615a9f2f41SSatish Balay ilink = ilink->next; 5621b9fc7fcSBarry Smith } 5636bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5641b9fc7fcSBarry Smith } 565c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 5660971522cSBarry Smith PetscFunctionReturn(0); 5670971522cSBarry Smith } 5680971522cSBarry Smith 5695a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 570ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 571ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5725a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 573ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 574ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 57579416396SBarry Smith 5760971522cSBarry Smith #undef __FUNCT__ 5773b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5783b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5793b224e63SBarry Smith { 5803b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5813b224e63SBarry Smith PetscErrorCode ierr; 5823b224e63SBarry Smith KSP ksp; 5833b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5843b224e63SBarry Smith 5853b224e63SBarry Smith PetscFunctionBegin; 5863b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5873b224e63SBarry Smith 588c5d2311dSJed Brown switch (jac->schurfactorization) { 589c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 590a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 591c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 592c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 593c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 594c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 595c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 596c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 598c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 599c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 600c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 601c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 602c5d2311dSJed Brown break; 603c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 604a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 605c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 606c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 607c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 608c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 609c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 610c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 611c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 612c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 613c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 614c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 615c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 616c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 617c5d2311dSJed Brown break; 618c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 619a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 620c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 621c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 622c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 623c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 624c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 625c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 626c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 627c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 628c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 629c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 630c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 631c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 632c5d2311dSJed Brown break; 633c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 634a04f6461SBarry 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 */ 6353b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6363b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6373b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6383b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6393b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6403b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6413b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6423b224e63SBarry Smith 6433b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6443b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6453b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6463b224e63SBarry Smith 6473b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6483b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6493b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6503b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6513b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 652c5d2311dSJed Brown } 6533b224e63SBarry Smith PetscFunctionReturn(0); 6543b224e63SBarry Smith } 6553b224e63SBarry Smith 6563b224e63SBarry Smith #undef __FUNCT__ 6570971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6580971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6590971522cSBarry Smith { 6600971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6610971522cSBarry Smith PetscErrorCode ierr; 6625a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 663af93645bSJed Brown PetscInt cnt; 6640971522cSBarry Smith 6650971522cSBarry Smith PetscFunctionBegin; 66651f519a2SBarry Smith CHKMEMQ; 66751f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 66851f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 66951f519a2SBarry Smith 67079416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6711b9fc7fcSBarry Smith if (jac->defaultsplit) { 6720971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6735a9f2f41SSatish Balay while (ilink) { 6745a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6755a9f2f41SSatish Balay ilink = ilink->next; 6760971522cSBarry Smith } 6770971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6781b9fc7fcSBarry Smith } else { 679efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6805a9f2f41SSatish Balay while (ilink) { 6815a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6825a9f2f41SSatish Balay ilink = ilink->next; 6831b9fc7fcSBarry Smith } 6841b9fc7fcSBarry Smith } 68516913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 68679416396SBarry Smith if (!jac->w1) { 68779416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 68879416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 68979416396SBarry Smith } 690efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6915a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6923e197d65SBarry Smith cnt = 1; 6935a9f2f41SSatish Balay while (ilink->next) { 6945a9f2f41SSatish Balay ilink = ilink->next; 6953e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6963e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6973e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6983e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6993e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7003e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7013e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7023e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7033e197d65SBarry Smith } 70451f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 70511755939SBarry Smith cnt -= 2; 70651f519a2SBarry Smith while (ilink->previous) { 70751f519a2SBarry Smith ilink = ilink->previous; 70811755939SBarry Smith /* compute the residual only over the part of the vector needed */ 70911755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 71011755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 71111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 71211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 71311755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 71411755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 71511755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 71651f519a2SBarry Smith } 71711755939SBarry Smith } 71865e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 71951f519a2SBarry Smith CHKMEMQ; 7200971522cSBarry Smith PetscFunctionReturn(0); 7210971522cSBarry Smith } 7220971522cSBarry Smith 723421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 724ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 725ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 726421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 727ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 728ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 729421e10b8SBarry Smith 730421e10b8SBarry Smith #undef __FUNCT__ 7318c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 732421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 733421e10b8SBarry Smith { 734421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 735421e10b8SBarry Smith PetscErrorCode ierr; 736421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 737421e10b8SBarry Smith 738421e10b8SBarry Smith PetscFunctionBegin; 739421e10b8SBarry Smith CHKMEMQ; 740421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 741421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 742421e10b8SBarry Smith 743421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 744421e10b8SBarry Smith if (jac->defaultsplit) { 745421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 746421e10b8SBarry Smith while (ilink) { 747421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 748421e10b8SBarry Smith ilink = ilink->next; 749421e10b8SBarry Smith } 750421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 751421e10b8SBarry Smith } else { 752421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 753421e10b8SBarry Smith while (ilink) { 754421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 755421e10b8SBarry Smith ilink = ilink->next; 756421e10b8SBarry Smith } 757421e10b8SBarry Smith } 758421e10b8SBarry Smith } else { 759421e10b8SBarry Smith if (!jac->w1) { 760421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 761421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 762421e10b8SBarry Smith } 763421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 764421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 765421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 766421e10b8SBarry Smith while (ilink->next) { 767421e10b8SBarry Smith ilink = ilink->next; 7689989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 769421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 770421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 771421e10b8SBarry Smith } 772421e10b8SBarry Smith while (ilink->previous) { 773421e10b8SBarry Smith ilink = ilink->previous; 7749989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 775421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 776421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 777421e10b8SBarry Smith } 778421e10b8SBarry Smith } else { 779421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 780421e10b8SBarry Smith ilink = ilink->next; 781421e10b8SBarry Smith } 782421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 783421e10b8SBarry Smith while (ilink->previous) { 784421e10b8SBarry Smith ilink = ilink->previous; 7859989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 786421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 787421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 788421e10b8SBarry Smith } 789421e10b8SBarry Smith } 790421e10b8SBarry Smith } 791421e10b8SBarry Smith CHKMEMQ; 792421e10b8SBarry Smith PetscFunctionReturn(0); 793421e10b8SBarry Smith } 794421e10b8SBarry Smith 7950971522cSBarry Smith #undef __FUNCT__ 796574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 797574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 7980971522cSBarry Smith { 7990971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8000971522cSBarry Smith PetscErrorCode ierr; 8015a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 8020971522cSBarry Smith 8030971522cSBarry Smith PetscFunctionBegin; 8045a9f2f41SSatish Balay while (ilink) { 805574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 806fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 807fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 808fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 809fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 8105f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); } 8115a9f2f41SSatish Balay next = ilink->next; 8125a9f2f41SSatish Balay ilink = next; 8130971522cSBarry Smith } 81405b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 815574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 816574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 817574deadeSBarry Smith } else if (jac->mat) { 818574deadeSBarry Smith jac->mat = PETSC_NULL; 819574deadeSBarry Smith } 82097bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 82168dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8226bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8236bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8246bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8256bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8266bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8276bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8286bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 82963ec74ffSBarry Smith jac->reset = PETSC_TRUE; 830574deadeSBarry Smith PetscFunctionReturn(0); 831574deadeSBarry Smith } 832574deadeSBarry Smith 833574deadeSBarry Smith #undef __FUNCT__ 834574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 835574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 836574deadeSBarry Smith { 837574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 838574deadeSBarry Smith PetscErrorCode ierr; 839574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 840574deadeSBarry Smith 841574deadeSBarry Smith PetscFunctionBegin; 842574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 843574deadeSBarry Smith while (ilink) { 8446bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 845574deadeSBarry Smith next = ilink->next; 846574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 847574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 8485d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 849574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 850574deadeSBarry Smith ilink = next; 851574deadeSBarry Smith } 852574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 853c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8540971522cSBarry Smith PetscFunctionReturn(0); 8550971522cSBarry Smith } 8560971522cSBarry Smith 8570971522cSBarry Smith #undef __FUNCT__ 8580971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8590971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8600971522cSBarry Smith { 8611b9fc7fcSBarry Smith PetscErrorCode ierr; 8626c924f48SJed Brown PetscInt bs; 863bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8649dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8653b224e63SBarry Smith PCCompositeType ctype; 8661b9fc7fcSBarry Smith 8670971522cSBarry Smith PetscFunctionBegin; 8681b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 869acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 87051f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 87151f519a2SBarry Smith if (flg) { 87251f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 87351f519a2SBarry Smith } 874704ba839SBarry Smith 875435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 876c0adfefeSBarry Smith if (stokes) { 877c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 878c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 879c0adfefeSBarry Smith } 880c0adfefeSBarry Smith 8813b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 8823b224e63SBarry Smith if (flg) { 8833b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 8843b224e63SBarry Smith } 885d32f9abdSBarry Smith 886c30613efSMatthew Knepley /* Only setup fields once */ 887c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 888d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 889d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 8906c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 8916c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 892d32f9abdSBarry Smith } 893c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 894c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 895084e4875SJed 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); 896c5d2311dSJed Brown } 8971b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8980971522cSBarry Smith PetscFunctionReturn(0); 8990971522cSBarry Smith } 9000971522cSBarry Smith 9010971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 9020971522cSBarry Smith 9030971522cSBarry Smith EXTERN_C_BEGIN 9040971522cSBarry Smith #undef __FUNCT__ 9050971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 9065d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 9070971522cSBarry Smith { 90897bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9090971522cSBarry Smith PetscErrorCode ierr; 9105a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 91169a612a9SBarry Smith char prefix[128]; 9125d4c12cdSJungho Lee PetscInt i; 9130971522cSBarry Smith 9140971522cSBarry Smith PetscFunctionBegin; 9156c924f48SJed Brown if (jac->splitdefined) { 9166c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9176c924f48SJed Brown PetscFunctionReturn(0); 9186c924f48SJed Brown } 91951f519a2SBarry Smith for (i=0; i<n; i++) { 920e32f2f54SBarry 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); 921e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 92251f519a2SBarry Smith } 923704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 924a04f6461SBarry Smith if (splitname) { 925db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 926a04f6461SBarry Smith } else { 927a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 928a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 929a04f6461SBarry Smith } 930704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 9315a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 9325d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 9335d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 9345a9f2f41SSatish Balay ilink->nfields = n; 9355a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9367adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9371cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9385a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9399005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 94069a612a9SBarry Smith 941a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9425a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9430971522cSBarry Smith 9440971522cSBarry Smith if (!next) { 9455a9f2f41SSatish Balay jac->head = ilink; 94651f519a2SBarry Smith ilink->previous = PETSC_NULL; 9470971522cSBarry Smith } else { 9480971522cSBarry Smith while (next->next) { 9490971522cSBarry Smith next = next->next; 9500971522cSBarry Smith } 9515a9f2f41SSatish Balay next->next = ilink; 95251f519a2SBarry Smith ilink->previous = next; 9530971522cSBarry Smith } 9540971522cSBarry Smith jac->nsplits++; 9550971522cSBarry Smith PetscFunctionReturn(0); 9560971522cSBarry Smith } 9570971522cSBarry Smith EXTERN_C_END 9580971522cSBarry Smith 959e69d4d44SBarry Smith EXTERN_C_BEGIN 960e69d4d44SBarry Smith #undef __FUNCT__ 961e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9627087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 963e69d4d44SBarry Smith { 964e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 965e69d4d44SBarry Smith PetscErrorCode ierr; 966e69d4d44SBarry Smith 967e69d4d44SBarry Smith PetscFunctionBegin; 968519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 969e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 970e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 97113e0d083SBarry Smith if (n) *n = jac->nsplits; 972e69d4d44SBarry Smith PetscFunctionReturn(0); 973e69d4d44SBarry Smith } 974e69d4d44SBarry Smith EXTERN_C_END 9750971522cSBarry Smith 9760971522cSBarry Smith EXTERN_C_BEGIN 9770971522cSBarry Smith #undef __FUNCT__ 97869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 9797087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 9800971522cSBarry Smith { 9810971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9820971522cSBarry Smith PetscErrorCode ierr; 9830971522cSBarry Smith PetscInt cnt = 0; 9845a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 9850971522cSBarry Smith 9860971522cSBarry Smith PetscFunctionBegin; 9875d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 9885a9f2f41SSatish Balay while (ilink) { 9895a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 9905a9f2f41SSatish Balay ilink = ilink->next; 9910971522cSBarry Smith } 9925d480477SMatthew 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); 99313e0d083SBarry Smith if (n) *n = jac->nsplits; 9940971522cSBarry Smith PetscFunctionReturn(0); 9950971522cSBarry Smith } 9960971522cSBarry Smith EXTERN_C_END 9970971522cSBarry Smith 998704ba839SBarry Smith EXTERN_C_BEGIN 999704ba839SBarry Smith #undef __FUNCT__ 1000704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 10017087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1002704ba839SBarry Smith { 1003704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1004704ba839SBarry Smith PetscErrorCode ierr; 1005704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1006704ba839SBarry Smith char prefix[128]; 1007704ba839SBarry Smith 1008704ba839SBarry Smith PetscFunctionBegin; 10096c924f48SJed Brown if (jac->splitdefined) { 10106c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10116c924f48SJed Brown PetscFunctionReturn(0); 10126c924f48SJed Brown } 101316913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1014a04f6461SBarry Smith if (splitname) { 1015db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1016a04f6461SBarry Smith } else { 1017a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 10185d4c12cdSJungho Lee ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1019a04f6461SBarry Smith } 10201ab39975SBarry Smith ilink->is = is; 10215f4ab4e1SJungho Lee ilink->is_col = is; 10221ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1023704ba839SBarry Smith ilink->next = PETSC_NULL; 1024704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10251cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1026704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10279005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1028704ba839SBarry Smith 1029a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1030704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1031704ba839SBarry Smith 1032704ba839SBarry Smith if (!next) { 1033704ba839SBarry Smith jac->head = ilink; 1034704ba839SBarry Smith ilink->previous = PETSC_NULL; 1035704ba839SBarry Smith } else { 1036704ba839SBarry Smith while (next->next) { 1037704ba839SBarry Smith next = next->next; 1038704ba839SBarry Smith } 1039704ba839SBarry Smith next->next = ilink; 1040704ba839SBarry Smith ilink->previous = next; 1041704ba839SBarry Smith } 1042704ba839SBarry Smith jac->nsplits++; 1043704ba839SBarry Smith 1044704ba839SBarry Smith PetscFunctionReturn(0); 1045704ba839SBarry Smith } 1046704ba839SBarry Smith EXTERN_C_END 1047704ba839SBarry Smith 10480971522cSBarry Smith #undef __FUNCT__ 10490971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10500971522cSBarry Smith /*@ 10510971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10520971522cSBarry Smith 1053ad4df100SBarry Smith Logically Collective on PC 10540971522cSBarry Smith 10550971522cSBarry Smith Input Parameters: 10560971522cSBarry Smith + pc - the preconditioner context 1057a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10580971522cSBarry Smith . n - the number of fields in this split 1059db4c96c1SJed Brown - fields - the fields in this split 10600971522cSBarry Smith 10610971522cSBarry Smith Level: intermediate 10620971522cSBarry Smith 1063d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1064d32f9abdSBarry Smith 1065d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1066d32f9abdSBarry 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 1067d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1068d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1069d32f9abdSBarry Smith 1070db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1071db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1072db4c96c1SJed Brown 10735d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 10745d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 10755d4c12cdSJungho Lee available when this routine is called. 10765d4c12cdSJungho Lee 1077d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 10780971522cSBarry Smith 10790971522cSBarry Smith @*/ 10805d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10810971522cSBarry Smith { 10824ac538c5SBarry Smith PetscErrorCode ierr; 10830971522cSBarry Smith 10840971522cSBarry Smith PetscFunctionBegin; 10850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1086db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1087db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1088db4c96c1SJed Brown PetscValidIntPointer(fields,3); 10895d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 10900971522cSBarry Smith PetscFunctionReturn(0); 10910971522cSBarry Smith } 10920971522cSBarry Smith 10930971522cSBarry Smith #undef __FUNCT__ 1094704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1095704ba839SBarry Smith /*@ 1096704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1097704ba839SBarry Smith 1098ad4df100SBarry Smith Logically Collective on PC 1099704ba839SBarry Smith 1100704ba839SBarry Smith Input Parameters: 1101704ba839SBarry Smith + pc - the preconditioner context 1102a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1103db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1104704ba839SBarry Smith 1105d32f9abdSBarry Smith 1106a6ffb8dbSJed Brown Notes: 1107a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1108a6ffb8dbSJed Brown 1109db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1110db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1111d32f9abdSBarry Smith 1112704ba839SBarry Smith Level: intermediate 1113704ba839SBarry Smith 1114704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1115704ba839SBarry Smith 1116704ba839SBarry Smith @*/ 11177087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1118704ba839SBarry Smith { 11194ac538c5SBarry Smith PetscErrorCode ierr; 1120704ba839SBarry Smith 1121704ba839SBarry Smith PetscFunctionBegin; 11220700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11235d4c12cdSJungho Lee PetscValidCharPointer(splitname,2); 1124db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11254ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1126704ba839SBarry Smith PetscFunctionReturn(0); 1127704ba839SBarry Smith } 1128704ba839SBarry Smith 1129704ba839SBarry Smith #undef __FUNCT__ 113057a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 113157a9adfeSMatthew G Knepley /*@ 113257a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 113357a9adfeSMatthew G Knepley 113457a9adfeSMatthew G Knepley Logically Collective on PC 113557a9adfeSMatthew G Knepley 113657a9adfeSMatthew G Knepley Input Parameters: 113757a9adfeSMatthew G Knepley + pc - the preconditioner context 113857a9adfeSMatthew G Knepley - splitname - name of this split 113957a9adfeSMatthew G Knepley 114057a9adfeSMatthew G Knepley Output Parameter: 114157a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 114257a9adfeSMatthew G Knepley 114357a9adfeSMatthew G Knepley Level: intermediate 114457a9adfeSMatthew G Knepley 114557a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 114657a9adfeSMatthew G Knepley 114757a9adfeSMatthew G Knepley @*/ 114857a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 114957a9adfeSMatthew G Knepley { 115057a9adfeSMatthew G Knepley PetscErrorCode ierr; 115157a9adfeSMatthew G Knepley 115257a9adfeSMatthew G Knepley PetscFunctionBegin; 115357a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 115457a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 115557a9adfeSMatthew G Knepley PetscValidPointer(is,3); 115657a9adfeSMatthew G Knepley { 115757a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 115857a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 115957a9adfeSMatthew G Knepley PetscBool found; 116057a9adfeSMatthew G Knepley 116157a9adfeSMatthew G Knepley *is = PETSC_NULL; 116257a9adfeSMatthew G Knepley while(ilink) { 116357a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 116457a9adfeSMatthew G Knepley if (found) { 116557a9adfeSMatthew G Knepley *is = ilink->is; 116657a9adfeSMatthew G Knepley break; 116757a9adfeSMatthew G Knepley } 116857a9adfeSMatthew G Knepley ilink = ilink->next; 116957a9adfeSMatthew G Knepley } 117057a9adfeSMatthew G Knepley } 117157a9adfeSMatthew G Knepley PetscFunctionReturn(0); 117257a9adfeSMatthew G Knepley } 117357a9adfeSMatthew G Knepley 117457a9adfeSMatthew G Knepley #undef __FUNCT__ 117551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 117651f519a2SBarry Smith /*@ 117751f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 117851f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 117951f519a2SBarry Smith 1180ad4df100SBarry Smith Logically Collective on PC 118151f519a2SBarry Smith 118251f519a2SBarry Smith Input Parameters: 118351f519a2SBarry Smith + pc - the preconditioner context 118451f519a2SBarry Smith - bs - the block size 118551f519a2SBarry Smith 118651f519a2SBarry Smith Level: intermediate 118751f519a2SBarry Smith 118851f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 118951f519a2SBarry Smith 119051f519a2SBarry Smith @*/ 11917087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 119251f519a2SBarry Smith { 11934ac538c5SBarry Smith PetscErrorCode ierr; 119451f519a2SBarry Smith 119551f519a2SBarry Smith PetscFunctionBegin; 11960700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1197c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 11984ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 119951f519a2SBarry Smith PetscFunctionReturn(0); 120051f519a2SBarry Smith } 120151f519a2SBarry Smith 120251f519a2SBarry Smith #undef __FUNCT__ 120369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 12040971522cSBarry Smith /*@C 120569a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 12060971522cSBarry Smith 120769a612a9SBarry Smith Collective on KSP 12080971522cSBarry Smith 12090971522cSBarry Smith Input Parameter: 12100971522cSBarry Smith . pc - the preconditioner context 12110971522cSBarry Smith 12120971522cSBarry Smith Output Parameters: 121313e0d083SBarry Smith + n - the number of splits 121469a612a9SBarry Smith - pc - the array of KSP contexts 12150971522cSBarry Smith 12160971522cSBarry Smith Note: 1217d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1218d32f9abdSBarry Smith (not the KSP just the array that contains them). 12190971522cSBarry Smith 122069a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12210971522cSBarry Smith 12220971522cSBarry Smith Level: advanced 12230971522cSBarry Smith 12240971522cSBarry Smith .seealso: PCFIELDSPLIT 12250971522cSBarry Smith @*/ 12267087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12270971522cSBarry Smith { 12284ac538c5SBarry Smith PetscErrorCode ierr; 12290971522cSBarry Smith 12300971522cSBarry Smith PetscFunctionBegin; 12310700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 123213e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 12334ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 12340971522cSBarry Smith PetscFunctionReturn(0); 12350971522cSBarry Smith } 12360971522cSBarry Smith 1237e69d4d44SBarry Smith #undef __FUNCT__ 1238e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1239e69d4d44SBarry Smith /*@ 1240e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1241a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1242e69d4d44SBarry Smith 1243e69d4d44SBarry Smith Collective on PC 1244e69d4d44SBarry Smith 1245e69d4d44SBarry Smith Input Parameters: 1246e69d4d44SBarry Smith + pc - the preconditioner context 1247fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1248084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1249084e4875SJed Brown 1250e69d4d44SBarry Smith Options Database: 1251084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1252e69d4d44SBarry Smith 1253fd1303e9SJungho Lee Notes: 1254fd1303e9SJungho Lee $ If ptype is 1255fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1256fd1303e9SJungho Lee $ to this function). 1257fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1258fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1259fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1260fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1261fd1303e9SJungho Lee $ preconditioner 1262fd1303e9SJungho Lee 1263fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1264fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1265fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1266fd1303e9SJungho Lee 1267fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1268fd1303e9SJungho Lee the name since it sets a proceedure to use. 1269fd1303e9SJungho Lee 1270e69d4d44SBarry Smith Level: intermediate 1271e69d4d44SBarry Smith 1272fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1273e69d4d44SBarry Smith 1274e69d4d44SBarry Smith @*/ 12757087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1276e69d4d44SBarry Smith { 12774ac538c5SBarry Smith PetscErrorCode ierr; 1278e69d4d44SBarry Smith 1279e69d4d44SBarry Smith PetscFunctionBegin; 12800700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12814ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1282e69d4d44SBarry Smith PetscFunctionReturn(0); 1283e69d4d44SBarry Smith } 1284e69d4d44SBarry Smith 1285e69d4d44SBarry Smith EXTERN_C_BEGIN 1286e69d4d44SBarry Smith #undef __FUNCT__ 1287e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 12887087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1289e69d4d44SBarry Smith { 1290e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1291084e4875SJed Brown PetscErrorCode ierr; 1292e69d4d44SBarry Smith 1293e69d4d44SBarry Smith PetscFunctionBegin; 1294084e4875SJed Brown jac->schurpre = ptype; 1295084e4875SJed Brown if (pre) { 12966bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1297084e4875SJed Brown jac->schur_user = pre; 1298084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1299084e4875SJed Brown } 1300e69d4d44SBarry Smith PetscFunctionReturn(0); 1301e69d4d44SBarry Smith } 1302e69d4d44SBarry Smith EXTERN_C_END 1303e69d4d44SBarry Smith 130430ad9308SMatthew Knepley #undef __FUNCT__ 130530ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 130630ad9308SMatthew Knepley /*@C 13078c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 130830ad9308SMatthew Knepley 130930ad9308SMatthew Knepley Collective on KSP 131030ad9308SMatthew Knepley 131130ad9308SMatthew Knepley Input Parameter: 131230ad9308SMatthew Knepley . pc - the preconditioner context 131330ad9308SMatthew Knepley 131430ad9308SMatthew Knepley Output Parameters: 1315a04f6461SBarry Smith + A00 - the (0,0) block 1316a04f6461SBarry Smith . A01 - the (0,1) block 1317a04f6461SBarry Smith . A10 - the (1,0) block 1318a04f6461SBarry Smith - A11 - the (1,1) block 131930ad9308SMatthew Knepley 132030ad9308SMatthew Knepley Level: advanced 132130ad9308SMatthew Knepley 132230ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 132330ad9308SMatthew Knepley @*/ 1324a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 132530ad9308SMatthew Knepley { 132630ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 132730ad9308SMatthew Knepley 132830ad9308SMatthew Knepley PetscFunctionBegin; 13290700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1330c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1331a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1332a04f6461SBarry Smith if (A01) *A01 = jac->B; 1333a04f6461SBarry Smith if (A10) *A10 = jac->C; 1334a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 133530ad9308SMatthew Knepley PetscFunctionReturn(0); 133630ad9308SMatthew Knepley } 133730ad9308SMatthew Knepley 133879416396SBarry Smith EXTERN_C_BEGIN 133979416396SBarry Smith #undef __FUNCT__ 134079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 13417087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 134279416396SBarry Smith { 134379416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1344e69d4d44SBarry Smith PetscErrorCode ierr; 134579416396SBarry Smith 134679416396SBarry Smith PetscFunctionBegin; 134779416396SBarry Smith jac->type = type; 13483b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 13493b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 13503b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1351e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1352e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1353e69d4d44SBarry Smith 13543b224e63SBarry Smith } else { 13553b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 13563b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1357e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13589e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 13593b224e63SBarry Smith } 136079416396SBarry Smith PetscFunctionReturn(0); 136179416396SBarry Smith } 136279416396SBarry Smith EXTERN_C_END 136379416396SBarry Smith 136451f519a2SBarry Smith EXTERN_C_BEGIN 136551f519a2SBarry Smith #undef __FUNCT__ 136651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 13677087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 136851f519a2SBarry Smith { 136951f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 137051f519a2SBarry Smith 137151f519a2SBarry Smith PetscFunctionBegin; 137265e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 137365e19b50SBarry 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); 137451f519a2SBarry Smith jac->bs = bs; 137551f519a2SBarry Smith PetscFunctionReturn(0); 137651f519a2SBarry Smith } 137751f519a2SBarry Smith EXTERN_C_END 137851f519a2SBarry Smith 137979416396SBarry Smith #undef __FUNCT__ 138079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1381bc08b0f1SBarry Smith /*@ 138279416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 138379416396SBarry Smith 138479416396SBarry Smith Collective on PC 138579416396SBarry Smith 138679416396SBarry Smith Input Parameter: 138779416396SBarry Smith . pc - the preconditioner context 138881540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 138979416396SBarry Smith 139079416396SBarry Smith Options Database Key: 1391a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 139279416396SBarry Smith 139379416396SBarry Smith Level: Developer 139479416396SBarry Smith 139579416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 139679416396SBarry Smith 139779416396SBarry Smith .seealso: PCCompositeSetType() 139879416396SBarry Smith 139979416396SBarry Smith @*/ 14007087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 140179416396SBarry Smith { 14024ac538c5SBarry Smith PetscErrorCode ierr; 140379416396SBarry Smith 140479416396SBarry Smith PetscFunctionBegin; 14050700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14064ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 140779416396SBarry Smith PetscFunctionReturn(0); 140879416396SBarry Smith } 140979416396SBarry Smith 14100971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 14110971522cSBarry Smith /*MC 1412a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1413a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 14140971522cSBarry Smith 1415edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1416edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 14170971522cSBarry Smith 1418a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 141969a612a9SBarry Smith and set the options directly on the resulting KSP object 14200971522cSBarry Smith 14210971522cSBarry Smith Level: intermediate 14220971522cSBarry Smith 142379416396SBarry Smith Options Database Keys: 142481540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 142581540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 142681540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 142781540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 142881540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1429e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1430435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 143179416396SBarry Smith 14325d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 14335d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 14345d4c12cdSJungho Lee 1435c8a0d604SMatthew G Knepley Notes: 1436c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1437d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1438d32f9abdSBarry Smith 1439d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1440d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1441d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1442d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1443d32f9abdSBarry Smith 1444c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1445c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1446c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1447c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1448c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1449a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1450c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1451c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1452c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1453a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1454c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1455c8a0d604SMatthew G Knepley factorization type is set using -pc_fieldsplit_schur_factorization_type <diag, lower, upper, full>. The full is shown above, 1456c8a0d604SMatthew G Knepley but diag gives 1457c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1458c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 1459c8a0d604SMatthew G Knepley so that the preconditioner is positive definite. The lower factorization is the inverse of 1460c8a0d604SMatthew G Knepley $ ( A00 0 ) 1461c8a0d604SMatthew G Knepley $ ( A10 S ) 1462c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1463c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1464c8a0d604SMatthew G Knepley $ ( 0 S ) 1465c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1466e69d4d44SBarry Smith 1467edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1468edf189efSBarry Smith is used automatically for a second block. 1469edf189efSBarry Smith 1470ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1471ff218e97SBarry Smith Generally it should be used with the AIJ format. 1472ff218e97SBarry Smith 1473ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1474ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1475ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 14760716a85fSBarry Smith 1477a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 14780971522cSBarry Smith 14797e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1480e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 14810971522cSBarry Smith M*/ 14820971522cSBarry Smith 14830971522cSBarry Smith EXTERN_C_BEGIN 14840971522cSBarry Smith #undef __FUNCT__ 14850971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 14867087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 14870971522cSBarry Smith { 14880971522cSBarry Smith PetscErrorCode ierr; 14890971522cSBarry Smith PC_FieldSplit *jac; 14900971522cSBarry Smith 14910971522cSBarry Smith PetscFunctionBegin; 149238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 14930971522cSBarry Smith jac->bs = -1; 14940971522cSBarry Smith jac->nsplits = 0; 14953e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1496e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1497c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 149851f519a2SBarry Smith 14990971522cSBarry Smith pc->data = (void*)jac; 15000971522cSBarry Smith 15010971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1502421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 15030971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1504574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 15050971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 15060971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 15070971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 15080971522cSBarry Smith pc->ops->applyrichardson = 0; 15090971522cSBarry Smith 151069a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 151169a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 15120971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 15130971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1514704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1515704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 151679416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 151779416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 151851f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 151951f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 15200971522cSBarry Smith PetscFunctionReturn(0); 15210971522cSBarry Smith } 15220971522cSBarry Smith EXTERN_C_END 15230971522cSBarry Smith 15240971522cSBarry Smith 1525a541d17aSBarry Smith 1526