1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 33c48a1e8SJed Brown #include <petscdmcomposite.h> /*I "petscdmcomposite.h" I*/ 40971522cSBarry Smith 5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0}; 7c5d2311dSJed Brown 80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 90971522cSBarry Smith struct _PC_FieldSplitLink { 1069a612a9SBarry Smith KSP ksp; 110971522cSBarry Smith Vec x,y; 12db4c96c1SJed Brown char *splitname; 130971522cSBarry Smith PetscInt nfields; 145d4c12cdSJungho Lee PetscInt *fields,*fields_col; 151b9fc7fcSBarry Smith VecScatter sctx; 165d4c12cdSJungho Lee IS is,is_col; 1751f519a2SBarry Smith PC_FieldSplitLink next,previous; 180971522cSBarry Smith }; 190971522cSBarry Smith 200971522cSBarry Smith typedef struct { 2168dd23aaSBarry Smith PCCompositeType type; 22ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 23ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 24ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 2530ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2630ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2779416396SBarry Smith Vec *x,*y,w1,w2; 28519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 29519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3030ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 31ace3abfcSBarry Smith PetscBool issetup; 3230ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3330ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3430ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 35a04f6461SBarry Smith Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01 */ 36084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 37084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 38c9c6ffaaSJed Brown PCFieldSplitSchurFactType schurfactorization; 3930ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 4097bbdb24SBarry Smith PC_FieldSplitLink head; 4163ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 42c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 430971522cSBarry Smith } PC_FieldSplit; 440971522cSBarry Smith 4516913363SBarry Smith /* 4616913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4716913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4816913363SBarry Smith PC you could change this. 4916913363SBarry Smith */ 50084e4875SJed Brown 51e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 52084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 53084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 54084e4875SJed Brown { 55084e4875SJed Brown switch (jac->schurpre) { 56084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 57084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 58084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 59084e4875SJed Brown default: 60084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 61084e4875SJed Brown } 62084e4875SJed Brown } 63084e4875SJed Brown 64084e4875SJed Brown 650971522cSBarry Smith #undef __FUNCT__ 660971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 670971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 680971522cSBarry Smith { 690971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 700971522cSBarry Smith PetscErrorCode ierr; 71ace3abfcSBarry Smith PetscBool iascii; 720971522cSBarry Smith PetscInt i,j; 735a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 740971522cSBarry Smith 750971522cSBarry Smith PetscFunctionBegin; 76251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 770971522cSBarry Smith if (iascii) { 782c9966d7SBarry Smith if (jac->bs > 0) { 7951f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 802c9966d7SBarry Smith } else { 812c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 822c9966d7SBarry Smith } 8369a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 840971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 850971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 861ab39975SBarry Smith if (ilink->fields) { 870971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 8879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 895a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9079416396SBarry Smith if (j > 0) { 9179416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9279416396SBarry Smith } 935a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 940971522cSBarry Smith } 950971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9679416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 971ab39975SBarry Smith } else { 981ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 991ab39975SBarry Smith } 1005a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1015a9f2f41SSatish Balay ilink = ilink->next; 1020971522cSBarry Smith } 1030971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1040971522cSBarry Smith } else { 10565e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1060971522cSBarry Smith } 1070971522cSBarry Smith PetscFunctionReturn(0); 1080971522cSBarry Smith } 1090971522cSBarry Smith 1100971522cSBarry Smith #undef __FUNCT__ 1113b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1123b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1133b224e63SBarry Smith { 1143b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1153b224e63SBarry Smith PetscErrorCode ierr; 116ace3abfcSBarry Smith PetscBool iascii; 1173b224e63SBarry Smith PetscInt i,j; 1183b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1193b224e63SBarry Smith KSP ksp; 1203b224e63SBarry Smith 1213b224e63SBarry Smith PetscFunctionBegin; 122251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1233b224e63SBarry Smith if (iascii) { 1242c9966d7SBarry Smith if (jac->bs > 0) { 125c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1262c9966d7SBarry Smith } else { 1272c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1282c9966d7SBarry Smith } 129*3e8b8b31SMatthew G Knepley switch(jac->schurpre) { 130*3e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 131*3e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 132*3e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 133*3e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break; 134*3e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 135*3e8b8b31SMatthew G Knepley if (jac->schur_user) { 136*3e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 137*3e8b8b31SMatthew G Knepley } else { 138*3e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr); 139*3e8b8b31SMatthew G Knepley } 140*3e8b8b31SMatthew G Knepley break; 141*3e8b8b31SMatthew G Knepley default: 142*3e8b8b31SMatthew G Knepley SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 143*3e8b8b31SMatthew G Knepley } 1443b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1453b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1463b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1473b224e63SBarry Smith if (ilink->fields) { 1483b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1493b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1503b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1513b224e63SBarry Smith if (j > 0) { 1523b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1533b224e63SBarry Smith } 1543b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1553b224e63SBarry Smith } 1563b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1573b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1583b224e63SBarry Smith } else { 1593b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1603b224e63SBarry Smith } 1613b224e63SBarry Smith ilink = ilink->next; 1623b224e63SBarry Smith } 163435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1643b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 16512cae6f2SJed Brown if (jac->schur) { 16612cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1673b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 16812cae6f2SJed Brown } else { 16912cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 17012cae6f2SJed Brown } 1713b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 172435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1733b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 17412cae6f2SJed Brown if (jac->kspschur) { 1753b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 17612cae6f2SJed Brown } else { 17712cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 17812cae6f2SJed Brown } 1793b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1803b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1813b224e63SBarry Smith } else { 18265e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1833b224e63SBarry Smith } 1843b224e63SBarry Smith PetscFunctionReturn(0); 1853b224e63SBarry Smith } 1863b224e63SBarry Smith 1873b224e63SBarry Smith #undef __FUNCT__ 1886c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1896c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1906c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1916c924f48SJed Brown { 1926c924f48SJed Brown PetscErrorCode ierr; 1936c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1945d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 1955d4c12cdSJungho Lee PetscBool flg,flg_col; 1965d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 1976c924f48SJed Brown 1986c924f48SJed Brown PetscFunctionBegin; 1996c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2005d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2016c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2026c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2036c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2045d4c12cdSJungho Lee ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2056c924f48SJed Brown nfields = jac->bs; 20629499fbbSJungho Lee nfields_col = jac->bs; 2076c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2085d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2096c924f48SJed Brown if (!flg) break; 2105d4c12cdSJungho Lee else if (flg && !flg_col) { 2115d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2125d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2135d4c12cdSJungho Lee } 2145d4c12cdSJungho Lee else { 2155d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2165d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2175d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2185d4c12cdSJungho Lee } 2196c924f48SJed Brown } 2206c924f48SJed Brown if (i > 0) { 2216c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2226c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2236c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2246c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2256c924f48SJed Brown } 2266c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2275d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2286c924f48SJed Brown PetscFunctionReturn(0); 2296c924f48SJed Brown } 2306c924f48SJed Brown 2316c924f48SJed Brown #undef __FUNCT__ 23269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 23369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2340971522cSBarry Smith { 2350971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2360971522cSBarry Smith PetscErrorCode ierr; 2375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2386ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2396c924f48SJed Brown PetscInt i; 2400971522cSBarry Smith 2410971522cSBarry Smith PetscFunctionBegin; 242d32f9abdSBarry Smith if (!ilink) { 243d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 244d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2457b62db95SJungho Lee PetscInt numFields, f; 2460784a22cSJed Brown char **fieldNames; 2477b62db95SJungho Lee IS *fields; 248e7c4fc90SDmitry Karpeev DM *dms; 249e7c4fc90SDmitry Karpeev ierr = DMCreateDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms); CHKERRQ(ierr); 2507b62db95SJungho Lee for(f = 0; f < numFields; ++f) { 2517b62db95SJungho Lee ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 2527b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 2537b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 2547b62db95SJungho Lee } 2557b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 2567b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 257e7c4fc90SDmitry Karpeev if(dms) { 2588b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2597b62db95SJungho Lee for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) { 2607b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 2617b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 262e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 2632fa5ba8aSJed Brown } 2647b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 2658b8307b2SJed Brown } 26666ffff09SJed Brown } else { 267521d7252SBarry Smith if (jac->bs <= 0) { 268704ba839SBarry Smith if (pc->pmat) { 269521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 270704ba839SBarry Smith } else { 271704ba839SBarry Smith jac->bs = 1; 272704ba839SBarry Smith } 273521d7252SBarry Smith } 274d32f9abdSBarry Smith 275acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2766ce1633cSBarry Smith if (stokes) { 2776ce1633cSBarry Smith IS zerodiags,rest; 2786ce1633cSBarry Smith PetscInt nmin,nmax; 2796ce1633cSBarry Smith 2806ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2816ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2826ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2836ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2846ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 285fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 286fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2876ce1633cSBarry Smith } else { 288d32f9abdSBarry Smith if (!flg) { 289d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 290d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2916c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2926c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 293d32f9abdSBarry Smith } 2946c924f48SJed Brown if (flg || !jac->splitdefined) { 295d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 296db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2976c924f48SJed Brown char splitname[8]; 2986c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2995d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 30079416396SBarry Smith } 3015d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 302521d7252SBarry Smith } 30366ffff09SJed Brown } 3046ce1633cSBarry Smith } 305edf189efSBarry Smith } else if (jac->nsplits == 1) { 306edf189efSBarry Smith if (ilink->is) { 307edf189efSBarry Smith IS is2; 308edf189efSBarry Smith PetscInt nmin,nmax; 309edf189efSBarry Smith 310edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 311edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 312db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 313fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 3147b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 31563ec74ffSBarry Smith } else if (jac->reset) { 31663ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 317d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 318d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 319d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 320d0af7cd3SBarry Smith if (pc->dm && !stokes) { 321d0af7cd3SBarry Smith PetscBool dmcomposite; 322251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 323d0af7cd3SBarry Smith if (dmcomposite) { 324d0af7cd3SBarry Smith PetscInt nDM; 325d0af7cd3SBarry Smith IS *fields; 3267b62db95SJungho Lee DM *dms; 327d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 328d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 329d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 3307b62db95SJungho Lee ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr); 3317b62db95SJungho Lee ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr); 332d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 3337b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr); 3347b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr); 335d0af7cd3SBarry Smith ilink->is = fields[i]; 336d0af7cd3SBarry Smith ilink = ilink->next; 337edf189efSBarry Smith } 338d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 3397b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 340d0af7cd3SBarry Smith } 341d0af7cd3SBarry Smith } else { 342d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 343d0af7cd3SBarry Smith if (stokes) { 344d0af7cd3SBarry Smith IS zerodiags,rest; 345d0af7cd3SBarry Smith PetscInt nmin,nmax; 346d0af7cd3SBarry Smith 347d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 348d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 349d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 35020b26d62SBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 35120b26d62SBarry Smith ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr); 352d0af7cd3SBarry Smith ilink->is = rest; 353d0af7cd3SBarry Smith ilink->next->is = zerodiags; 35420b26d62SBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 355d0af7cd3SBarry Smith } 356d0af7cd3SBarry Smith } 357d0af7cd3SBarry Smith 3587b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 35969a612a9SBarry Smith PetscFunctionReturn(0); 36069a612a9SBarry Smith } 36169a612a9SBarry Smith 36269a612a9SBarry Smith #undef __FUNCT__ 36369a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 36469a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 36569a612a9SBarry Smith { 36669a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 36769a612a9SBarry Smith PetscErrorCode ierr; 3685a9f2f41SSatish Balay PC_FieldSplitLink ilink; 3692c9966d7SBarry Smith PetscInt i,nsplit; 37069a612a9SBarry Smith MatStructure flag = pc->flag; 3715f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 37269a612a9SBarry Smith 37369a612a9SBarry Smith PetscFunctionBegin; 37469a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 37597bbdb24SBarry Smith nsplit = jac->nsplits; 3765a9f2f41SSatish Balay ilink = jac->head; 37797bbdb24SBarry Smith 37897bbdb24SBarry Smith /* get the matrices for each split */ 379704ba839SBarry Smith if (!jac->issetup) { 3801b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 38197bbdb24SBarry Smith 382704ba839SBarry Smith jac->issetup = PETSC_TRUE; 383704ba839SBarry Smith 3845d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 3852c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 3862c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 3872c9966d7SBarry Smith } 38851f519a2SBarry Smith bs = jac->bs; 38997bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 3901b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3911b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3921b9fc7fcSBarry Smith if (jac->defaultsplit) { 393704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 3945f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 395704ba839SBarry Smith } else if (!ilink->is) { 396ccb205f8SBarry Smith if (ilink->nfields > 1) { 3975f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 3985a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3995f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4001b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4011b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4021b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4035f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 40497bbdb24SBarry Smith } 40597bbdb24SBarry Smith } 406d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4075f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 408ccb205f8SBarry Smith } else { 409704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4105f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 411ccb205f8SBarry Smith } 4123e197d65SBarry Smith } 413704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4145f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4155f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 416704ba839SBarry Smith ilink = ilink->next; 4171b9fc7fcSBarry Smith } 4181b9fc7fcSBarry Smith } 4191b9fc7fcSBarry Smith 420704ba839SBarry Smith ilink = jac->head; 42197bbdb24SBarry Smith if (!jac->pmat) { 422bdddcaaaSMatthew G Knepley Vec xtmp; 423bdddcaaaSMatthew G Knepley 424bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 425cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 426bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 427cf502942SBarry Smith for (i=0; i<nsplit; i++) { 428bdddcaaaSMatthew G Knepley MatNullSpace sp; 429bdddcaaaSMatthew G Knepley 4305f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 431bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 432bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 433bdddcaaaSMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; 434bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 435bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 436bdddcaaaSMatthew G Knepley /* HACK: Check for the constant null space */ 437bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(pc->pmat, &sp);CHKERRQ(ierr); 438bdddcaaaSMatthew G Knepley if (sp) { 439bdddcaaaSMatthew G Knepley MatNullSpace subsp; 440bdddcaaaSMatthew G Knepley Vec ftmp, gtmp; 4419583d628SBarry Smith PetscReal norm; 442bdddcaaaSMatthew G Knepley PetscInt N; 443bdddcaaaSMatthew G Knepley 444bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat, >mp, PETSC_NULL);CHKERRQ(ierr); 445bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i], &ftmp, PETSC_NULL);CHKERRQ(ierr); 446bdddcaaaSMatthew G Knepley ierr = VecGetSize(ftmp, &N);CHKERRQ(ierr); 447bdddcaaaSMatthew G Knepley ierr = VecSet(ftmp, 1.0/N);CHKERRQ(ierr); 448bdddcaaaSMatthew G Knepley ierr = VecScatterBegin(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 449bdddcaaaSMatthew G Knepley ierr = VecScatterEnd(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 450bdddcaaaSMatthew G Knepley ierr = MatNullSpaceRemove(sp, gtmp, PETSC_NULL);CHKERRQ(ierr); 451bdddcaaaSMatthew G Knepley ierr = VecNorm(gtmp, NORM_2, &norm);CHKERRQ(ierr); 452bdddcaaaSMatthew G Knepley if (norm < 1.0e-10) { 453bdddcaaaSMatthew G Knepley ierr = MatNullSpaceCreate(((PetscObject)pc)->comm, PETSC_TRUE, 0, PETSC_NULL, &subsp);CHKERRQ(ierr); 454bdddcaaaSMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], subsp);CHKERRQ(ierr); 455bdddcaaaSMatthew G Knepley ierr = MatNullSpaceDestroy(&subsp);CHKERRQ(ierr); 456bdddcaaaSMatthew G Knepley } 457bdddcaaaSMatthew G Knepley ierr = VecDestroy(&ftmp);CHKERRQ(ierr); 458bdddcaaaSMatthew G Knepley ierr = VecDestroy(>mp);CHKERRQ(ierr); 459bdddcaaaSMatthew G Knepley } 460ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 461ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 462ed1f0337SMatthew G Knepley if (sp) { 463ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 464ed1f0337SMatthew G Knepley } 465704ba839SBarry Smith ilink = ilink->next; 466cf502942SBarry Smith } 467bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 46897bbdb24SBarry Smith } else { 469cf502942SBarry Smith for (i=0; i<nsplit; i++) { 4705f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 471704ba839SBarry Smith ilink = ilink->next; 472cf502942SBarry Smith } 47397bbdb24SBarry Smith } 474519d70e2SJed Brown if (jac->realdiagonal) { 475519d70e2SJed Brown ilink = jac->head; 476519d70e2SJed Brown if (!jac->mat) { 477519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 478519d70e2SJed Brown for (i=0; i<nsplit; i++) { 4795f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 480519d70e2SJed Brown ilink = ilink->next; 481519d70e2SJed Brown } 482519d70e2SJed Brown } else { 483519d70e2SJed Brown for (i=0; i<nsplit; i++) { 4845f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 485519d70e2SJed Brown ilink = ilink->next; 486519d70e2SJed Brown } 487519d70e2SJed Brown } 488519d70e2SJed Brown } else { 489519d70e2SJed Brown jac->mat = jac->pmat; 490519d70e2SJed Brown } 49197bbdb24SBarry Smith 4926c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 49368dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 49468dd23aaSBarry Smith ilink = jac->head; 49568dd23aaSBarry Smith if (!jac->Afield) { 49668dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 49768dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4984aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 49968dd23aaSBarry Smith ilink = ilink->next; 50068dd23aaSBarry Smith } 50168dd23aaSBarry Smith } else { 50268dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5034aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 50468dd23aaSBarry Smith ilink = ilink->next; 50568dd23aaSBarry Smith } 50668dd23aaSBarry Smith } 50768dd23aaSBarry Smith } 50868dd23aaSBarry Smith 5093b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5103b224e63SBarry Smith IS ccis; 5114aa3045dSJed Brown PetscInt rstart,rend; 512093c86ffSJed Brown char lscname[256]; 513093c86ffSJed Brown PetscObject LSC_L; 514e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 51568dd23aaSBarry Smith 516e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 517e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 518e6cab6aaSJed Brown 5193b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5203b224e63SBarry Smith if (jac->schur) { 5213b224e63SBarry Smith ilink = jac->head; 52249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5234aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 524fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5253b224e63SBarry Smith ilink = ilink->next; 52649bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5274aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 528fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 529519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 530084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 5313b224e63SBarry Smith 5323b224e63SBarry Smith } else { 5331cee3971SBarry Smith KSP ksp; 5346c924f48SJed Brown char schurprefix[256]; 535bdddcaaaSMatthew G Knepley MatNullSpace sp; 5363b224e63SBarry Smith 537a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 5383b224e63SBarry Smith ilink = jac->head; 53949bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5404aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 541fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5423b224e63SBarry Smith ilink = ilink->next; 54349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5444aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 545fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5467d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 547519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 548bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 549bdddcaaaSMatthew G Knepley if (sp) {ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);} 550a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 5511cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 552e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 553a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 554a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 55520b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 55620b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 5577b62db95SJungho Lee ierr = MatSetUp(jac->schur);CHKERRQ(ierr); 5583b224e63SBarry Smith 5593b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 5609005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 5611cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 562084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 563084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 564084e4875SJed Brown PC pc; 565cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 566084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 567084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 568e69d4d44SBarry Smith } 5696c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 5706c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 5713b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 57220b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 57320b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 5743b224e63SBarry Smith } 575093c86ffSJed Brown 576093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 577093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 578093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 579093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 580093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 581093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 582093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 583093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 584093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 5853b224e63SBarry Smith } else { 58697bbdb24SBarry Smith /* set up the individual PCs */ 58797bbdb24SBarry Smith i = 0; 5885a9f2f41SSatish Balay ilink = jac->head; 5895a9f2f41SSatish Balay while (ilink) { 590519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 5913b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 592c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 5935a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 59497bbdb24SBarry Smith i++; 5955a9f2f41SSatish Balay ilink = ilink->next; 5960971522cSBarry Smith } 5973b224e63SBarry Smith } 5983b224e63SBarry Smith 599c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 6000971522cSBarry Smith PetscFunctionReturn(0); 6010971522cSBarry Smith } 6020971522cSBarry Smith 6035a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 604ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 605ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 6065a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 607ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 608ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 60979416396SBarry Smith 6100971522cSBarry Smith #undef __FUNCT__ 6113b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 6123b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 6133b224e63SBarry Smith { 6143b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6153b224e63SBarry Smith PetscErrorCode ierr; 6163b224e63SBarry Smith KSP ksp; 6173b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 6183b224e63SBarry Smith 6193b224e63SBarry Smith PetscFunctionBegin; 6203b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 6213b224e63SBarry Smith 622c5d2311dSJed Brown switch (jac->schurfactorization) { 623c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 624a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 625c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 626c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 627c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 628c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 629c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 630c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 631c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 632c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 633c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 634c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 635c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 636c5d2311dSJed Brown break; 637c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 638a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 639c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 640c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 641c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 642c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 643c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 644c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 645c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 646c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 647c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 648c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 649c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 650c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 651c5d2311dSJed Brown break; 652c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 653a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 654c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 655c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 656c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 657c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 658c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 659c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 660c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 661c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 662c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 663c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 664c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 665c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 666c5d2311dSJed Brown break; 667c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 668a04f6461SBarry 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 */ 6693b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6703b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6713b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6723b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6733b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6743b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6753b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6763b224e63SBarry Smith 6773b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6783b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6793b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6803b224e63SBarry Smith 6813b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6823b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6833b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6843b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6853b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 686c5d2311dSJed Brown } 6873b224e63SBarry Smith PetscFunctionReturn(0); 6883b224e63SBarry Smith } 6893b224e63SBarry Smith 6903b224e63SBarry Smith #undef __FUNCT__ 6910971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6920971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6930971522cSBarry Smith { 6940971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6950971522cSBarry Smith PetscErrorCode ierr; 6965a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 697939b8a20SBarry Smith PetscInt cnt,bs; 6980971522cSBarry Smith 6990971522cSBarry Smith PetscFunctionBegin; 70051f519a2SBarry Smith CHKMEMQ; 70179416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 7021b9fc7fcSBarry Smith if (jac->defaultsplit) { 703939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 704939b8a20SBarry 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); 705939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 706939b8a20SBarry 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); 7070971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 7085a9f2f41SSatish Balay while (ilink) { 7095a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7105a9f2f41SSatish Balay ilink = ilink->next; 7110971522cSBarry Smith } 7120971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 7131b9fc7fcSBarry Smith } else { 714efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7155a9f2f41SSatish Balay while (ilink) { 7165a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 7175a9f2f41SSatish Balay ilink = ilink->next; 7181b9fc7fcSBarry Smith } 7191b9fc7fcSBarry Smith } 72016913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 72179416396SBarry Smith if (!jac->w1) { 72279416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 72379416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 72479416396SBarry Smith } 725efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7265a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 7273e197d65SBarry Smith cnt = 1; 7285a9f2f41SSatish Balay while (ilink->next) { 7295a9f2f41SSatish Balay ilink = ilink->next; 7303e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 7313e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 7323e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 7333e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7343e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7353e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7363e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7373e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7383e197d65SBarry Smith } 73951f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 74011755939SBarry Smith cnt -= 2; 74151f519a2SBarry Smith while (ilink->previous) { 74251f519a2SBarry Smith ilink = ilink->previous; 74311755939SBarry Smith /* compute the residual only over the part of the vector needed */ 74411755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 74511755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 74611755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 74711755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 74811755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 74911755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 75011755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 75151f519a2SBarry Smith } 75211755939SBarry Smith } 75365e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 75451f519a2SBarry Smith CHKMEMQ; 7550971522cSBarry Smith PetscFunctionReturn(0); 7560971522cSBarry Smith } 7570971522cSBarry Smith 758421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 759ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 760ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 761421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 762ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 763ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 764421e10b8SBarry Smith 765421e10b8SBarry Smith #undef __FUNCT__ 7668c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 767421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 768421e10b8SBarry Smith { 769421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 770421e10b8SBarry Smith PetscErrorCode ierr; 771421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 772939b8a20SBarry Smith PetscInt bs; 773421e10b8SBarry Smith 774421e10b8SBarry Smith PetscFunctionBegin; 775421e10b8SBarry Smith CHKMEMQ; 776421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 777421e10b8SBarry Smith if (jac->defaultsplit) { 778939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 779939b8a20SBarry 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); 780939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 781939b8a20SBarry 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); 782421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 783421e10b8SBarry Smith while (ilink) { 784421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 785421e10b8SBarry Smith ilink = ilink->next; 786421e10b8SBarry Smith } 787421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 788421e10b8SBarry Smith } else { 789421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 790421e10b8SBarry Smith while (ilink) { 791421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 792421e10b8SBarry Smith ilink = ilink->next; 793421e10b8SBarry Smith } 794421e10b8SBarry Smith } 795421e10b8SBarry Smith } else { 796421e10b8SBarry Smith if (!jac->w1) { 797421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 798421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 799421e10b8SBarry Smith } 800421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 801421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 802421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 803421e10b8SBarry Smith while (ilink->next) { 804421e10b8SBarry Smith ilink = ilink->next; 8059989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 806421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 807421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 808421e10b8SBarry Smith } 809421e10b8SBarry Smith while (ilink->previous) { 810421e10b8SBarry Smith ilink = ilink->previous; 8119989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 812421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 813421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 814421e10b8SBarry Smith } 815421e10b8SBarry Smith } else { 816421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 817421e10b8SBarry Smith ilink = ilink->next; 818421e10b8SBarry Smith } 819421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 820421e10b8SBarry Smith while (ilink->previous) { 821421e10b8SBarry Smith ilink = ilink->previous; 8229989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 823421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 824421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 825421e10b8SBarry Smith } 826421e10b8SBarry Smith } 827421e10b8SBarry Smith } 828421e10b8SBarry Smith CHKMEMQ; 829421e10b8SBarry Smith PetscFunctionReturn(0); 830421e10b8SBarry Smith } 831421e10b8SBarry Smith 8320971522cSBarry Smith #undef __FUNCT__ 833574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 834574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 8350971522cSBarry Smith { 8360971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8370971522cSBarry Smith PetscErrorCode ierr; 8385a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 8390971522cSBarry Smith 8400971522cSBarry Smith PetscFunctionBegin; 8415a9f2f41SSatish Balay while (ilink) { 842574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 843fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 844fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 845fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 846fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 847b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 8485a9f2f41SSatish Balay next = ilink->next; 8495a9f2f41SSatish Balay ilink = next; 8500971522cSBarry Smith } 85105b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 852574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 853574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 854574deadeSBarry Smith } else if (jac->mat) { 855574deadeSBarry Smith jac->mat = PETSC_NULL; 856574deadeSBarry Smith } 85797bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 85868dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8596bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8606bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8616bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8626bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8636bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8646bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8656bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 86663ec74ffSBarry Smith jac->reset = PETSC_TRUE; 867574deadeSBarry Smith PetscFunctionReturn(0); 868574deadeSBarry Smith } 869574deadeSBarry Smith 870574deadeSBarry Smith #undef __FUNCT__ 871574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 872574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 873574deadeSBarry Smith { 874574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 875574deadeSBarry Smith PetscErrorCode ierr; 876574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 877574deadeSBarry Smith 878574deadeSBarry Smith PetscFunctionBegin; 879574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 880574deadeSBarry Smith while (ilink) { 8816bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 882574deadeSBarry Smith next = ilink->next; 883574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 884574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 8855d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 886574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 887574deadeSBarry Smith ilink = next; 888574deadeSBarry Smith } 889574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 890c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8917b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 8927b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 8937b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 8947b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 8957b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 896ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 897c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 8980971522cSBarry Smith PetscFunctionReturn(0); 8990971522cSBarry Smith } 9000971522cSBarry Smith 9010971522cSBarry Smith #undef __FUNCT__ 9020971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 9030971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 9040971522cSBarry Smith { 9051b9fc7fcSBarry Smith PetscErrorCode ierr; 9066c924f48SJed Brown PetscInt bs; 907bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 9089dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9093b224e63SBarry Smith PCCompositeType ctype; 910e7c4fc90SDmitry Karpeev DM ddm; 911e7c4fc90SDmitry Karpeev char ddm_name[1025]; 9121b9fc7fcSBarry Smith 9130971522cSBarry Smith PetscFunctionBegin; 9141b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 915e7c4fc90SDmitry Karpeev if(pc->dm) { 916e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 917e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 918e7c4fc90SDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition_dm", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 919e7c4fc90SDmitry Karpeev if(flg) { 920e7c4fc90SDmitry Karpeev ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 921e7c4fc90SDmitry Karpeev if(!ddm) { 922e7c4fc90SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name); 923e7c4fc90SDmitry Karpeev } 924e7c4fc90SDmitry Karpeev ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr); 925e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 926e7c4fc90SDmitry Karpeev } 927e7c4fc90SDmitry Karpeev } 928acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 92951f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 93051f519a2SBarry Smith if (flg) { 93151f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 93251f519a2SBarry Smith } 933704ba839SBarry Smith 934435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 935c0adfefeSBarry Smith if (stokes) { 936c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 937c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 938c0adfefeSBarry Smith } 939c0adfefeSBarry Smith 9403b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 9413b224e63SBarry Smith if (flg) { 9423b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 9433b224e63SBarry Smith } 944c30613efSMatthew Knepley /* Only setup fields once */ 945c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 946d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 947d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 9486c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 9496c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 950d32f9abdSBarry Smith } 951c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 952c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 953c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 954c9c6ffaaSJed 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); 955084e4875SJed 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); 956c5d2311dSJed Brown } 9571b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 9580971522cSBarry Smith PetscFunctionReturn(0); 9590971522cSBarry Smith } 9600971522cSBarry Smith 9610971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 9620971522cSBarry Smith 9630971522cSBarry Smith EXTERN_C_BEGIN 9640971522cSBarry Smith #undef __FUNCT__ 9650971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 9665d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 9670971522cSBarry Smith { 96897bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9690971522cSBarry Smith PetscErrorCode ierr; 9705a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 97169a612a9SBarry Smith char prefix[128]; 9725d4c12cdSJungho Lee PetscInt i; 9730971522cSBarry Smith 9740971522cSBarry Smith PetscFunctionBegin; 9756c924f48SJed Brown if (jac->splitdefined) { 9766c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9776c924f48SJed Brown PetscFunctionReturn(0); 9786c924f48SJed Brown } 97951f519a2SBarry Smith for (i=0; i<n; i++) { 980e32f2f54SBarry 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); 981e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 98251f519a2SBarry Smith } 983704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 984a04f6461SBarry Smith if (splitname) { 985db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 986a04f6461SBarry Smith } else { 987a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 988a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 989a04f6461SBarry Smith } 990704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 9915a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 9925d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 9935d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 9945a9f2f41SSatish Balay ilink->nfields = n; 9955a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9967adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9971cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9985a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9999005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 100069a612a9SBarry Smith 1001a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 10025a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 10030971522cSBarry Smith 10040971522cSBarry Smith if (!next) { 10055a9f2f41SSatish Balay jac->head = ilink; 100651f519a2SBarry Smith ilink->previous = PETSC_NULL; 10070971522cSBarry Smith } else { 10080971522cSBarry Smith while (next->next) { 10090971522cSBarry Smith next = next->next; 10100971522cSBarry Smith } 10115a9f2f41SSatish Balay next->next = ilink; 101251f519a2SBarry Smith ilink->previous = next; 10130971522cSBarry Smith } 10140971522cSBarry Smith jac->nsplits++; 10150971522cSBarry Smith PetscFunctionReturn(0); 10160971522cSBarry Smith } 10170971522cSBarry Smith EXTERN_C_END 10180971522cSBarry Smith 1019e69d4d44SBarry Smith EXTERN_C_BEGIN 1020e69d4d44SBarry Smith #undef __FUNCT__ 1021e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 10227087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1023e69d4d44SBarry Smith { 1024e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1025e69d4d44SBarry Smith PetscErrorCode ierr; 1026e69d4d44SBarry Smith 1027e69d4d44SBarry Smith PetscFunctionBegin; 1028519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1029e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1030e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 103113e0d083SBarry Smith if (n) *n = jac->nsplits; 1032e69d4d44SBarry Smith PetscFunctionReturn(0); 1033e69d4d44SBarry Smith } 1034e69d4d44SBarry Smith EXTERN_C_END 10350971522cSBarry Smith 10360971522cSBarry Smith EXTERN_C_BEGIN 10370971522cSBarry Smith #undef __FUNCT__ 103869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 10397087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 10400971522cSBarry Smith { 10410971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10420971522cSBarry Smith PetscErrorCode ierr; 10430971522cSBarry Smith PetscInt cnt = 0; 10445a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 10450971522cSBarry Smith 10460971522cSBarry Smith PetscFunctionBegin; 10475d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 10485a9f2f41SSatish Balay while (ilink) { 10495a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 10505a9f2f41SSatish Balay ilink = ilink->next; 10510971522cSBarry Smith } 10525d480477SMatthew 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); 105313e0d083SBarry Smith if (n) *n = jac->nsplits; 10540971522cSBarry Smith PetscFunctionReturn(0); 10550971522cSBarry Smith } 10560971522cSBarry Smith EXTERN_C_END 10570971522cSBarry Smith 1058704ba839SBarry Smith EXTERN_C_BEGIN 1059704ba839SBarry Smith #undef __FUNCT__ 1060704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 10617087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1062704ba839SBarry Smith { 1063704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1064704ba839SBarry Smith PetscErrorCode ierr; 1065704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1066704ba839SBarry Smith char prefix[128]; 1067704ba839SBarry Smith 1068704ba839SBarry Smith PetscFunctionBegin; 10696c924f48SJed Brown if (jac->splitdefined) { 10706c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10716c924f48SJed Brown PetscFunctionReturn(0); 10726c924f48SJed Brown } 107316913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1074a04f6461SBarry Smith if (splitname) { 1075db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1076a04f6461SBarry Smith } else { 1077b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1078b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1079a04f6461SBarry Smith } 10801ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1081b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1082b5787286SJed Brown ilink->is = is; 1083b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1084b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1085b5787286SJed Brown ilink->is_col = is; 1086704ba839SBarry Smith ilink->next = PETSC_NULL; 1087704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10881cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1089704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10909005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1091704ba839SBarry Smith 1092a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1093704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1094704ba839SBarry Smith 1095704ba839SBarry Smith if (!next) { 1096704ba839SBarry Smith jac->head = ilink; 1097704ba839SBarry Smith ilink->previous = PETSC_NULL; 1098704ba839SBarry Smith } else { 1099704ba839SBarry Smith while (next->next) { 1100704ba839SBarry Smith next = next->next; 1101704ba839SBarry Smith } 1102704ba839SBarry Smith next->next = ilink; 1103704ba839SBarry Smith ilink->previous = next; 1104704ba839SBarry Smith } 1105704ba839SBarry Smith jac->nsplits++; 1106704ba839SBarry Smith 1107704ba839SBarry Smith PetscFunctionReturn(0); 1108704ba839SBarry Smith } 1109704ba839SBarry Smith EXTERN_C_END 1110704ba839SBarry Smith 11110971522cSBarry Smith #undef __FUNCT__ 11120971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 11130971522cSBarry Smith /*@ 11140971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 11150971522cSBarry Smith 1116ad4df100SBarry Smith Logically Collective on PC 11170971522cSBarry Smith 11180971522cSBarry Smith Input Parameters: 11190971522cSBarry Smith + pc - the preconditioner context 1120a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 11210971522cSBarry Smith . n - the number of fields in this split 1122db4c96c1SJed Brown - fields - the fields in this split 11230971522cSBarry Smith 11240971522cSBarry Smith Level: intermediate 11250971522cSBarry Smith 1126d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1127d32f9abdSBarry Smith 1128d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1129d32f9abdSBarry 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 1130d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1131d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1132d32f9abdSBarry Smith 1133db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1134db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1135db4c96c1SJed Brown 11365d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 11375d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 11385d4c12cdSJungho Lee available when this routine is called. 11395d4c12cdSJungho Lee 1140d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 11410971522cSBarry Smith 11420971522cSBarry Smith @*/ 11435d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11440971522cSBarry Smith { 11454ac538c5SBarry Smith PetscErrorCode ierr; 11460971522cSBarry Smith 11470971522cSBarry Smith PetscFunctionBegin; 11480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1149db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1150db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1151db4c96c1SJed Brown PetscValidIntPointer(fields,3); 11525d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 11530971522cSBarry Smith PetscFunctionReturn(0); 11540971522cSBarry Smith } 11550971522cSBarry Smith 11560971522cSBarry Smith #undef __FUNCT__ 1157704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1158704ba839SBarry Smith /*@ 1159704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1160704ba839SBarry Smith 1161ad4df100SBarry Smith Logically Collective on PC 1162704ba839SBarry Smith 1163704ba839SBarry Smith Input Parameters: 1164704ba839SBarry Smith + pc - the preconditioner context 1165a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1166db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1167704ba839SBarry Smith 1168d32f9abdSBarry Smith 1169a6ffb8dbSJed Brown Notes: 1170a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1171a6ffb8dbSJed Brown 1172db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1173db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1174d32f9abdSBarry Smith 1175704ba839SBarry Smith Level: intermediate 1176704ba839SBarry Smith 1177704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1178704ba839SBarry Smith 1179704ba839SBarry Smith @*/ 11807087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1181704ba839SBarry Smith { 11824ac538c5SBarry Smith PetscErrorCode ierr; 1183704ba839SBarry Smith 1184704ba839SBarry Smith PetscFunctionBegin; 11850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11867b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1187db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11884ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1189704ba839SBarry Smith PetscFunctionReturn(0); 1190704ba839SBarry Smith } 1191704ba839SBarry Smith 1192704ba839SBarry Smith #undef __FUNCT__ 119357a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 119457a9adfeSMatthew G Knepley /*@ 119557a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 119657a9adfeSMatthew G Knepley 119757a9adfeSMatthew G Knepley Logically Collective on PC 119857a9adfeSMatthew G Knepley 119957a9adfeSMatthew G Knepley Input Parameters: 120057a9adfeSMatthew G Knepley + pc - the preconditioner context 120157a9adfeSMatthew G Knepley - splitname - name of this split 120257a9adfeSMatthew G Knepley 120357a9adfeSMatthew G Knepley Output Parameter: 120457a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 120557a9adfeSMatthew G Knepley 120657a9adfeSMatthew G Knepley Level: intermediate 120757a9adfeSMatthew G Knepley 120857a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 120957a9adfeSMatthew G Knepley 121057a9adfeSMatthew G Knepley @*/ 121157a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 121257a9adfeSMatthew G Knepley { 121357a9adfeSMatthew G Knepley PetscErrorCode ierr; 121457a9adfeSMatthew G Knepley 121557a9adfeSMatthew G Knepley PetscFunctionBegin; 121657a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 121757a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 121857a9adfeSMatthew G Knepley PetscValidPointer(is,3); 121957a9adfeSMatthew G Knepley { 122057a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 122157a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 122257a9adfeSMatthew G Knepley PetscBool found; 122357a9adfeSMatthew G Knepley 122457a9adfeSMatthew G Knepley *is = PETSC_NULL; 122557a9adfeSMatthew G Knepley while(ilink) { 122657a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 122757a9adfeSMatthew G Knepley if (found) { 122857a9adfeSMatthew G Knepley *is = ilink->is; 122957a9adfeSMatthew G Knepley break; 123057a9adfeSMatthew G Knepley } 123157a9adfeSMatthew G Knepley ilink = ilink->next; 123257a9adfeSMatthew G Knepley } 123357a9adfeSMatthew G Knepley } 123457a9adfeSMatthew G Knepley PetscFunctionReturn(0); 123557a9adfeSMatthew G Knepley } 123657a9adfeSMatthew G Knepley 123757a9adfeSMatthew G Knepley #undef __FUNCT__ 123851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 123951f519a2SBarry Smith /*@ 124051f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 124151f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 124251f519a2SBarry Smith 1243ad4df100SBarry Smith Logically Collective on PC 124451f519a2SBarry Smith 124551f519a2SBarry Smith Input Parameters: 124651f519a2SBarry Smith + pc - the preconditioner context 124751f519a2SBarry Smith - bs - the block size 124851f519a2SBarry Smith 124951f519a2SBarry Smith Level: intermediate 125051f519a2SBarry Smith 125151f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 125251f519a2SBarry Smith 125351f519a2SBarry Smith @*/ 12547087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 125551f519a2SBarry Smith { 12564ac538c5SBarry Smith PetscErrorCode ierr; 125751f519a2SBarry Smith 125851f519a2SBarry Smith PetscFunctionBegin; 12590700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1260c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 12614ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 126251f519a2SBarry Smith PetscFunctionReturn(0); 126351f519a2SBarry Smith } 126451f519a2SBarry Smith 126551f519a2SBarry Smith #undef __FUNCT__ 126669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 12670971522cSBarry Smith /*@C 126869a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 12690971522cSBarry Smith 127069a612a9SBarry Smith Collective on KSP 12710971522cSBarry Smith 12720971522cSBarry Smith Input Parameter: 12730971522cSBarry Smith . pc - the preconditioner context 12740971522cSBarry Smith 12750971522cSBarry Smith Output Parameters: 127613e0d083SBarry Smith + n - the number of splits 127769a612a9SBarry Smith - pc - the array of KSP contexts 12780971522cSBarry Smith 12790971522cSBarry Smith Note: 1280d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1281d32f9abdSBarry Smith (not the KSP just the array that contains them). 12820971522cSBarry Smith 128369a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12840971522cSBarry Smith 12850971522cSBarry Smith Level: advanced 12860971522cSBarry Smith 12870971522cSBarry Smith .seealso: PCFIELDSPLIT 12880971522cSBarry Smith @*/ 12897087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12900971522cSBarry Smith { 12914ac538c5SBarry Smith PetscErrorCode ierr; 12920971522cSBarry Smith 12930971522cSBarry Smith PetscFunctionBegin; 12940700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 129513e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 12964ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 12970971522cSBarry Smith PetscFunctionReturn(0); 12980971522cSBarry Smith } 12990971522cSBarry Smith 1300e69d4d44SBarry Smith #undef __FUNCT__ 1301e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1302e69d4d44SBarry Smith /*@ 1303e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1304a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1305e69d4d44SBarry Smith 1306e69d4d44SBarry Smith Collective on PC 1307e69d4d44SBarry Smith 1308e69d4d44SBarry Smith Input Parameters: 1309e69d4d44SBarry Smith + pc - the preconditioner context 1310fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1311084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1312084e4875SJed Brown 1313e69d4d44SBarry Smith Options Database: 1314084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1315e69d4d44SBarry Smith 1316fd1303e9SJungho Lee Notes: 1317fd1303e9SJungho Lee $ If ptype is 1318fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1319fd1303e9SJungho Lee $ to this function). 1320fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1321fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1322fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1323fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1324fd1303e9SJungho Lee $ preconditioner 1325fd1303e9SJungho Lee 1326fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1327fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1328fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1329fd1303e9SJungho Lee 1330fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1331fd1303e9SJungho Lee the name since it sets a proceedure to use. 1332fd1303e9SJungho Lee 1333e69d4d44SBarry Smith Level: intermediate 1334e69d4d44SBarry Smith 1335fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1336e69d4d44SBarry Smith 1337e69d4d44SBarry Smith @*/ 13387087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1339e69d4d44SBarry Smith { 13404ac538c5SBarry Smith PetscErrorCode ierr; 1341e69d4d44SBarry Smith 1342e69d4d44SBarry Smith PetscFunctionBegin; 13430700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13444ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1345e69d4d44SBarry Smith PetscFunctionReturn(0); 1346e69d4d44SBarry Smith } 1347e69d4d44SBarry Smith 1348e69d4d44SBarry Smith EXTERN_C_BEGIN 1349e69d4d44SBarry Smith #undef __FUNCT__ 1350e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 13517087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1352e69d4d44SBarry Smith { 1353e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1354084e4875SJed Brown PetscErrorCode ierr; 1355e69d4d44SBarry Smith 1356e69d4d44SBarry Smith PetscFunctionBegin; 1357084e4875SJed Brown jac->schurpre = ptype; 1358084e4875SJed Brown if (pre) { 13596bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1360084e4875SJed Brown jac->schur_user = pre; 1361084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1362084e4875SJed Brown } 1363e69d4d44SBarry Smith PetscFunctionReturn(0); 1364e69d4d44SBarry Smith } 1365e69d4d44SBarry Smith EXTERN_C_END 1366e69d4d44SBarry Smith 136730ad9308SMatthew Knepley #undef __FUNCT__ 1368c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1369ab1df9f5SJed Brown /*@ 1370c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1371ab1df9f5SJed Brown 1372ab1df9f5SJed Brown Collective on PC 1373ab1df9f5SJed Brown 1374ab1df9f5SJed Brown Input Parameters: 1375ab1df9f5SJed Brown + pc - the preconditioner context 1376c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1377ab1df9f5SJed Brown 1378ab1df9f5SJed Brown Options Database: 1379c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1380ab1df9f5SJed Brown 1381ab1df9f5SJed Brown 1382ab1df9f5SJed Brown Level: intermediate 1383ab1df9f5SJed Brown 1384ab1df9f5SJed Brown Notes: 1385ab1df9f5SJed Brown The FULL factorization is 1386ab1df9f5SJed Brown 1387ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1388ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1389ab1df9f5SJed Brown 13906be4592eSBarry 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, 13916be4592eSBarry 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). 1392ab1df9f5SJed Brown 13936be4592eSBarry 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 13946be4592eSBarry 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 13956be4592eSBarry 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, 13966be4592eSBarry 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. 1397ab1df9f5SJed Brown 13986be4592eSBarry 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 13996be4592eSBarry 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). 1400ab1df9f5SJed Brown 1401ab1df9f5SJed Brown References: 1402ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1403ab1df9f5SJed Brown 1404ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1405ab1df9f5SJed Brown 1406ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1407ab1df9f5SJed Brown @*/ 1408c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1409ab1df9f5SJed Brown { 1410ab1df9f5SJed Brown PetscErrorCode ierr; 1411ab1df9f5SJed Brown 1412ab1df9f5SJed Brown PetscFunctionBegin; 1413ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1414c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1415ab1df9f5SJed Brown PetscFunctionReturn(0); 1416ab1df9f5SJed Brown } 1417ab1df9f5SJed Brown 1418ab1df9f5SJed Brown #undef __FUNCT__ 1419c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1420c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1421ab1df9f5SJed Brown { 1422ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1423ab1df9f5SJed Brown 1424ab1df9f5SJed Brown PetscFunctionBegin; 1425ab1df9f5SJed Brown jac->schurfactorization = ftype; 1426ab1df9f5SJed Brown PetscFunctionReturn(0); 1427ab1df9f5SJed Brown } 1428ab1df9f5SJed Brown 1429ab1df9f5SJed Brown #undef __FUNCT__ 143030ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 143130ad9308SMatthew Knepley /*@C 14328c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 143330ad9308SMatthew Knepley 143430ad9308SMatthew Knepley Collective on KSP 143530ad9308SMatthew Knepley 143630ad9308SMatthew Knepley Input Parameter: 143730ad9308SMatthew Knepley . pc - the preconditioner context 143830ad9308SMatthew Knepley 143930ad9308SMatthew Knepley Output Parameters: 1440a04f6461SBarry Smith + A00 - the (0,0) block 1441a04f6461SBarry Smith . A01 - the (0,1) block 1442a04f6461SBarry Smith . A10 - the (1,0) block 1443a04f6461SBarry Smith - A11 - the (1,1) block 144430ad9308SMatthew Knepley 144530ad9308SMatthew Knepley Level: advanced 144630ad9308SMatthew Knepley 144730ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 144830ad9308SMatthew Knepley @*/ 1449a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 145030ad9308SMatthew Knepley { 145130ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 145230ad9308SMatthew Knepley 145330ad9308SMatthew Knepley PetscFunctionBegin; 14540700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1455c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1456a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1457a04f6461SBarry Smith if (A01) *A01 = jac->B; 1458a04f6461SBarry Smith if (A10) *A10 = jac->C; 1459a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 146030ad9308SMatthew Knepley PetscFunctionReturn(0); 146130ad9308SMatthew Knepley } 146230ad9308SMatthew Knepley 146379416396SBarry Smith EXTERN_C_BEGIN 146479416396SBarry Smith #undef __FUNCT__ 146579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 14667087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 146779416396SBarry Smith { 146879416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1469e69d4d44SBarry Smith PetscErrorCode ierr; 147079416396SBarry Smith 147179416396SBarry Smith PetscFunctionBegin; 147279416396SBarry Smith jac->type = type; 14733b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 14743b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 14753b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1476e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1477e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1478c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1479e69d4d44SBarry Smith 14803b224e63SBarry Smith } else { 14813b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 14823b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1483e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14849e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1485c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 14863b224e63SBarry Smith } 148779416396SBarry Smith PetscFunctionReturn(0); 148879416396SBarry Smith } 148979416396SBarry Smith EXTERN_C_END 149079416396SBarry Smith 149151f519a2SBarry Smith EXTERN_C_BEGIN 149251f519a2SBarry Smith #undef __FUNCT__ 149351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 14947087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 149551f519a2SBarry Smith { 149651f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 149751f519a2SBarry Smith 149851f519a2SBarry Smith PetscFunctionBegin; 149965e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 150065e19b50SBarry 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); 150151f519a2SBarry Smith jac->bs = bs; 150251f519a2SBarry Smith PetscFunctionReturn(0); 150351f519a2SBarry Smith } 150451f519a2SBarry Smith EXTERN_C_END 150551f519a2SBarry Smith 150679416396SBarry Smith #undef __FUNCT__ 150779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1508bc08b0f1SBarry Smith /*@ 150979416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 151079416396SBarry Smith 151179416396SBarry Smith Collective on PC 151279416396SBarry Smith 151379416396SBarry Smith Input Parameter: 151479416396SBarry Smith . pc - the preconditioner context 151581540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 151679416396SBarry Smith 151779416396SBarry Smith Options Database Key: 1518a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 151979416396SBarry Smith 152079416396SBarry Smith Level: Developer 152179416396SBarry Smith 152279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 152379416396SBarry Smith 152479416396SBarry Smith .seealso: PCCompositeSetType() 152579416396SBarry Smith 152679416396SBarry Smith @*/ 15277087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 152879416396SBarry Smith { 15294ac538c5SBarry Smith PetscErrorCode ierr; 153079416396SBarry Smith 153179416396SBarry Smith PetscFunctionBegin; 15320700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15334ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 153479416396SBarry Smith PetscFunctionReturn(0); 153579416396SBarry Smith } 153679416396SBarry Smith 15370971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 15380971522cSBarry Smith /*MC 1539a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1540a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 15410971522cSBarry Smith 1542edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1543edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 15440971522cSBarry Smith 1545a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 154669a612a9SBarry Smith and set the options directly on the resulting KSP object 15470971522cSBarry Smith 15480971522cSBarry Smith Level: intermediate 15490971522cSBarry Smith 155079416396SBarry Smith Options Database Keys: 155181540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 155281540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 155381540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 155481540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 15550f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 15560f188ba9SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag 1557435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 155879416396SBarry Smith 15595d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 15605d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 15615d4c12cdSJungho Lee 1562c8a0d604SMatthew G Knepley Notes: 1563c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1564d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1565d32f9abdSBarry Smith 1566d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1567d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1568d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1569d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1570d32f9abdSBarry Smith 1571c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1572c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1573c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1574c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1575c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1576a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1577c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1578c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1579c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1580a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1581c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1582c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 15835668aaf4SBarry Smith diag gives 1584c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1585c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 15865668aaf4SBarry 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 1587c8a0d604SMatthew G Knepley $ ( A00 0 ) 1588c8a0d604SMatthew G Knepley $ ( A10 S ) 1589c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1590c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1591c8a0d604SMatthew G Knepley $ ( 0 S ) 1592c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1593e69d4d44SBarry Smith 1594edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1595edf189efSBarry Smith is used automatically for a second block. 1596edf189efSBarry Smith 1597ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1598ff218e97SBarry Smith Generally it should be used with the AIJ format. 1599ff218e97SBarry Smith 1600ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1601ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1602ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 16030716a85fSBarry Smith 1604a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 16050971522cSBarry Smith 16067e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1607e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 16080971522cSBarry Smith M*/ 16090971522cSBarry Smith 16100971522cSBarry Smith EXTERN_C_BEGIN 16110971522cSBarry Smith #undef __FUNCT__ 16120971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 16137087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 16140971522cSBarry Smith { 16150971522cSBarry Smith PetscErrorCode ierr; 16160971522cSBarry Smith PC_FieldSplit *jac; 16170971522cSBarry Smith 16180971522cSBarry Smith PetscFunctionBegin; 161938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 16200971522cSBarry Smith jac->bs = -1; 16210971522cSBarry Smith jac->nsplits = 0; 16223e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1623e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1624c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 162551f519a2SBarry Smith 16260971522cSBarry Smith pc->data = (void*)jac; 16270971522cSBarry Smith 16280971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1629421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 16300971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1631574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 16320971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 16330971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 16340971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 16350971522cSBarry Smith pc->ops->applyrichardson = 0; 16360971522cSBarry Smith 163769a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 163869a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 16390971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 16400971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1641704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1642704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 164379416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 164479416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 164551f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 164651f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 16470971522cSBarry Smith PetscFunctionReturn(0); 16480971522cSBarry Smith } 16490971522cSBarry Smith EXTERN_C_END 16500971522cSBarry Smith 16510971522cSBarry Smith 1652a541d17aSBarry Smith 1653