1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 33c48a1e8SJed Brown #include <petscdmcomposite.h> /*I "petscdmcomposite.h" I*/ 40971522cSBarry Smith 5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0}; 7c5d2311dSJed Brown 80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 90971522cSBarry Smith struct _PC_FieldSplitLink { 1069a612a9SBarry Smith KSP ksp; 110971522cSBarry Smith Vec x,y; 12db4c96c1SJed Brown char *splitname; 130971522cSBarry Smith PetscInt nfields; 145d4c12cdSJungho Lee PetscInt *fields,*fields_col; 151b9fc7fcSBarry Smith VecScatter sctx; 165d4c12cdSJungho Lee IS is,is_col; 177233a360SDmitry Karpeev DM dm; 1851f519a2SBarry Smith PC_FieldSplitLink next,previous; 190971522cSBarry Smith }; 200971522cSBarry Smith 210971522cSBarry Smith typedef struct { 2268dd23aaSBarry Smith PCCompositeType type; 23ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 24ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 25ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 2630ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2730ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2879416396SBarry Smith Vec *x,*y,w1,w2; 29519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 30519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3130ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 32ace3abfcSBarry Smith PetscBool issetup; 3330ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3430ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3530ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 36a04f6461SBarry Smith Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01 */ 37084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 38084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 39c9c6ffaaSJed Brown PCFieldSplitSchurFactType schurfactorization; 4030ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 4197bbdb24SBarry Smith PC_FieldSplitLink head; 4263ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 43c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 440971522cSBarry Smith } PC_FieldSplit; 450971522cSBarry Smith 4616913363SBarry Smith /* 4716913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4816913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4916913363SBarry Smith PC you could change this. 5016913363SBarry Smith */ 51084e4875SJed Brown 52e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 53084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 54084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 55084e4875SJed Brown { 56084e4875SJed Brown switch (jac->schurpre) { 57084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 58084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 59084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 60084e4875SJed Brown default: 61084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 62084e4875SJed Brown } 63084e4875SJed Brown } 64084e4875SJed Brown 65084e4875SJed Brown 660971522cSBarry Smith #undef __FUNCT__ 670971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 680971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 690971522cSBarry Smith { 700971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 710971522cSBarry Smith PetscErrorCode ierr; 72ace3abfcSBarry Smith PetscBool iascii; 730971522cSBarry Smith PetscInt i,j; 745a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 750971522cSBarry Smith 760971522cSBarry Smith PetscFunctionBegin; 77251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 780971522cSBarry Smith if (iascii) { 792c9966d7SBarry Smith if (jac->bs > 0) { 8051f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 812c9966d7SBarry Smith } else { 822c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 832c9966d7SBarry Smith } 84a3df900dSMatthew G Knepley if (jac->realdiagonal) { 85a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 86a3df900dSMatthew G Knepley } 8769a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 880971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 890971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 901ab39975SBarry Smith if (ilink->fields) { 910971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 9279416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 935a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9479416396SBarry Smith if (j > 0) { 9579416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9679416396SBarry Smith } 975a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 980971522cSBarry Smith } 990971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 10079416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1011ab39975SBarry Smith } else { 1021ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1031ab39975SBarry Smith } 1045a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1055a9f2f41SSatish Balay ilink = ilink->next; 1060971522cSBarry Smith } 1070971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1080971522cSBarry Smith } else { 10965e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1100971522cSBarry Smith } 1110971522cSBarry Smith PetscFunctionReturn(0); 1120971522cSBarry Smith } 1130971522cSBarry Smith 1140971522cSBarry Smith #undef __FUNCT__ 1153b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1163b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1173b224e63SBarry Smith { 1183b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1193b224e63SBarry Smith PetscErrorCode ierr; 120ace3abfcSBarry Smith PetscBool iascii; 1213b224e63SBarry Smith PetscInt i,j; 1223b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1233b224e63SBarry Smith KSP ksp; 1243b224e63SBarry Smith 1253b224e63SBarry Smith PetscFunctionBegin; 126251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1273b224e63SBarry Smith if (iascii) { 1282c9966d7SBarry Smith if (jac->bs > 0) { 129c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1302c9966d7SBarry Smith } else { 1312c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1322c9966d7SBarry Smith } 133a3df900dSMatthew G Knepley if (jac->realdiagonal) { 134a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 135a3df900dSMatthew G Knepley } 1363e8b8b31SMatthew G Knepley switch(jac->schurpre) { 1373e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1383e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 139b02e2d75SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_DIAG: 1403e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break; 1413e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1423e8b8b31SMatthew G Knepley if (jac->schur_user) { 1433e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1443e8b8b31SMatthew G Knepley } else { 1453e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr); 1463e8b8b31SMatthew G Knepley } 1473e8b8b31SMatthew G Knepley break; 1483e8b8b31SMatthew G Knepley default: 1493e8b8b31SMatthew G Knepley SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1503e8b8b31SMatthew G Knepley } 1513b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1523b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1533b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1543b224e63SBarry Smith if (ilink->fields) { 1553b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1563b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1573b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1583b224e63SBarry Smith if (j > 0) { 1593b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1603b224e63SBarry Smith } 1613b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1623b224e63SBarry Smith } 1633b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1643b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1653b224e63SBarry Smith } else { 1663b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1673b224e63SBarry Smith } 1683b224e63SBarry Smith ilink = ilink->next; 1693b224e63SBarry Smith } 170435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1713b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 17212cae6f2SJed Brown if (jac->schur) { 17312cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1743b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 17512cae6f2SJed Brown } else { 17612cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 17712cae6f2SJed Brown } 1783b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 179435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1803b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18112cae6f2SJed Brown if (jac->kspschur) { 1823b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 18312cae6f2SJed Brown } else { 18412cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 18512cae6f2SJed Brown } 1863b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1873b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1883b224e63SBarry Smith } else { 18965e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1903b224e63SBarry Smith } 1913b224e63SBarry Smith PetscFunctionReturn(0); 1923b224e63SBarry Smith } 1933b224e63SBarry Smith 1943b224e63SBarry Smith #undef __FUNCT__ 1956c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1966c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1976c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1986c924f48SJed Brown { 1996c924f48SJed Brown PetscErrorCode ierr; 2006c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2015d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2025d4c12cdSJungho Lee PetscBool flg,flg_col; 2035d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2046c924f48SJed Brown 2056c924f48SJed Brown PetscFunctionBegin; 2066c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2075d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2086c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2096c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2106c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2115d4c12cdSJungho Lee ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2126c924f48SJed Brown nfields = jac->bs; 21329499fbbSJungho Lee nfields_col = jac->bs; 2146c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2155d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2166c924f48SJed Brown if (!flg) break; 2175d4c12cdSJungho Lee else if (flg && !flg_col) { 2185d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2195d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2205d4c12cdSJungho Lee } 2215d4c12cdSJungho Lee else { 2225d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2235d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2245d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2255d4c12cdSJungho Lee } 2266c924f48SJed Brown } 2276c924f48SJed Brown if (i > 0) { 2286c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2296c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2306c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2316c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2326c924f48SJed Brown } 2336c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2345d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2356c924f48SJed Brown PetscFunctionReturn(0); 2366c924f48SJed Brown } 2376c924f48SJed Brown 2386c924f48SJed Brown #undef __FUNCT__ 23969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 24069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2410971522cSBarry Smith { 2420971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2430971522cSBarry Smith PetscErrorCode ierr; 2445a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2456ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2466c924f48SJed Brown PetscInt i; 2470971522cSBarry Smith 2480971522cSBarry Smith PetscFunctionBegin; 249d32f9abdSBarry Smith if (!ilink) { 250d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 251d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2527b62db95SJungho Lee PetscInt numFields, f; 2530784a22cSJed Brown char **fieldNames; 2547b62db95SJungho Lee IS *fields; 255e7c4fc90SDmitry Karpeev DM *dms; 25616621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms); CHKERRQ(ierr); 2577b62db95SJungho Lee for(f = 0; f < numFields; ++f) { 2587b62db95SJungho Lee ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 2597b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 2607b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 2617b62db95SJungho Lee } 2627b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 2637b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 264e7c4fc90SDmitry Karpeev if(dms) { 2658b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2667b62db95SJungho Lee for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) { 26737d9a391SDmitry Karpeev ierr = DMDestroy(&ilink->dm);CHKERRQ(ierr); 2687233a360SDmitry Karpeev ilink->dm = dms[i]; 26937d9a391SDmitry Karpeev if(ilink->dm) { 27037d9a391SDmitry Karpeev ierr = PetscObjectReference((PetscObject)ilink->dm);CHKERRQ(ierr); 27137d9a391SDmitry Karpeev } 2727b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 2737b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 274e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 2752fa5ba8aSJed Brown } 2767b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 2778b8307b2SJed Brown } 27866ffff09SJed Brown } else { 279521d7252SBarry Smith if (jac->bs <= 0) { 280704ba839SBarry Smith if (pc->pmat) { 281521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 282704ba839SBarry Smith } else { 283704ba839SBarry Smith jac->bs = 1; 284704ba839SBarry Smith } 285521d7252SBarry Smith } 286d32f9abdSBarry Smith 287acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2886ce1633cSBarry Smith if (stokes) { 2896ce1633cSBarry Smith IS zerodiags,rest; 2906ce1633cSBarry Smith PetscInt nmin,nmax; 2916ce1633cSBarry Smith 2926ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2936ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2946ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2956ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2966ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 297fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 298fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2996ce1633cSBarry Smith } else { 300d32f9abdSBarry Smith if (!flg) { 301d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 302d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 3036c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 3046c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 305d32f9abdSBarry Smith } 3066c924f48SJed Brown if (flg || !jac->splitdefined) { 307d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 308db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 3096c924f48SJed Brown char splitname[8]; 3106c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 3115d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 31279416396SBarry Smith } 3135d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 314521d7252SBarry Smith } 31566ffff09SJed Brown } 3166ce1633cSBarry Smith } 317edf189efSBarry Smith } else if (jac->nsplits == 1) { 318edf189efSBarry Smith if (ilink->is) { 319edf189efSBarry Smith IS is2; 320edf189efSBarry Smith PetscInt nmin,nmax; 321edf189efSBarry Smith 322edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 323edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 324db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 325fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 3267b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 32763ec74ffSBarry Smith } else if (jac->reset) { 32863ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 329d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 330d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 331d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 332d0af7cd3SBarry Smith if (pc->dm && !stokes) { 333d0af7cd3SBarry Smith PetscBool dmcomposite; 334251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 335d0af7cd3SBarry Smith if (dmcomposite) { 336d0af7cd3SBarry Smith PetscInt nDM; 337d0af7cd3SBarry Smith IS *fields; 3387b62db95SJungho Lee DM *dms; 339d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 340d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 341d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 3427b62db95SJungho Lee ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr); 3437b62db95SJungho Lee ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr); 344d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 34537d9a391SDmitry Karpeev ierr = DMDestroy(&ilink->dm);CHKERRQ(ierr); 34637d9a391SDmitry Karpeev ilink->dm = dms[i]; 34737d9a391SDmitry Karpeev if(ilink->dm) { 34837d9a391SDmitry Karpeev ierr = PetscObjectReference((PetscObject)ilink->dm);CHKERRQ(ierr); 34937d9a391SDmitry Karpeev } 3507b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr); 3517b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr); 352d0af7cd3SBarry Smith ilink->is = fields[i]; 353d0af7cd3SBarry Smith ilink = ilink->next; 354edf189efSBarry Smith } 355d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 3567b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 357d0af7cd3SBarry Smith } 358d0af7cd3SBarry Smith } else { 359d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 360d0af7cd3SBarry Smith if (stokes) { 361d0af7cd3SBarry Smith IS zerodiags,rest; 362d0af7cd3SBarry Smith PetscInt nmin,nmax; 363d0af7cd3SBarry Smith 364d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 365d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 366d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 36720b26d62SBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 36820b26d62SBarry Smith ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr); 369d0af7cd3SBarry Smith ilink->is = rest; 370d0af7cd3SBarry Smith ilink->next->is = zerodiags; 37120b26d62SBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 372d0af7cd3SBarry Smith } 373d0af7cd3SBarry Smith } 374d0af7cd3SBarry Smith 3757b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 37669a612a9SBarry Smith PetscFunctionReturn(0); 37769a612a9SBarry Smith } 37869a612a9SBarry Smith 37969a612a9SBarry Smith #undef __FUNCT__ 38069a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 38169a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 38269a612a9SBarry Smith { 38369a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 38469a612a9SBarry Smith PetscErrorCode ierr; 3855a9f2f41SSatish Balay PC_FieldSplitLink ilink; 3862c9966d7SBarry Smith PetscInt i,nsplit; 38769a612a9SBarry Smith MatStructure flag = pc->flag; 3885f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 38969a612a9SBarry Smith 39069a612a9SBarry Smith PetscFunctionBegin; 39169a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 39297bbdb24SBarry Smith nsplit = jac->nsplits; 3935a9f2f41SSatish Balay ilink = jac->head; 39497bbdb24SBarry Smith 39597bbdb24SBarry Smith /* get the matrices for each split */ 396704ba839SBarry Smith if (!jac->issetup) { 3971b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 39897bbdb24SBarry Smith 399704ba839SBarry Smith jac->issetup = PETSC_TRUE; 400704ba839SBarry Smith 4015d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4022c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4032c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4042c9966d7SBarry Smith } 40551f519a2SBarry Smith bs = jac->bs; 40697bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4071b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4081b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4091b9fc7fcSBarry Smith if (jac->defaultsplit) { 410704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4115f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 412704ba839SBarry Smith } else if (!ilink->is) { 413ccb205f8SBarry Smith if (ilink->nfields > 1) { 4145f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4155a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4165f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4171b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4181b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4191b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4205f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 42197bbdb24SBarry Smith } 42297bbdb24SBarry Smith } 423d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4245f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 425ccb205f8SBarry Smith } else { 426704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4275f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 428ccb205f8SBarry Smith } 4293e197d65SBarry Smith } 430704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4315f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4325f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 433704ba839SBarry Smith ilink = ilink->next; 4341b9fc7fcSBarry Smith } 4351b9fc7fcSBarry Smith } 4361b9fc7fcSBarry Smith 437704ba839SBarry Smith ilink = jac->head; 43897bbdb24SBarry Smith if (!jac->pmat) { 439bdddcaaaSMatthew G Knepley Vec xtmp; 440bdddcaaaSMatthew G Knepley 441bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 442cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 443bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 444cf502942SBarry Smith for (i=0; i<nsplit; i++) { 445bdddcaaaSMatthew G Knepley MatNullSpace sp; 446bdddcaaaSMatthew G Knepley 447a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 448a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 449a3df900dSMatthew G Knepley if (jac->pmat[i]) { 450a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 451a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 452a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 453a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 454a3df900dSMatthew G Knepley } 455a3df900dSMatthew G Knepley } else { 4565f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 457a3df900dSMatthew G Knepley } 458bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 459bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 460bdddcaaaSMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; 461bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 462bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 463bdddcaaaSMatthew G Knepley /* HACK: Check for the constant null space */ 464f1dc53b5SDmitry Karpeev ierr = MatGetNullSpace(pc->pmat, &sp);CHKERRQ(ierr); 465bdddcaaaSMatthew G Knepley if (sp) { 466bdddcaaaSMatthew G Knepley MatNullSpace subsp; 467bdddcaaaSMatthew G Knepley Vec ftmp, gtmp; 4689583d628SBarry Smith PetscReal norm; 469bdddcaaaSMatthew G Knepley PetscInt N; 470f1dc53b5SDmitry Karpeev 471f1dc53b5SDmitry Karpeev ierr = MatGetVecs(pc->pmat, >mp, PETSC_NULL);CHKERRQ(ierr); 472f1dc53b5SDmitry Karpeev ierr = MatGetVecs(jac->pmat[i], &ftmp, PETSC_NULL);CHKERRQ(ierr); 473bdddcaaaSMatthew G Knepley ierr = VecGetSize(ftmp, &N);CHKERRQ(ierr); 474bdddcaaaSMatthew G Knepley ierr = VecSet(ftmp, 1.0/N);CHKERRQ(ierr); 475f1dc53b5SDmitry Karpeev ierr = VecScatterBegin(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 476f1dc53b5SDmitry Karpeev ierr = VecScatterEnd(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 477f1dc53b5SDmitry Karpeev ierr = MatNullSpaceRemove(sp, gtmp, PETSC_NULL);CHKERRQ(ierr); 478bdddcaaaSMatthew G Knepley ierr = VecNorm(gtmp, NORM_2, &norm);CHKERRQ(ierr); 479bdddcaaaSMatthew G Knepley if (norm < 1.0e-10) { 480bdddcaaaSMatthew G Knepley ierr = MatNullSpaceCreate(((PetscObject)pc)->comm, PETSC_TRUE, 0, PETSC_NULL, &subsp);CHKERRQ(ierr); 481bdddcaaaSMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], subsp);CHKERRQ(ierr); 482bdddcaaaSMatthew G Knepley ierr = MatNullSpaceDestroy(&subsp);CHKERRQ(ierr); 483bdddcaaaSMatthew G Knepley } 484bdddcaaaSMatthew G Knepley ierr = VecDestroy(&ftmp);CHKERRQ(ierr); 485bdddcaaaSMatthew G Knepley ierr = VecDestroy(>mp);CHKERRQ(ierr); 486bdddcaaaSMatthew G Knepley } 487ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 488ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 489ed1f0337SMatthew G Knepley if (sp) { 490ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 491ed1f0337SMatthew G Knepley } 492704ba839SBarry Smith ilink = ilink->next; 493cf502942SBarry Smith } 494bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 49597bbdb24SBarry Smith } else { 496cf502942SBarry Smith for (i=0; i<nsplit; i++) { 497a3df900dSMatthew G Knepley Mat pmat; 498a3df900dSMatthew G Knepley 499a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 500a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 501a3df900dSMatthew G Knepley if (!pmat) { 5025f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 503a3df900dSMatthew G Knepley } 504704ba839SBarry Smith ilink = ilink->next; 505cf502942SBarry Smith } 50697bbdb24SBarry Smith } 507519d70e2SJed Brown if (jac->realdiagonal) { 508519d70e2SJed Brown ilink = jac->head; 509519d70e2SJed Brown if (!jac->mat) { 510519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 511519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5125f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 513519d70e2SJed Brown ilink = ilink->next; 514519d70e2SJed Brown } 515519d70e2SJed Brown } else { 516519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5175f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 518519d70e2SJed Brown ilink = ilink->next; 519519d70e2SJed Brown } 520519d70e2SJed Brown } 521519d70e2SJed Brown } else { 522519d70e2SJed Brown jac->mat = jac->pmat; 523519d70e2SJed Brown } 52497bbdb24SBarry Smith 5256c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 52668dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 52768dd23aaSBarry Smith ilink = jac->head; 52868dd23aaSBarry Smith if (!jac->Afield) { 52968dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 53068dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5314aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 53268dd23aaSBarry Smith ilink = ilink->next; 53368dd23aaSBarry Smith } 53468dd23aaSBarry Smith } else { 53568dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5364aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 53768dd23aaSBarry Smith ilink = ilink->next; 53868dd23aaSBarry Smith } 53968dd23aaSBarry Smith } 54068dd23aaSBarry Smith } 54168dd23aaSBarry Smith 5423b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5433b224e63SBarry Smith IS ccis; 5444aa3045dSJed Brown PetscInt rstart,rend; 545093c86ffSJed Brown char lscname[256]; 546093c86ffSJed Brown PetscObject LSC_L; 547e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 54868dd23aaSBarry Smith 549e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 550e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 551e6cab6aaSJed Brown 552*68bd789dSDmitry Karpeev /* 553*68bd789dSDmitry Karpeev Set from options only the A00 split. The other split's solver won't be used with Schur. 554*68bd789dSDmitry Karpeev Should it be destroyed? Should KSPCreate() be moved here from PCFieldSplitSetIS() and invoked 555*68bd789dSDmitry Karpeev only when necessary? 556*68bd789dSDmitry Karpeev */ 557*68bd789dSDmitry Karpeev ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 558*68bd789dSDmitry Karpeev 5593b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5603b224e63SBarry Smith if (jac->schur) { 5613b224e63SBarry Smith ilink = jac->head; 56249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5634aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 564fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5653b224e63SBarry Smith ilink = ilink->next; 56649bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5674aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 568fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 569a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 570084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 5713b224e63SBarry Smith 5723b224e63SBarry Smith } else { 5731cee3971SBarry Smith KSP ksp; 5746c924f48SJed Brown char schurprefix[256]; 575bdddcaaaSMatthew G Knepley MatNullSpace sp; 5763b224e63SBarry Smith 577a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 5783b224e63SBarry Smith ilink = jac->head; 57949bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5804aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 581fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5823b224e63SBarry Smith ilink = ilink->next; 58349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5844aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 585fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5867d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 587519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 588bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 5891d33e650SDmitry Karpeev /* 5901d33e650SDmitry Karpeev Do we really want to attach the A11-block's nullspace to S? Note that this may generate 5911d33e650SDmitry Karpeev "false positives", when the A11's nullspace isn't S's: Stokes or A = [1, 1; 1, 0]. 5921d33e650SDmitry Karpeev Using a false nullspace may prevent KSP(S) from converging, since it might force an inconsistent rhs. 5931d33e650SDmitry Karpeev */ 5941d33e650SDmitry Karpeev /* if (sp) {ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);} */ 5951d33e650SDmitry Karpeev /* set tabbing, options prefix and DM of KSP inside the MatSchurComplement */ 5961cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 597e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 5987233a360SDmitry Karpeev { 5997233a360SDmitry Karpeev PC pcinner; 6007233a360SDmitry Karpeev ierr = KSPGetPC(ksp, &pcinner); CHKERRQ(ierr); 6017233a360SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)pcinner,(PetscObject)pc,2);CHKERRQ(ierr); 6027233a360SDmitry Karpeev } 603a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 604a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 6057233a360SDmitry Karpeev /* Can't do KSPGetDM(jac->head->ksp,&dminner); KSPSetDM(ksp,dminner): KSPGetDM() will create DMShell, if the DM hasn't been set - not what we want. */ 6067233a360SDmitry Karpeev ierr = KSPSetDM(ksp,jac->head->dm); CHKERRQ(ierr); 6077233a360SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE); CHKERRQ(ierr); 6087233a360SDmitry Karpeev ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 60920b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 61020b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 6117b62db95SJungho Lee ierr = MatSetUp(jac->schur);CHKERRQ(ierr); 6127233a360SDmitry Karpeev /* Create and set up the KSP for the Schur complement; forward the second split's DM and set up tabbing, including for the contained PC. */ 6133b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 6149005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 6151cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 6167233a360SDmitry Karpeev { 6177233a360SDmitry Karpeev PC pcschur; 6187233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur, &pcschur); CHKERRQ(ierr); 6197233a360SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)pcschur,(PetscObject)pc,1);CHKERRQ(ierr); 6207233a360SDmitry Karpeev } 6217233a360SDmitry Karpeev /* Can't do KSPGetDM(ilink->ksp,&dmschur); KSPSetDM(kspshur,dmschur): KSPGetDM() will create DMShell, if the DM hasn't been set - not what we want. */ 6227233a360SDmitry Karpeev ierr = KSPSetDM(jac->kspschur,ilink->dm); CHKERRQ(ierr); 6237233a360SDmitry Karpeev ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE); CHKERRQ(ierr); 6247233a360SDmitry Karpeev 625084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 626084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 6277233a360SDmitry Karpeev PC pcschur; 6287233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 6297233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 630084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 631e69d4d44SBarry Smith } 6323b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 63320b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 6347233a360SDmitry Karpeev ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 6357233a360SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 63620b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur); CHKERRQ(ierr); 6373b224e63SBarry Smith } 638093c86ffSJed Brown 639093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 640093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 641093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 642093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 643093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 644093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 645093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 646093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 647093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 6487233a360SDmitry Karpeev } 649*68bd789dSDmitry Karpeev else { 650*68bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 65197bbdb24SBarry Smith i = 0; 6525a9f2f41SSatish Balay ilink = jac->head; 6535a9f2f41SSatish Balay while (ilink) { 654519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag); CHKERRQ(ierr); 6553b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 656*68bd789dSDmitry Karpeev if (!jac->suboptionsset) { 657*68bd789dSDmitry Karpeev ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 658*68bd789dSDmitry Karpeev } 65997bbdb24SBarry Smith i++; 6605a9f2f41SSatish Balay ilink = ilink->next; 6610971522cSBarry Smith } 662*68bd789dSDmitry Karpeev } 6633b224e63SBarry Smith 664c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 6650971522cSBarry Smith PetscFunctionReturn(0); 6660971522cSBarry Smith } 6670971522cSBarry Smith 6685a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 669ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 670ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 6715a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 672ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 673ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 67479416396SBarry Smith 6750971522cSBarry Smith #undef __FUNCT__ 6763b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 6773b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 6783b224e63SBarry Smith { 6793b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6803b224e63SBarry Smith PetscErrorCode ierr; 6813b224e63SBarry Smith KSP ksp; 6823b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 6833b224e63SBarry Smith 6843b224e63SBarry Smith PetscFunctionBegin; 6853b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 6863b224e63SBarry Smith 687c5d2311dSJed Brown switch (jac->schurfactorization) { 688c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 689a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 690c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 691c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 692c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 693c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 694c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 695c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 696c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 697c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 698c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 699c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 700c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 701c5d2311dSJed Brown break; 702c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 703a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 704c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 705c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 706c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 707c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 708c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 709c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 710c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 711c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 712c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 713c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 714c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 715c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 716c5d2311dSJed Brown break; 717c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 718a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 719c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 720c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 721c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 722c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 723c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 724c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 725c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 726c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 727c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 728c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 729c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 730c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 731c5d2311dSJed Brown break; 732c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 733a04f6461SBarry 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 */ 7343b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7353b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7363b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 7373b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 7383b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 7393b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7403b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7413b224e63SBarry Smith 7423b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 7433b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7443b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7453b224e63SBarry Smith 7463b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 7473b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 7483b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 7493b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7503b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 751c5d2311dSJed Brown } 7523b224e63SBarry Smith PetscFunctionReturn(0); 7533b224e63SBarry Smith } 7543b224e63SBarry Smith 7553b224e63SBarry Smith #undef __FUNCT__ 7560971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 7570971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 7580971522cSBarry Smith { 7590971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7600971522cSBarry Smith PetscErrorCode ierr; 7615a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 762939b8a20SBarry Smith PetscInt cnt,bs; 7630971522cSBarry Smith 7640971522cSBarry Smith PetscFunctionBegin; 76551f519a2SBarry Smith CHKMEMQ; 76679416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 7671b9fc7fcSBarry Smith if (jac->defaultsplit) { 768939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 769939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 770939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 771939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 7720971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 7735a9f2f41SSatish Balay while (ilink) { 7745a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7755a9f2f41SSatish Balay ilink = ilink->next; 7760971522cSBarry Smith } 7770971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 7781b9fc7fcSBarry Smith } else { 779efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7805a9f2f41SSatish Balay while (ilink) { 7815a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 7825a9f2f41SSatish Balay ilink = ilink->next; 7831b9fc7fcSBarry Smith } 7841b9fc7fcSBarry Smith } 78516913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 78679416396SBarry Smith if (!jac->w1) { 78779416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 78879416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 78979416396SBarry Smith } 790efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7915a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 7923e197d65SBarry Smith cnt = 1; 7935a9f2f41SSatish Balay while (ilink->next) { 7945a9f2f41SSatish Balay ilink = ilink->next; 7953e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 7963e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 7973e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 7983e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7993e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8003e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8013e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8023e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8033e197d65SBarry Smith } 80451f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 80511755939SBarry Smith cnt -= 2; 80651f519a2SBarry Smith while (ilink->previous) { 80751f519a2SBarry Smith ilink = ilink->previous; 80811755939SBarry Smith /* compute the residual only over the part of the vector needed */ 80911755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 81011755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 81111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81311755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 81411755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 81511755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 81651f519a2SBarry Smith } 81711755939SBarry Smith } 81865e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 81951f519a2SBarry Smith CHKMEMQ; 8200971522cSBarry Smith PetscFunctionReturn(0); 8210971522cSBarry Smith } 8220971522cSBarry Smith 823421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 824ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 825ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 826421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 827ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 828ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 829421e10b8SBarry Smith 830421e10b8SBarry Smith #undef __FUNCT__ 8318c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 832421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 833421e10b8SBarry Smith { 834421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 835421e10b8SBarry Smith PetscErrorCode ierr; 836421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 837939b8a20SBarry Smith PetscInt bs; 838421e10b8SBarry Smith 839421e10b8SBarry Smith PetscFunctionBegin; 840421e10b8SBarry Smith CHKMEMQ; 841421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 842421e10b8SBarry Smith if (jac->defaultsplit) { 843939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 844939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 845939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 846939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 847421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 848421e10b8SBarry Smith while (ilink) { 849421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 850421e10b8SBarry Smith ilink = ilink->next; 851421e10b8SBarry Smith } 852421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 853421e10b8SBarry Smith } else { 854421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 855421e10b8SBarry Smith while (ilink) { 856421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 857421e10b8SBarry Smith ilink = ilink->next; 858421e10b8SBarry Smith } 859421e10b8SBarry Smith } 860421e10b8SBarry Smith } else { 861421e10b8SBarry Smith if (!jac->w1) { 862421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 863421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 864421e10b8SBarry Smith } 865421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 866421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 867421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 868421e10b8SBarry Smith while (ilink->next) { 869421e10b8SBarry Smith ilink = ilink->next; 8709989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 871421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 872421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 873421e10b8SBarry Smith } 874421e10b8SBarry Smith while (ilink->previous) { 875421e10b8SBarry Smith ilink = ilink->previous; 8769989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 877421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 878421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 879421e10b8SBarry Smith } 880421e10b8SBarry Smith } else { 881421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 882421e10b8SBarry Smith ilink = ilink->next; 883421e10b8SBarry Smith } 884421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 885421e10b8SBarry Smith while (ilink->previous) { 886421e10b8SBarry Smith ilink = ilink->previous; 8879989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 888421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 889421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 890421e10b8SBarry Smith } 891421e10b8SBarry Smith } 892421e10b8SBarry Smith } 893421e10b8SBarry Smith CHKMEMQ; 894421e10b8SBarry Smith PetscFunctionReturn(0); 895421e10b8SBarry Smith } 896421e10b8SBarry Smith 8970971522cSBarry Smith #undef __FUNCT__ 898574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 899574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9000971522cSBarry Smith { 9010971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9020971522cSBarry Smith PetscErrorCode ierr; 9035a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9040971522cSBarry Smith 9050971522cSBarry Smith PetscFunctionBegin; 9065a9f2f41SSatish Balay while (ilink) { 907574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 90837d9a391SDmitry Karpeev ierr = DMDestroy(&ilink->dm);CHKERRQ(ierr); 909fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 910fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 911fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 912fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 913b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9145a9f2f41SSatish Balay next = ilink->next; 9155a9f2f41SSatish Balay ilink = next; 9160971522cSBarry Smith } 91705b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 918574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 919574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 920574deadeSBarry Smith } else if (jac->mat) { 921574deadeSBarry Smith jac->mat = PETSC_NULL; 922574deadeSBarry Smith } 92397bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 92468dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 9256bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 9266bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 9276bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 9286bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 9296bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 9306bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 9316bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 93263ec74ffSBarry Smith jac->reset = PETSC_TRUE; 933574deadeSBarry Smith PetscFunctionReturn(0); 934574deadeSBarry Smith } 935574deadeSBarry Smith 936574deadeSBarry Smith #undef __FUNCT__ 937574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 938574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 939574deadeSBarry Smith { 940574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 941574deadeSBarry Smith PetscErrorCode ierr; 942574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 943574deadeSBarry Smith 944574deadeSBarry Smith PetscFunctionBegin; 945574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 946574deadeSBarry Smith while (ilink) { 9476bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 948574deadeSBarry Smith next = ilink->next; 949574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 950574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 9515d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 952574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 953574deadeSBarry Smith ilink = next; 954574deadeSBarry Smith } 955574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 956c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 9577b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 9587b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 9597b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 9607b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 9617b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 962ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 963c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 9640971522cSBarry Smith PetscFunctionReturn(0); 9650971522cSBarry Smith } 9660971522cSBarry Smith 9670971522cSBarry Smith #undef __FUNCT__ 9680971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 9690971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 9700971522cSBarry Smith { 9711b9fc7fcSBarry Smith PetscErrorCode ierr; 9726c924f48SJed Brown PetscInt bs; 973bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 9749dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9753b224e63SBarry Smith PCCompositeType ctype; 976e7c4fc90SDmitry Karpeev DM ddm; 977e7c4fc90SDmitry Karpeev char ddm_name[1025]; 9781b9fc7fcSBarry Smith 9790971522cSBarry Smith PetscFunctionBegin; 9801b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 981e7c4fc90SDmitry Karpeev if(pc->dm) { 982e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 983e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 984731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 985e7c4fc90SDmitry Karpeev if(flg) { 98616621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 987e7c4fc90SDmitry Karpeev if(!ddm) { 98816621825SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 989e7c4fc90SDmitry Karpeev } 99016621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 991e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 992e7c4fc90SDmitry Karpeev } 993e7c4fc90SDmitry Karpeev } 994acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 99551f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 99651f519a2SBarry Smith if (flg) { 99751f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 99851f519a2SBarry Smith } 999704ba839SBarry Smith 1000435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 1001c0adfefeSBarry Smith if (stokes) { 1002c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1003c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1004c0adfefeSBarry Smith } 1005c0adfefeSBarry Smith 10063b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10073b224e63SBarry Smith if (flg) { 10083b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10093b224e63SBarry Smith } 1010c30613efSMatthew Knepley /* Only setup fields once */ 1011c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1012d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1013d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10146c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10156c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1016d32f9abdSBarry Smith } 1017c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1018c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1019c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 1020c9c6ffaaSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_fact_type","Which off-diagonal parts of the block factorization to use","PCFieldSplitSetSchurFactType",PCFieldSplitSchurFactTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 1021084e4875SJed 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); 1022c5d2311dSJed Brown } 10231b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10240971522cSBarry Smith PetscFunctionReturn(0); 10250971522cSBarry Smith } 10260971522cSBarry Smith 10270971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10280971522cSBarry Smith 10290971522cSBarry Smith EXTERN_C_BEGIN 10300971522cSBarry Smith #undef __FUNCT__ 10310971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 10325d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10330971522cSBarry Smith { 103497bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10350971522cSBarry Smith PetscErrorCode ierr; 10365a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 103769a612a9SBarry Smith char prefix[128]; 10385d4c12cdSJungho Lee PetscInt i; 10390971522cSBarry Smith 10400971522cSBarry Smith PetscFunctionBegin; 10416c924f48SJed Brown if (jac->splitdefined) { 10426c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10436c924f48SJed Brown PetscFunctionReturn(0); 10446c924f48SJed Brown } 104551f519a2SBarry Smith for (i=0; i<n; i++) { 1046e32f2f54SBarry 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); 1047e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 104851f519a2SBarry Smith } 1049704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1050a04f6461SBarry Smith if (splitname) { 1051db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1052a04f6461SBarry Smith } else { 1053a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1054a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1055a04f6461SBarry Smith } 1056704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 10575a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 10585d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 10595d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 10605a9f2f41SSatish Balay ilink->nfields = n; 10615a9f2f41SSatish Balay ilink->next = PETSC_NULL; 10627adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10631cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 10645a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10659005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 106669a612a9SBarry Smith 1067a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 10685a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 10690971522cSBarry Smith 10700971522cSBarry Smith if (!next) { 10715a9f2f41SSatish Balay jac->head = ilink; 107251f519a2SBarry Smith ilink->previous = PETSC_NULL; 10730971522cSBarry Smith } else { 10740971522cSBarry Smith while (next->next) { 10750971522cSBarry Smith next = next->next; 10760971522cSBarry Smith } 10775a9f2f41SSatish Balay next->next = ilink; 107851f519a2SBarry Smith ilink->previous = next; 10790971522cSBarry Smith } 10800971522cSBarry Smith jac->nsplits++; 10810971522cSBarry Smith PetscFunctionReturn(0); 10820971522cSBarry Smith } 10830971522cSBarry Smith EXTERN_C_END 10840971522cSBarry Smith 1085e69d4d44SBarry Smith EXTERN_C_BEGIN 1086e69d4d44SBarry Smith #undef __FUNCT__ 1087e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 10887087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1089e69d4d44SBarry Smith { 1090e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1091e69d4d44SBarry Smith PetscErrorCode ierr; 1092e69d4d44SBarry Smith 1093e69d4d44SBarry Smith PetscFunctionBegin; 1094519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1095e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1096e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 109713e0d083SBarry Smith if (n) *n = jac->nsplits; 1098e69d4d44SBarry Smith PetscFunctionReturn(0); 1099e69d4d44SBarry Smith } 1100e69d4d44SBarry Smith EXTERN_C_END 11010971522cSBarry Smith 11020971522cSBarry Smith EXTERN_C_BEGIN 11030971522cSBarry Smith #undef __FUNCT__ 110469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11057087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11060971522cSBarry Smith { 11070971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11080971522cSBarry Smith PetscErrorCode ierr; 11090971522cSBarry Smith PetscInt cnt = 0; 11105a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11110971522cSBarry Smith 11120971522cSBarry Smith PetscFunctionBegin; 11135d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11145a9f2f41SSatish Balay while (ilink) { 11155a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11165a9f2f41SSatish Balay ilink = ilink->next; 11170971522cSBarry Smith } 11185d480477SMatthew 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); 111913e0d083SBarry Smith if (n) *n = jac->nsplits; 11200971522cSBarry Smith PetscFunctionReturn(0); 11210971522cSBarry Smith } 11220971522cSBarry Smith EXTERN_C_END 11230971522cSBarry Smith 1124704ba839SBarry Smith EXTERN_C_BEGIN 1125704ba839SBarry Smith #undef __FUNCT__ 1126704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11277087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1128704ba839SBarry Smith { 1129704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1130704ba839SBarry Smith PetscErrorCode ierr; 1131704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1132704ba839SBarry Smith char prefix[128]; 1133704ba839SBarry Smith 1134704ba839SBarry Smith PetscFunctionBegin; 11356c924f48SJed Brown if (jac->splitdefined) { 11366c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11376c924f48SJed Brown PetscFunctionReturn(0); 11386c924f48SJed Brown } 113916913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1140a04f6461SBarry Smith if (splitname) { 1141db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1142a04f6461SBarry Smith } else { 1143b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1144b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1145a04f6461SBarry Smith } 11461ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1147b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1148b5787286SJed Brown ilink->is = is; 1149b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1150b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1151b5787286SJed Brown ilink->is_col = is; 1152704ba839SBarry Smith ilink->next = PETSC_NULL; 1153704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 11541cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11557233a360SDmitry Karpeev { 11567233a360SDmitry Karpeev PC ilinkpc; 11577233a360SDmitry Karpeev ierr = KSPGetPC(ilink->ksp, &ilinkpc); CHKERRQ(ierr); 11587233a360SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)ilinkpc,(PetscObject)pc,1);CHKERRQ(ierr); 11597233a360SDmitry Karpeev } 1160704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11619005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1162704ba839SBarry Smith 1163a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1164704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1165704ba839SBarry Smith 1166704ba839SBarry Smith if (!next) { 1167704ba839SBarry Smith jac->head = ilink; 1168704ba839SBarry Smith ilink->previous = PETSC_NULL; 1169704ba839SBarry Smith } else { 1170704ba839SBarry Smith while (next->next) { 1171704ba839SBarry Smith next = next->next; 1172704ba839SBarry Smith } 1173704ba839SBarry Smith next->next = ilink; 1174704ba839SBarry Smith ilink->previous = next; 1175704ba839SBarry Smith } 1176704ba839SBarry Smith jac->nsplits++; 1177704ba839SBarry Smith 1178704ba839SBarry Smith PetscFunctionReturn(0); 1179704ba839SBarry Smith } 1180704ba839SBarry Smith EXTERN_C_END 1181704ba839SBarry Smith 11820971522cSBarry Smith #undef __FUNCT__ 11830971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 11840971522cSBarry Smith /*@ 11850971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 11860971522cSBarry Smith 1187ad4df100SBarry Smith Logically Collective on PC 11880971522cSBarry Smith 11890971522cSBarry Smith Input Parameters: 11900971522cSBarry Smith + pc - the preconditioner context 1191a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 11920971522cSBarry Smith . n - the number of fields in this split 1193db4c96c1SJed Brown - fields - the fields in this split 11940971522cSBarry Smith 11950971522cSBarry Smith Level: intermediate 11960971522cSBarry Smith 1197d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1198d32f9abdSBarry Smith 1199d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1200d32f9abdSBarry 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 1201d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1202d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1203d32f9abdSBarry Smith 1204db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1205db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1206db4c96c1SJed Brown 12075d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12085d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12095d4c12cdSJungho Lee available when this routine is called. 12105d4c12cdSJungho Lee 1211d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12120971522cSBarry Smith 12130971522cSBarry Smith @*/ 12145d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12150971522cSBarry Smith { 12164ac538c5SBarry Smith PetscErrorCode ierr; 12170971522cSBarry Smith 12180971522cSBarry Smith PetscFunctionBegin; 12190700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1220db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1221db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1222db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12235d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12240971522cSBarry Smith PetscFunctionReturn(0); 12250971522cSBarry Smith } 12260971522cSBarry Smith 12270971522cSBarry Smith #undef __FUNCT__ 1228704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1229704ba839SBarry Smith /*@ 1230704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1231704ba839SBarry Smith 1232ad4df100SBarry Smith Logically Collective on PC 1233704ba839SBarry Smith 1234704ba839SBarry Smith Input Parameters: 1235704ba839SBarry Smith + pc - the preconditioner context 1236a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1237db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1238704ba839SBarry Smith 1239d32f9abdSBarry Smith 1240a6ffb8dbSJed Brown Notes: 1241a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1242a6ffb8dbSJed Brown 1243db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1244db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1245d32f9abdSBarry Smith 1246704ba839SBarry Smith Level: intermediate 1247704ba839SBarry Smith 1248704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1249704ba839SBarry Smith 1250704ba839SBarry Smith @*/ 12517087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1252704ba839SBarry Smith { 12534ac538c5SBarry Smith PetscErrorCode ierr; 1254704ba839SBarry Smith 1255704ba839SBarry Smith PetscFunctionBegin; 12560700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12577b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1258db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 12594ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1260704ba839SBarry Smith PetscFunctionReturn(0); 1261704ba839SBarry Smith } 1262704ba839SBarry Smith 1263704ba839SBarry Smith #undef __FUNCT__ 126457a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 126557a9adfeSMatthew G Knepley /*@ 126657a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 126757a9adfeSMatthew G Knepley 126857a9adfeSMatthew G Knepley Logically Collective on PC 126957a9adfeSMatthew G Knepley 127057a9adfeSMatthew G Knepley Input Parameters: 127157a9adfeSMatthew G Knepley + pc - the preconditioner context 127257a9adfeSMatthew G Knepley - splitname - name of this split 127357a9adfeSMatthew G Knepley 127457a9adfeSMatthew G Knepley Output Parameter: 127557a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 127657a9adfeSMatthew G Knepley 127757a9adfeSMatthew G Knepley Level: intermediate 127857a9adfeSMatthew G Knepley 127957a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 128057a9adfeSMatthew G Knepley 128157a9adfeSMatthew G Knepley @*/ 128257a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 128357a9adfeSMatthew G Knepley { 128457a9adfeSMatthew G Knepley PetscErrorCode ierr; 128557a9adfeSMatthew G Knepley 128657a9adfeSMatthew G Knepley PetscFunctionBegin; 128757a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 128857a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 128957a9adfeSMatthew G Knepley PetscValidPointer(is,3); 129057a9adfeSMatthew G Knepley { 129157a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 129257a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 129357a9adfeSMatthew G Knepley PetscBool found; 129457a9adfeSMatthew G Knepley 129557a9adfeSMatthew G Knepley *is = PETSC_NULL; 129657a9adfeSMatthew G Knepley while(ilink) { 129757a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 129857a9adfeSMatthew G Knepley if (found) { 129957a9adfeSMatthew G Knepley *is = ilink->is; 130057a9adfeSMatthew G Knepley break; 130157a9adfeSMatthew G Knepley } 130257a9adfeSMatthew G Knepley ilink = ilink->next; 130357a9adfeSMatthew G Knepley } 130457a9adfeSMatthew G Knepley } 130557a9adfeSMatthew G Knepley PetscFunctionReturn(0); 130657a9adfeSMatthew G Knepley } 130757a9adfeSMatthew G Knepley 130857a9adfeSMatthew G Knepley #undef __FUNCT__ 130951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 131051f519a2SBarry Smith /*@ 131151f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 131251f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 131351f519a2SBarry Smith 1314ad4df100SBarry Smith Logically Collective on PC 131551f519a2SBarry Smith 131651f519a2SBarry Smith Input Parameters: 131751f519a2SBarry Smith + pc - the preconditioner context 131851f519a2SBarry Smith - bs - the block size 131951f519a2SBarry Smith 132051f519a2SBarry Smith Level: intermediate 132151f519a2SBarry Smith 132251f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 132351f519a2SBarry Smith 132451f519a2SBarry Smith @*/ 13257087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 132651f519a2SBarry Smith { 13274ac538c5SBarry Smith PetscErrorCode ierr; 132851f519a2SBarry Smith 132951f519a2SBarry Smith PetscFunctionBegin; 13300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1331c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13324ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 133351f519a2SBarry Smith PetscFunctionReturn(0); 133451f519a2SBarry Smith } 133551f519a2SBarry Smith 133651f519a2SBarry Smith #undef __FUNCT__ 133769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 13380971522cSBarry Smith /*@C 133969a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 13400971522cSBarry Smith 134169a612a9SBarry Smith Collective on KSP 13420971522cSBarry Smith 13430971522cSBarry Smith Input Parameter: 13440971522cSBarry Smith . pc - the preconditioner context 13450971522cSBarry Smith 13460971522cSBarry Smith Output Parameters: 134713e0d083SBarry Smith + n - the number of splits 134869a612a9SBarry Smith - pc - the array of KSP contexts 13490971522cSBarry Smith 13500971522cSBarry Smith Note: 1351d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1352d32f9abdSBarry Smith (not the KSP just the array that contains them). 13530971522cSBarry Smith 135469a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 13550971522cSBarry Smith 13560971522cSBarry Smith Level: advanced 13570971522cSBarry Smith 13580971522cSBarry Smith .seealso: PCFIELDSPLIT 13590971522cSBarry Smith @*/ 13607087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 13610971522cSBarry Smith { 13624ac538c5SBarry Smith PetscErrorCode ierr; 13630971522cSBarry Smith 13640971522cSBarry Smith PetscFunctionBegin; 13650700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 136613e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 13674ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 13680971522cSBarry Smith PetscFunctionReturn(0); 13690971522cSBarry Smith } 13700971522cSBarry Smith 1371e69d4d44SBarry Smith #undef __FUNCT__ 1372e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1373e69d4d44SBarry Smith /*@ 1374e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1375a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1376e69d4d44SBarry Smith 1377e69d4d44SBarry Smith Collective on PC 1378e69d4d44SBarry Smith 1379e69d4d44SBarry Smith Input Parameters: 1380e69d4d44SBarry Smith + pc - the preconditioner context 1381fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1382084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1383084e4875SJed Brown 1384e69d4d44SBarry Smith Options Database: 1385084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1386e69d4d44SBarry Smith 1387fd1303e9SJungho Lee Notes: 1388fd1303e9SJungho Lee $ If ptype is 1389fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1390fd1303e9SJungho Lee $ to this function). 1391fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1392fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1393fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1394fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1395fd1303e9SJungho Lee $ preconditioner 1396fd1303e9SJungho Lee 1397fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1398fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1399fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1400fd1303e9SJungho Lee 1401fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1402fd1303e9SJungho Lee the name since it sets a proceedure to use. 1403fd1303e9SJungho Lee 1404e69d4d44SBarry Smith Level: intermediate 1405e69d4d44SBarry Smith 1406fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1407e69d4d44SBarry Smith 1408e69d4d44SBarry Smith @*/ 14097087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1410e69d4d44SBarry Smith { 14114ac538c5SBarry Smith PetscErrorCode ierr; 1412e69d4d44SBarry Smith 1413e69d4d44SBarry Smith PetscFunctionBegin; 14140700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14154ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1416e69d4d44SBarry Smith PetscFunctionReturn(0); 1417e69d4d44SBarry Smith } 1418e69d4d44SBarry Smith 1419e69d4d44SBarry Smith EXTERN_C_BEGIN 1420e69d4d44SBarry Smith #undef __FUNCT__ 1421e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14227087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1423e69d4d44SBarry Smith { 1424e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1425084e4875SJed Brown PetscErrorCode ierr; 1426e69d4d44SBarry Smith 1427e69d4d44SBarry Smith PetscFunctionBegin; 1428084e4875SJed Brown jac->schurpre = ptype; 1429084e4875SJed Brown if (pre) { 14306bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1431084e4875SJed Brown jac->schur_user = pre; 1432084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1433084e4875SJed Brown } 1434e69d4d44SBarry Smith PetscFunctionReturn(0); 1435e69d4d44SBarry Smith } 1436e69d4d44SBarry Smith EXTERN_C_END 1437e69d4d44SBarry Smith 143830ad9308SMatthew Knepley #undef __FUNCT__ 1439c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1440ab1df9f5SJed Brown /*@ 1441c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1442ab1df9f5SJed Brown 1443ab1df9f5SJed Brown Collective on PC 1444ab1df9f5SJed Brown 1445ab1df9f5SJed Brown Input Parameters: 1446ab1df9f5SJed Brown + pc - the preconditioner context 1447c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1448ab1df9f5SJed Brown 1449ab1df9f5SJed Brown Options Database: 1450c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1451ab1df9f5SJed Brown 1452ab1df9f5SJed Brown 1453ab1df9f5SJed Brown Level: intermediate 1454ab1df9f5SJed Brown 1455ab1df9f5SJed Brown Notes: 1456ab1df9f5SJed Brown The FULL factorization is 1457ab1df9f5SJed Brown 1458ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1459ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1460ab1df9f5SJed Brown 14616be4592eSBarry Smith where S = D - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D, 14626be4592eSBarry Smith and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES). 1463ab1df9f5SJed Brown 14646be4592eSBarry Smith If applied exactly, FULL factorization is a direct solver. The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial 14656be4592eSBarry Smith of degree 2, so KSPGMRES converges in 2 iterations. If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner 14666be4592eSBarry Smith application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice. With DIAG, 14676be4592eSBarry Smith the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations. 1468ab1df9f5SJed Brown 14696be4592eSBarry Smith For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES. Note that a flexible method like KSPFGMRES 14706be4592eSBarry Smith or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used inside a split). 1471ab1df9f5SJed Brown 1472ab1df9f5SJed Brown References: 1473ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1474ab1df9f5SJed Brown 1475ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1476ab1df9f5SJed Brown 1477ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1478ab1df9f5SJed Brown @*/ 1479c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1480ab1df9f5SJed Brown { 1481ab1df9f5SJed Brown PetscErrorCode ierr; 1482ab1df9f5SJed Brown 1483ab1df9f5SJed Brown PetscFunctionBegin; 1484ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1485c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1486ab1df9f5SJed Brown PetscFunctionReturn(0); 1487ab1df9f5SJed Brown } 1488ab1df9f5SJed Brown 1489ab1df9f5SJed Brown #undef __FUNCT__ 1490c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1491c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1492ab1df9f5SJed Brown { 1493ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1494ab1df9f5SJed Brown 1495ab1df9f5SJed Brown PetscFunctionBegin; 1496ab1df9f5SJed Brown jac->schurfactorization = ftype; 1497ab1df9f5SJed Brown PetscFunctionReturn(0); 1498ab1df9f5SJed Brown } 1499ab1df9f5SJed Brown 1500ab1df9f5SJed Brown #undef __FUNCT__ 150130ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 150230ad9308SMatthew Knepley /*@C 15038c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 150430ad9308SMatthew Knepley 150530ad9308SMatthew Knepley Collective on KSP 150630ad9308SMatthew Knepley 150730ad9308SMatthew Knepley Input Parameter: 150830ad9308SMatthew Knepley . pc - the preconditioner context 150930ad9308SMatthew Knepley 151030ad9308SMatthew Knepley Output Parameters: 1511a04f6461SBarry Smith + A00 - the (0,0) block 1512a04f6461SBarry Smith . A01 - the (0,1) block 1513a04f6461SBarry Smith . A10 - the (1,0) block 1514a04f6461SBarry Smith - A11 - the (1,1) block 151530ad9308SMatthew Knepley 151630ad9308SMatthew Knepley Level: advanced 151730ad9308SMatthew Knepley 151830ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 151930ad9308SMatthew Knepley @*/ 1520a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 152130ad9308SMatthew Knepley { 152230ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 152330ad9308SMatthew Knepley 152430ad9308SMatthew Knepley PetscFunctionBegin; 15250700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1526c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1527a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1528a04f6461SBarry Smith if (A01) *A01 = jac->B; 1529a04f6461SBarry Smith if (A10) *A10 = jac->C; 1530a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 153130ad9308SMatthew Knepley PetscFunctionReturn(0); 153230ad9308SMatthew Knepley } 153330ad9308SMatthew Knepley 153479416396SBarry Smith EXTERN_C_BEGIN 153579416396SBarry Smith #undef __FUNCT__ 153679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 15377087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 153879416396SBarry Smith { 153979416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1540e69d4d44SBarry Smith PetscErrorCode ierr; 154179416396SBarry Smith 154279416396SBarry Smith PetscFunctionBegin; 154379416396SBarry Smith jac->type = type; 15443b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 15453b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 15463b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1547e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1548e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1549c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1550e69d4d44SBarry Smith 15513b224e63SBarry Smith } else { 15523b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 15533b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1554e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 15559e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1556c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 15573b224e63SBarry Smith } 155879416396SBarry Smith PetscFunctionReturn(0); 155979416396SBarry Smith } 156079416396SBarry Smith EXTERN_C_END 156179416396SBarry Smith 156251f519a2SBarry Smith EXTERN_C_BEGIN 156351f519a2SBarry Smith #undef __FUNCT__ 156451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 15657087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 156651f519a2SBarry Smith { 156751f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 156851f519a2SBarry Smith 156951f519a2SBarry Smith PetscFunctionBegin; 157065e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 157165e19b50SBarry 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); 157251f519a2SBarry Smith jac->bs = bs; 157351f519a2SBarry Smith PetscFunctionReturn(0); 157451f519a2SBarry Smith } 157551f519a2SBarry Smith EXTERN_C_END 157651f519a2SBarry Smith 157779416396SBarry Smith #undef __FUNCT__ 157879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1579bc08b0f1SBarry Smith /*@ 158079416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 158179416396SBarry Smith 158279416396SBarry Smith Collective on PC 158379416396SBarry Smith 158479416396SBarry Smith Input Parameter: 158579416396SBarry Smith . pc - the preconditioner context 158681540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 158779416396SBarry Smith 158879416396SBarry Smith Options Database Key: 1589a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 159079416396SBarry Smith 1591b02e2d75SMatthew G Knepley Level: Intermediate 159279416396SBarry Smith 159379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 159479416396SBarry Smith 159579416396SBarry Smith .seealso: PCCompositeSetType() 159679416396SBarry Smith 159779416396SBarry Smith @*/ 15987087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 159979416396SBarry Smith { 16004ac538c5SBarry Smith PetscErrorCode ierr; 160179416396SBarry Smith 160279416396SBarry Smith PetscFunctionBegin; 16030700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16044ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 160579416396SBarry Smith PetscFunctionReturn(0); 160679416396SBarry Smith } 160779416396SBarry Smith 1608b02e2d75SMatthew G Knepley #undef __FUNCT__ 1609b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1610b02e2d75SMatthew G Knepley /*@ 1611b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1612b02e2d75SMatthew G Knepley 1613b02e2d75SMatthew G Knepley Not collective 1614b02e2d75SMatthew G Knepley 1615b02e2d75SMatthew G Knepley Input Parameter: 1616b02e2d75SMatthew G Knepley . pc - the preconditioner context 1617b02e2d75SMatthew G Knepley 1618b02e2d75SMatthew G Knepley Output Parameter: 1619b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1620b02e2d75SMatthew G Knepley 1621b02e2d75SMatthew G Knepley Level: Intermediate 1622b02e2d75SMatthew G Knepley 1623b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1624b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1625b02e2d75SMatthew G Knepley @*/ 1626b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1627b02e2d75SMatthew G Knepley { 1628b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1629b02e2d75SMatthew G Knepley 1630b02e2d75SMatthew G Knepley PetscFunctionBegin; 1631b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1632b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1633b02e2d75SMatthew G Knepley *type = jac->type; 1634b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1635b02e2d75SMatthew G Knepley } 1636b02e2d75SMatthew G Knepley 16370971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 16380971522cSBarry Smith /*MC 1639a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1640a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 16410971522cSBarry Smith 1642edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1643edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 16440971522cSBarry Smith 1645a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 164669a612a9SBarry Smith and set the options directly on the resulting KSP object 16470971522cSBarry Smith 16480971522cSBarry Smith Level: intermediate 16490971522cSBarry Smith 165079416396SBarry Smith Options Database Keys: 165181540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 165281540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 165381540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 165481540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 16550f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 16560f188ba9SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag 1657435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 165879416396SBarry Smith 16595d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 16605d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 16615d4c12cdSJungho Lee 1662c8a0d604SMatthew G Knepley Notes: 1663c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1664d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1665d32f9abdSBarry Smith 1666d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1667d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1668d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1669d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1670d32f9abdSBarry Smith 1671c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1672c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1673c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1674c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1675c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1676a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1677c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1678c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1679c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1680a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1681c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1682c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 16835668aaf4SBarry Smith diag gives 1684c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1685c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 16865668aaf4SBarry Smith note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of 1687c8a0d604SMatthew G Knepley $ ( A00 0 ) 1688c8a0d604SMatthew G Knepley $ ( A10 S ) 1689c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1690c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1691c8a0d604SMatthew G Knepley $ ( 0 S ) 1692c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1693e69d4d44SBarry Smith 1694edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1695edf189efSBarry Smith is used automatically for a second block. 1696edf189efSBarry Smith 1697ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1698ff218e97SBarry Smith Generally it should be used with the AIJ format. 1699ff218e97SBarry Smith 1700ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1701ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1702ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 17030716a85fSBarry Smith 1704a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 17050971522cSBarry Smith 17067e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1707e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 17080971522cSBarry Smith M*/ 17090971522cSBarry Smith 17100971522cSBarry Smith EXTERN_C_BEGIN 17110971522cSBarry Smith #undef __FUNCT__ 17120971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 17137087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 17140971522cSBarry Smith { 17150971522cSBarry Smith PetscErrorCode ierr; 17160971522cSBarry Smith PC_FieldSplit *jac; 17170971522cSBarry Smith 17180971522cSBarry Smith PetscFunctionBegin; 171938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 17200971522cSBarry Smith jac->bs = -1; 17210971522cSBarry Smith jac->nsplits = 0; 17223e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1723e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1724c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 172551f519a2SBarry Smith 17260971522cSBarry Smith pc->data = (void*)jac; 17270971522cSBarry Smith 17280971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1729421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 17300971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1731574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 17320971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 17330971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 17340971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 17350971522cSBarry Smith pc->ops->applyrichardson = 0; 17360971522cSBarry Smith 173769a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 173869a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 17390971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 17400971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1741704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1742704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 174379416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 174479416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 174551f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 174651f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 17470971522cSBarry Smith PetscFunctionReturn(0); 17480971522cSBarry Smith } 17490971522cSBarry Smith EXTERN_C_END 17500971522cSBarry Smith 17510971522cSBarry Smith 1752a541d17aSBarry Smith 1753