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 } 83a3df900dSMatthew G Knepley if (jac->realdiagonal) { 84a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 85a3df900dSMatthew G Knepley } 8669a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 870971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 880971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 891ab39975SBarry Smith if (ilink->fields) { 900971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 9179416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 925a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9379416396SBarry Smith if (j > 0) { 9479416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9579416396SBarry Smith } 965a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 970971522cSBarry Smith } 980971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1001ab39975SBarry Smith } else { 1011ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1021ab39975SBarry Smith } 1035a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1045a9f2f41SSatish Balay ilink = ilink->next; 1050971522cSBarry Smith } 1060971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1070971522cSBarry Smith } else { 10865e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1090971522cSBarry Smith } 1100971522cSBarry Smith PetscFunctionReturn(0); 1110971522cSBarry Smith } 1120971522cSBarry Smith 1130971522cSBarry Smith #undef __FUNCT__ 1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1163b224e63SBarry Smith { 1173b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1183b224e63SBarry Smith PetscErrorCode ierr; 119ace3abfcSBarry Smith PetscBool iascii; 1203b224e63SBarry Smith PetscInt i,j; 1213b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1223b224e63SBarry Smith KSP ksp; 1233b224e63SBarry Smith 1243b224e63SBarry Smith PetscFunctionBegin; 125251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1263b224e63SBarry Smith if (iascii) { 1272c9966d7SBarry Smith if (jac->bs > 0) { 128c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1292c9966d7SBarry Smith } else { 1302c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1312c9966d7SBarry Smith } 132a3df900dSMatthew G Knepley if (jac->realdiagonal) { 133a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 134a3df900dSMatthew G Knepley } 1353e8b8b31SMatthew G Knepley switch(jac->schurpre) { 1363e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1373e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 138b02e2d75SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_DIAG: 1393e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break; 1403e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1413e8b8b31SMatthew G Knepley if (jac->schur_user) { 1423e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1433e8b8b31SMatthew G Knepley } else { 1443e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr); 1453e8b8b31SMatthew G Knepley } 1463e8b8b31SMatthew G Knepley break; 1473e8b8b31SMatthew G Knepley default: 1483e8b8b31SMatthew G Knepley SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1493e8b8b31SMatthew G Knepley } 1503b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1513b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1523b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1533b224e63SBarry Smith if (ilink->fields) { 1543b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1553b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1563b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1573b224e63SBarry Smith if (j > 0) { 1583b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1593b224e63SBarry Smith } 1603b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1613b224e63SBarry Smith } 1623b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1633b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1643b224e63SBarry Smith } else { 1653b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1663b224e63SBarry Smith } 1673b224e63SBarry Smith ilink = ilink->next; 1683b224e63SBarry Smith } 169435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1703b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 17112cae6f2SJed Brown if (jac->schur) { 17212cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1733b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 17412cae6f2SJed Brown } else { 17512cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 17612cae6f2SJed Brown } 1773b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 178435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1793b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 18012cae6f2SJed Brown if (jac->kspschur) { 1813b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 18212cae6f2SJed Brown } else { 18312cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 18412cae6f2SJed Brown } 1853b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1863b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1873b224e63SBarry Smith } else { 18865e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1893b224e63SBarry Smith } 1903b224e63SBarry Smith PetscFunctionReturn(0); 1913b224e63SBarry Smith } 1923b224e63SBarry Smith 1933b224e63SBarry Smith #undef __FUNCT__ 1946c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1956c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1966c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1976c924f48SJed Brown { 1986c924f48SJed Brown PetscErrorCode ierr; 1996c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2005d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2015d4c12cdSJungho Lee PetscBool flg,flg_col; 2025d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2036c924f48SJed Brown 2046c924f48SJed Brown PetscFunctionBegin; 2056c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2065d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2076c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2086c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2096c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2105d4c12cdSJungho Lee ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2116c924f48SJed Brown nfields = jac->bs; 21229499fbbSJungho Lee nfields_col = jac->bs; 2136c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2145d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2156c924f48SJed Brown if (!flg) break; 2165d4c12cdSJungho Lee else if (flg && !flg_col) { 2175d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2185d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2195d4c12cdSJungho Lee } 2205d4c12cdSJungho Lee else { 2215d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2225d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2235d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2245d4c12cdSJungho Lee } 2256c924f48SJed Brown } 2266c924f48SJed Brown if (i > 0) { 2276c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2286c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2296c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2306c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2316c924f48SJed Brown } 2326c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2335d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2346c924f48SJed Brown PetscFunctionReturn(0); 2356c924f48SJed Brown } 2366c924f48SJed Brown 2376c924f48SJed Brown #undef __FUNCT__ 23869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 23969a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2400971522cSBarry Smith { 2410971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2420971522cSBarry Smith PetscErrorCode ierr; 2435a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2447287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 2456c924f48SJed Brown PetscInt i; 2460971522cSBarry Smith 2470971522cSBarry Smith PetscFunctionBegin; 2487287d2a3SDmitry Karpeev /* 2497287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 2507287d2a3SDmitry Karpeev Should probably be rewritten. 2517287d2a3SDmitry Karpeev */ 2527287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 253d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 254d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2557b62db95SJungho Lee PetscInt numFields, f; 2560784a22cSJed Brown char **fieldNames; 2577b62db95SJungho Lee IS *fields; 258e7c4fc90SDmitry Karpeev DM *dms; 25916621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms); CHKERRQ(ierr); 2607b62db95SJungho Lee for(f = 0; f < numFields; ++f) { 2617287d2a3SDmitry Karpeev if(jac->reset) { 2627287d2a3SDmitry Karpeev ilink->is = fields[i]; 2637287d2a3SDmitry Karpeev ilink = ilink->next; 2647287d2a3SDmitry Karpeev } 2657287d2a3SDmitry Karpeev else { 2667b62db95SJungho Lee ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 2677287d2a3SDmitry Karpeev } 2687b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 2697b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 2707b62db95SJungho Lee } 2717b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 2727b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 273e7c4fc90SDmitry Karpeev if(dms) { 2748b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2757b62db95SJungho Lee for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) { 2767287d2a3SDmitry Karpeev const char* prefix; 2777287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr); 2787287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix); CHKERRQ(ierr); 2797b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]); CHKERRQ(ierr); 2807b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE); CHKERRQ(ierr); 2817287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr); 282e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 2832fa5ba8aSJed Brown } 2847b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 2858b8307b2SJed Brown } 28666ffff09SJed Brown } else { 287521d7252SBarry Smith if (jac->bs <= 0) { 288704ba839SBarry Smith if (pc->pmat) { 289521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 290704ba839SBarry Smith } else { 291704ba839SBarry Smith jac->bs = 1; 292704ba839SBarry Smith } 293521d7252SBarry Smith } 294d32f9abdSBarry Smith 2957287d2a3SDmitry Karpeev ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr); 2966ce1633cSBarry Smith if (stokes) { 2976ce1633cSBarry Smith IS zerodiags,rest; 2986ce1633cSBarry Smith PetscInt nmin,nmax; 2996ce1633cSBarry Smith 3006ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3016ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3026ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 3037287d2a3SDmitry Karpeev if(jac->reset) { 3047287d2a3SDmitry Karpeev jac->head->is = rest; 3057287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 3067287d2a3SDmitry Karpeev } 3077287d2a3SDmitry Karpeev else { 3086ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 3096ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 3107287d2a3SDmitry Karpeev } 311fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 312fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 3136ce1633cSBarry Smith } else { 3147287d2a3SDmitry Karpeev if(jac->reset) 3157287d2a3SDmitry Karpeev SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 3167287d2a3SDmitry Karpeev if (!fieldsplit_default) { 317d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 318d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 3196c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 3206c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 321d32f9abdSBarry Smith } 3227287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 323d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 324db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 3256c924f48SJed Brown char splitname[8]; 3266c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 3275d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 32879416396SBarry Smith } 3295d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 330521d7252SBarry Smith } 33166ffff09SJed Brown } 3326ce1633cSBarry Smith } 333edf189efSBarry Smith } else if (jac->nsplits == 1) { 334edf189efSBarry Smith if (ilink->is) { 335edf189efSBarry Smith IS is2; 336edf189efSBarry Smith PetscInt nmin,nmax; 337edf189efSBarry Smith 338edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 339edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 340db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 341fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 3427b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 343edf189efSBarry Smith } 344d0af7cd3SBarry Smith 345d0af7cd3SBarry Smith 3467b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 34769a612a9SBarry Smith PetscFunctionReturn(0); 34869a612a9SBarry Smith } 34969a612a9SBarry Smith 35069a612a9SBarry Smith #undef __FUNCT__ 35169a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 35269a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 35369a612a9SBarry Smith { 35469a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 35569a612a9SBarry Smith PetscErrorCode ierr; 3565a9f2f41SSatish Balay PC_FieldSplitLink ilink; 3572c9966d7SBarry Smith PetscInt i,nsplit; 35869a612a9SBarry Smith MatStructure flag = pc->flag; 3595f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 36069a612a9SBarry Smith 36169a612a9SBarry Smith PetscFunctionBegin; 36269a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 36397bbdb24SBarry Smith nsplit = jac->nsplits; 3645a9f2f41SSatish Balay ilink = jac->head; 36597bbdb24SBarry Smith 36697bbdb24SBarry Smith /* get the matrices for each split */ 367704ba839SBarry Smith if (!jac->issetup) { 3681b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 36997bbdb24SBarry Smith 370704ba839SBarry Smith jac->issetup = PETSC_TRUE; 371704ba839SBarry Smith 3725d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 3732c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 3742c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 3752c9966d7SBarry Smith } 37651f519a2SBarry Smith bs = jac->bs; 37797bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 3781b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3791b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3801b9fc7fcSBarry Smith if (jac->defaultsplit) { 381704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 3825f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 383704ba839SBarry Smith } else if (!ilink->is) { 384ccb205f8SBarry Smith if (ilink->nfields > 1) { 3855f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 3865a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3875f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 3881b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3891b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3901b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 3915f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 39297bbdb24SBarry Smith } 39397bbdb24SBarry Smith } 394d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 3955f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 396ccb205f8SBarry Smith } else { 397704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 3985f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 399ccb205f8SBarry Smith } 4003e197d65SBarry Smith } 401704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4025f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4035f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 404704ba839SBarry Smith ilink = ilink->next; 4051b9fc7fcSBarry Smith } 4061b9fc7fcSBarry Smith } 4071b9fc7fcSBarry Smith 408704ba839SBarry Smith ilink = jac->head; 40997bbdb24SBarry Smith if (!jac->pmat) { 410bdddcaaaSMatthew G Knepley Vec xtmp; 411bdddcaaaSMatthew G Knepley 412bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 413cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 414bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 415cf502942SBarry Smith for (i=0; i<nsplit; i++) { 416bdddcaaaSMatthew G Knepley MatNullSpace sp; 417bdddcaaaSMatthew G Knepley 418a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 419a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 420a3df900dSMatthew G Knepley if (jac->pmat[i]) { 421a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 422a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 423a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 424a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 425a3df900dSMatthew G Knepley } 426a3df900dSMatthew G Knepley } else { 4275f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 428a3df900dSMatthew G Knepley } 429bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 430bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 431bdddcaaaSMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; 432bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 433bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 434bdddcaaaSMatthew G Knepley /* HACK: Check for the constant null space */ 435bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(pc->pmat, &sp);CHKERRQ(ierr); 436bdddcaaaSMatthew G Knepley if (sp) { 437bdddcaaaSMatthew G Knepley MatNullSpace subsp; 438bdddcaaaSMatthew G Knepley Vec ftmp, gtmp; 4399583d628SBarry Smith PetscReal norm; 440bdddcaaaSMatthew G Knepley PetscInt N; 441bdddcaaaSMatthew G Knepley 442bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat, >mp, PETSC_NULL);CHKERRQ(ierr); 443bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i], &ftmp, PETSC_NULL);CHKERRQ(ierr); 444bdddcaaaSMatthew G Knepley ierr = VecGetSize(ftmp, &N);CHKERRQ(ierr); 445bdddcaaaSMatthew G Knepley ierr = VecSet(ftmp, 1.0/N);CHKERRQ(ierr); 446bdddcaaaSMatthew G Knepley ierr = VecScatterBegin(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 447bdddcaaaSMatthew G Knepley ierr = VecScatterEnd(ilink->sctx, ftmp, gtmp, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 448bdddcaaaSMatthew G Knepley ierr = MatNullSpaceRemove(sp, gtmp, PETSC_NULL);CHKERRQ(ierr); 449bdddcaaaSMatthew G Knepley ierr = VecNorm(gtmp, NORM_2, &norm);CHKERRQ(ierr); 450bdddcaaaSMatthew G Knepley if (norm < 1.0e-10) { 451bdddcaaaSMatthew G Knepley ierr = MatNullSpaceCreate(((PetscObject)pc)->comm, PETSC_TRUE, 0, PETSC_NULL, &subsp);CHKERRQ(ierr); 452bdddcaaaSMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], subsp);CHKERRQ(ierr); 453bdddcaaaSMatthew G Knepley ierr = MatNullSpaceDestroy(&subsp);CHKERRQ(ierr); 454bdddcaaaSMatthew G Knepley } 455bdddcaaaSMatthew G Knepley ierr = VecDestroy(&ftmp);CHKERRQ(ierr); 456bdddcaaaSMatthew G Knepley ierr = VecDestroy(>mp);CHKERRQ(ierr); 457bdddcaaaSMatthew G Knepley } 458ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 459ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 460ed1f0337SMatthew G Knepley if (sp) { 461ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 462ed1f0337SMatthew G Knepley } 463704ba839SBarry Smith ilink = ilink->next; 464cf502942SBarry Smith } 465bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 46697bbdb24SBarry Smith } else { 467cf502942SBarry Smith for (i=0; i<nsplit; i++) { 468a3df900dSMatthew G Knepley Mat pmat; 469a3df900dSMatthew G Knepley 470a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 471a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 472a3df900dSMatthew G Knepley if (!pmat) { 4735f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 474a3df900dSMatthew G Knepley } 475704ba839SBarry Smith ilink = ilink->next; 476cf502942SBarry Smith } 47797bbdb24SBarry Smith } 478519d70e2SJed Brown if (jac->realdiagonal) { 479519d70e2SJed Brown ilink = jac->head; 480519d70e2SJed Brown if (!jac->mat) { 481519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 482519d70e2SJed Brown for (i=0; i<nsplit; i++) { 4835f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 484519d70e2SJed Brown ilink = ilink->next; 485519d70e2SJed Brown } 486519d70e2SJed Brown } else { 487519d70e2SJed Brown for (i=0; i<nsplit; i++) { 4885f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 489519d70e2SJed Brown ilink = ilink->next; 490519d70e2SJed Brown } 491519d70e2SJed Brown } 492519d70e2SJed Brown } else { 493519d70e2SJed Brown jac->mat = jac->pmat; 494519d70e2SJed Brown } 49597bbdb24SBarry Smith 4966c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 49768dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 49868dd23aaSBarry Smith ilink = jac->head; 49968dd23aaSBarry Smith if (!jac->Afield) { 50068dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 50168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5024aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 50368dd23aaSBarry Smith ilink = ilink->next; 50468dd23aaSBarry Smith } 50568dd23aaSBarry Smith } else { 50668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5074aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 50868dd23aaSBarry Smith ilink = ilink->next; 50968dd23aaSBarry Smith } 51068dd23aaSBarry Smith } 51168dd23aaSBarry Smith } 51268dd23aaSBarry Smith 5133b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5143b224e63SBarry Smith IS ccis; 5154aa3045dSJed Brown PetscInt rstart,rend; 516093c86ffSJed Brown char lscname[256]; 517093c86ffSJed Brown PetscObject LSC_L; 518e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 51968dd23aaSBarry Smith 520e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 521e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 522e6cab6aaSJed Brown 5233b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5243b224e63SBarry Smith if (jac->schur) { 5253b224e63SBarry Smith ilink = jac->head; 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->B);CHKERRQ(ierr); 528fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5293b224e63SBarry Smith ilink = ilink->next; 53049bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5314aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 532fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 533a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 534084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 5353b224e63SBarry Smith } else { 5361cee3971SBarry Smith KSP ksp; 5377287d2a3SDmitry Karpeev const char *Aprefix, *Dprefix; 538bdddcaaaSMatthew G Knepley MatNullSpace sp; 5397287d2a3SDmitry Karpeev DM Adm; 5403b224e63SBarry Smith 541a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 5423b224e63SBarry Smith ilink = jac->head; 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->B);CHKERRQ(ierr); 545fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5463b224e63SBarry Smith ilink = ilink->next; 54749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5484aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 549fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5507d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 551519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 552bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 553bdddcaaaSMatthew G Knepley if (sp) {ierr = MatSetNullSpace(jac->schur, sp); CHKERRQ(ierr);} 5547287d2a3SDmitry Karpeev /* set tabbing, options prefix and DM of KSP inside the MatSchur (inherited from the split) */ 5551cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp); CHKERRQ(ierr); 5567287d2a3SDmitry Karpeev ierr = KSPGetDM(jac->head->ksp,&Adm); CHKERRQ(ierr); 5577287d2a3SDmitry Karpeev ierr = KSPSetDM(ksp,Adm); CHKERRQ(ierr); 5587287d2a3SDmitry Karpeev ierr = KSPSetDMActive(ksp,PETSC_FALSE); CHKERRQ(ierr); 5597287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->ksp,&Aprefix); CHKERRQ(ierr); 5607287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(ksp,Aprefix); CHKERRQ(ierr); 5617287d2a3SDmitry Karpeev ierr = KSPIncrementTabLevel(ksp,(PetscObject)pc,2); CHKERRQ(ierr); 5627287d2a3SDmitry Karpeev ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 56320b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 56420b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 5657b62db95SJungho Lee ierr = MatSetUp(jac->schur); CHKERRQ(ierr); 5663b224e63SBarry Smith 5677287d2a3SDmitry Karpeev 5683b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur); CHKERRQ(ierr); 5699005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 5707287d2a3SDmitry Karpeev ierr = KSPIncrementTabLevel(jac->kspschur,(PetscObject)pc,1); CHKERRQ(ierr); 571084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 572084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 573084e4875SJed Brown PC pc; 574cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 575084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 576084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 577e69d4d44SBarry Smith } 5787287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr); 5797287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix); CHKERRQ(ierr); 5803b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 58120b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 58220b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 5833b224e63SBarry Smith } 584093c86ffSJed Brown 585093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 586093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 587093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 588093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 589093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 590093c86ffSJed Brown ierr = PetscSNPrintf(lscname,sizeof lscname,"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 591093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 592093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 593093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 5943b224e63SBarry Smith } else { 59597bbdb24SBarry Smith /* set up the individual PCs */ 59697bbdb24SBarry Smith i = 0; 5975a9f2f41SSatish Balay ilink = jac->head; 5985a9f2f41SSatish Balay while (ilink) { 599519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 6003b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 601c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 6025a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 60397bbdb24SBarry Smith i++; 6045a9f2f41SSatish Balay ilink = ilink->next; 6050971522cSBarry Smith } 6063b224e63SBarry Smith } 6073b224e63SBarry Smith 608c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 6090971522cSBarry Smith PetscFunctionReturn(0); 6100971522cSBarry Smith } 6110971522cSBarry Smith 6125a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 613ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 614ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 6155a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 616ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 617ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 61879416396SBarry Smith 6190971522cSBarry Smith #undef __FUNCT__ 6203b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 6213b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 6223b224e63SBarry Smith { 6233b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6243b224e63SBarry Smith PetscErrorCode ierr; 6253b224e63SBarry Smith KSP ksp; 6263b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 6273b224e63SBarry Smith 6283b224e63SBarry Smith PetscFunctionBegin; 6293b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 6303b224e63SBarry Smith 631c5d2311dSJed Brown switch (jac->schurfactorization) { 632c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 633a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 634c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 635c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 636c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 637c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 638c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 639c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 640c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 641c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 642c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 643c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 644c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 645c5d2311dSJed Brown break; 646c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 647a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 648c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 649c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 650c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 651c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 652c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 653c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 654c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 655c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 656c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 657c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 658c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 659c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 660c5d2311dSJed Brown break; 661c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 662a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 663c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 664c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 665c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 666c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 667c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 668c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 669c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 670c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 671c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 672c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 673c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 674c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 675c5d2311dSJed Brown break; 676c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 677a04f6461SBarry 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 */ 6783b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6793b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6803b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6813b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6823b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6833b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6843b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6853b224e63SBarry Smith 6863b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6873b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6883b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6893b224e63SBarry Smith 6903b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6913b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6923b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6933b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6943b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 695c5d2311dSJed Brown } 6963b224e63SBarry Smith PetscFunctionReturn(0); 6973b224e63SBarry Smith } 6983b224e63SBarry Smith 6993b224e63SBarry Smith #undef __FUNCT__ 7000971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 7010971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 7020971522cSBarry Smith { 7030971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7040971522cSBarry Smith PetscErrorCode ierr; 7055a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 706939b8a20SBarry Smith PetscInt cnt,bs; 7070971522cSBarry Smith 7080971522cSBarry Smith PetscFunctionBegin; 70951f519a2SBarry Smith CHKMEMQ; 71079416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 7111b9fc7fcSBarry Smith if (jac->defaultsplit) { 712939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 713939b8a20SBarry 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); 714939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 715939b8a20SBarry 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); 7160971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 7175a9f2f41SSatish Balay while (ilink) { 7185a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7195a9f2f41SSatish Balay ilink = ilink->next; 7200971522cSBarry Smith } 7210971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 7221b9fc7fcSBarry Smith } else { 723efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7245a9f2f41SSatish Balay while (ilink) { 7255a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 7265a9f2f41SSatish Balay ilink = ilink->next; 7271b9fc7fcSBarry Smith } 7281b9fc7fcSBarry Smith } 72916913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 73079416396SBarry Smith if (!jac->w1) { 73179416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 73279416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 73379416396SBarry Smith } 734efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 7355a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 7363e197d65SBarry Smith cnt = 1; 7375a9f2f41SSatish Balay while (ilink->next) { 7385a9f2f41SSatish Balay ilink = ilink->next; 7393e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 7403e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 7413e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 7423e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7433e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7443e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7453e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7463e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7473e197d65SBarry Smith } 74851f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 74911755939SBarry Smith cnt -= 2; 75051f519a2SBarry Smith while (ilink->previous) { 75151f519a2SBarry Smith ilink = ilink->previous; 75211755939SBarry Smith /* compute the residual only over the part of the vector needed */ 75311755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 75411755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 75511755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 75611755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 75711755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 75811755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 75911755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 76051f519a2SBarry Smith } 76111755939SBarry Smith } 76265e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 76351f519a2SBarry Smith CHKMEMQ; 7640971522cSBarry Smith PetscFunctionReturn(0); 7650971522cSBarry Smith } 7660971522cSBarry Smith 767421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 768ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 769ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 770421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 771ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 772ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 773421e10b8SBarry Smith 774421e10b8SBarry Smith #undef __FUNCT__ 7758c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 776421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 777421e10b8SBarry Smith { 778421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 779421e10b8SBarry Smith PetscErrorCode ierr; 780421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 781939b8a20SBarry Smith PetscInt bs; 782421e10b8SBarry Smith 783421e10b8SBarry Smith PetscFunctionBegin; 784421e10b8SBarry Smith CHKMEMQ; 785421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 786421e10b8SBarry Smith if (jac->defaultsplit) { 787939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 788939b8a20SBarry 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); 789939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 790939b8a20SBarry 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); 791421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 792421e10b8SBarry Smith while (ilink) { 793421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 794421e10b8SBarry Smith ilink = ilink->next; 795421e10b8SBarry Smith } 796421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 797421e10b8SBarry Smith } else { 798421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 799421e10b8SBarry Smith while (ilink) { 800421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 801421e10b8SBarry Smith ilink = ilink->next; 802421e10b8SBarry Smith } 803421e10b8SBarry Smith } 804421e10b8SBarry Smith } else { 805421e10b8SBarry Smith if (!jac->w1) { 806421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 807421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 808421e10b8SBarry Smith } 809421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 810421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 811421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 812421e10b8SBarry Smith while (ilink->next) { 813421e10b8SBarry Smith ilink = ilink->next; 8149989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 815421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 816421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 817421e10b8SBarry Smith } 818421e10b8SBarry Smith while (ilink->previous) { 819421e10b8SBarry Smith ilink = ilink->previous; 8209989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 821421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 822421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 823421e10b8SBarry Smith } 824421e10b8SBarry Smith } else { 825421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 826421e10b8SBarry Smith ilink = ilink->next; 827421e10b8SBarry Smith } 828421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 829421e10b8SBarry Smith while (ilink->previous) { 830421e10b8SBarry Smith ilink = ilink->previous; 8319989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 832421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 833421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 834421e10b8SBarry Smith } 835421e10b8SBarry Smith } 836421e10b8SBarry Smith } 837421e10b8SBarry Smith CHKMEMQ; 838421e10b8SBarry Smith PetscFunctionReturn(0); 839421e10b8SBarry Smith } 840421e10b8SBarry Smith 8410971522cSBarry Smith #undef __FUNCT__ 842574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 843574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 8440971522cSBarry Smith { 8450971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8460971522cSBarry Smith PetscErrorCode ierr; 8475a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 8480971522cSBarry Smith 8490971522cSBarry Smith PetscFunctionBegin; 8505a9f2f41SSatish Balay while (ilink) { 851574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 852fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 853fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 854fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 855fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 856b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 8575a9f2f41SSatish Balay next = ilink->next; 8585a9f2f41SSatish Balay ilink = next; 8590971522cSBarry Smith } 86005b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 861574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 862574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 863574deadeSBarry Smith } else if (jac->mat) { 864574deadeSBarry Smith jac->mat = PETSC_NULL; 865574deadeSBarry Smith } 86697bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 86768dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8686bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8696bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8706bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8716bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8726bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8736bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8746bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 87563ec74ffSBarry Smith jac->reset = PETSC_TRUE; 876574deadeSBarry Smith PetscFunctionReturn(0); 877574deadeSBarry Smith } 878574deadeSBarry Smith 879574deadeSBarry Smith #undef __FUNCT__ 880574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 881574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 882574deadeSBarry Smith { 883574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 884574deadeSBarry Smith PetscErrorCode ierr; 885574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 886574deadeSBarry Smith 887574deadeSBarry Smith PetscFunctionBegin; 888574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 889574deadeSBarry Smith while (ilink) { 8906bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 891574deadeSBarry Smith next = ilink->next; 892574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 893574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 8945d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 895574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 896574deadeSBarry Smith ilink = next; 897574deadeSBarry Smith } 898574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 899c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 9007b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 9017b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 9027b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 9037b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 9047b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 905ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 906c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 9070971522cSBarry Smith PetscFunctionReturn(0); 9080971522cSBarry Smith } 9090971522cSBarry Smith 9100971522cSBarry Smith #undef __FUNCT__ 9110971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 9120971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 9130971522cSBarry Smith { 9141b9fc7fcSBarry Smith PetscErrorCode ierr; 9156c924f48SJed Brown PetscInt bs; 916bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 9179dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9183b224e63SBarry Smith PCCompositeType ctype; 919e7c4fc90SDmitry Karpeev DM ddm; 920e7c4fc90SDmitry Karpeev char ddm_name[1025]; 9211b9fc7fcSBarry Smith 9220971522cSBarry Smith PetscFunctionBegin; 9231b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 924e7c4fc90SDmitry Karpeev if(pc->dm) { 925e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 926e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 927*731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 928e7c4fc90SDmitry Karpeev if(flg) { 92916621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 930e7c4fc90SDmitry Karpeev if(!ddm) { 93116621825SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 932e7c4fc90SDmitry Karpeev } 93316621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 934e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 935e7c4fc90SDmitry Karpeev } 936e7c4fc90SDmitry Karpeev } 937acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 93851f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 93951f519a2SBarry Smith if (flg) { 94051f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 94151f519a2SBarry Smith } 942704ba839SBarry Smith 943435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 944c0adfefeSBarry Smith if (stokes) { 945c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 946c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 947c0adfefeSBarry Smith } 948c0adfefeSBarry Smith 9493b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 9503b224e63SBarry Smith if (flg) { 9513b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 9523b224e63SBarry Smith } 953c30613efSMatthew Knepley /* Only setup fields once */ 954c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 955d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 956d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 9576c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 9586c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 959d32f9abdSBarry Smith } 960c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 961c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 962c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 963c9c6ffaaSJed 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); 964084e4875SJed 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); 965c5d2311dSJed Brown } 9661b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 9670971522cSBarry Smith PetscFunctionReturn(0); 9680971522cSBarry Smith } 9690971522cSBarry Smith 9700971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 9710971522cSBarry Smith 9720971522cSBarry Smith EXTERN_C_BEGIN 9730971522cSBarry Smith #undef __FUNCT__ 9740971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 9755d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 9760971522cSBarry Smith { 97797bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9780971522cSBarry Smith PetscErrorCode ierr; 9795a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 98069a612a9SBarry Smith char prefix[128]; 9815d4c12cdSJungho Lee PetscInt i; 9820971522cSBarry Smith 9830971522cSBarry Smith PetscFunctionBegin; 9846c924f48SJed Brown if (jac->splitdefined) { 9856c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9866c924f48SJed Brown PetscFunctionReturn(0); 9876c924f48SJed Brown } 98851f519a2SBarry Smith for (i=0; i<n; i++) { 989e32f2f54SBarry 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); 990e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 99151f519a2SBarry Smith } 992704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 993a04f6461SBarry Smith if (splitname) { 994db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 995a04f6461SBarry Smith } else { 996a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 997a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 998a04f6461SBarry Smith } 999704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 10005a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 10015d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 10025d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 10035a9f2f41SSatish Balay ilink->nfields = n; 10045a9f2f41SSatish Balay ilink->next = PETSC_NULL; 10057adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10067287d2a3SDmitry Karpeev ierr = KSPIncrementTabLevel(ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 10075a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10089005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 100969a612a9SBarry Smith 1010a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 10115a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 10120971522cSBarry Smith 10130971522cSBarry Smith if (!next) { 10145a9f2f41SSatish Balay jac->head = ilink; 101551f519a2SBarry Smith ilink->previous = PETSC_NULL; 10160971522cSBarry Smith } else { 10170971522cSBarry Smith while (next->next) { 10180971522cSBarry Smith next = next->next; 10190971522cSBarry Smith } 10205a9f2f41SSatish Balay next->next = ilink; 102151f519a2SBarry Smith ilink->previous = next; 10220971522cSBarry Smith } 10230971522cSBarry Smith jac->nsplits++; 10240971522cSBarry Smith PetscFunctionReturn(0); 10250971522cSBarry Smith } 10260971522cSBarry Smith EXTERN_C_END 10270971522cSBarry Smith 1028e69d4d44SBarry Smith EXTERN_C_BEGIN 1029e69d4d44SBarry Smith #undef __FUNCT__ 1030e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 10317087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1032e69d4d44SBarry Smith { 1033e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1034e69d4d44SBarry Smith PetscErrorCode ierr; 1035e69d4d44SBarry Smith 1036e69d4d44SBarry Smith PetscFunctionBegin; 1037519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1038e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1039e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 104013e0d083SBarry Smith if (n) *n = jac->nsplits; 1041e69d4d44SBarry Smith PetscFunctionReturn(0); 1042e69d4d44SBarry Smith } 1043e69d4d44SBarry Smith EXTERN_C_END 10440971522cSBarry Smith 10450971522cSBarry Smith EXTERN_C_BEGIN 10460971522cSBarry Smith #undef __FUNCT__ 104769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 10487087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 10490971522cSBarry Smith { 10500971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10510971522cSBarry Smith PetscErrorCode ierr; 10520971522cSBarry Smith PetscInt cnt = 0; 10535a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 10540971522cSBarry Smith 10550971522cSBarry Smith PetscFunctionBegin; 10565d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 10575a9f2f41SSatish Balay while (ilink) { 10585a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 10595a9f2f41SSatish Balay ilink = ilink->next; 10600971522cSBarry Smith } 10615d480477SMatthew 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); 106213e0d083SBarry Smith if (n) *n = jac->nsplits; 10630971522cSBarry Smith PetscFunctionReturn(0); 10640971522cSBarry Smith } 10650971522cSBarry Smith EXTERN_C_END 10660971522cSBarry Smith 1067704ba839SBarry Smith EXTERN_C_BEGIN 1068704ba839SBarry Smith #undef __FUNCT__ 1069704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 10707087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1071704ba839SBarry Smith { 1072704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1073704ba839SBarry Smith PetscErrorCode ierr; 1074704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1075704ba839SBarry Smith char prefix[128]; 1076704ba839SBarry Smith 1077704ba839SBarry Smith PetscFunctionBegin; 10786c924f48SJed Brown if (jac->splitdefined) { 10796c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10806c924f48SJed Brown PetscFunctionReturn(0); 10816c924f48SJed Brown } 108216913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1083a04f6461SBarry Smith if (splitname) { 1084db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1085a04f6461SBarry Smith } else { 1086b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1087b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1088a04f6461SBarry Smith } 10891ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1090b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1091b5787286SJed Brown ilink->is = is; 1092b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1093b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1094b5787286SJed Brown ilink->is_col = is; 1095704ba839SBarry Smith ilink->next = PETSC_NULL; 1096704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10977287d2a3SDmitry Karpeev ierr = KSPIncrementTabLevel(ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1098704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10999005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1100704ba839SBarry Smith 1101a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1102704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1103704ba839SBarry Smith 1104704ba839SBarry Smith if (!next) { 1105704ba839SBarry Smith jac->head = ilink; 1106704ba839SBarry Smith ilink->previous = PETSC_NULL; 1107704ba839SBarry Smith } else { 1108704ba839SBarry Smith while (next->next) { 1109704ba839SBarry Smith next = next->next; 1110704ba839SBarry Smith } 1111704ba839SBarry Smith next->next = ilink; 1112704ba839SBarry Smith ilink->previous = next; 1113704ba839SBarry Smith } 1114704ba839SBarry Smith jac->nsplits++; 1115704ba839SBarry Smith 1116704ba839SBarry Smith PetscFunctionReturn(0); 1117704ba839SBarry Smith } 1118704ba839SBarry Smith EXTERN_C_END 1119704ba839SBarry Smith 11200971522cSBarry Smith #undef __FUNCT__ 11210971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 11220971522cSBarry Smith /*@ 11230971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 11240971522cSBarry Smith 1125ad4df100SBarry Smith Logically Collective on PC 11260971522cSBarry Smith 11270971522cSBarry Smith Input Parameters: 11280971522cSBarry Smith + pc - the preconditioner context 1129a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 11300971522cSBarry Smith . n - the number of fields in this split 1131db4c96c1SJed Brown - fields - the fields in this split 11320971522cSBarry Smith 11330971522cSBarry Smith Level: intermediate 11340971522cSBarry Smith 1135d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1136d32f9abdSBarry Smith 11377287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1138d32f9abdSBarry 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 1139d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1140d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1141d32f9abdSBarry Smith 1142db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1143db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1144db4c96c1SJed Brown 11455d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 11465d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 11475d4c12cdSJungho Lee available when this routine is called. 11485d4c12cdSJungho Lee 1149d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 11500971522cSBarry Smith 11510971522cSBarry Smith @*/ 11525d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11530971522cSBarry Smith { 11544ac538c5SBarry Smith PetscErrorCode ierr; 11550971522cSBarry Smith 11560971522cSBarry Smith PetscFunctionBegin; 11570700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1158db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1159db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1160db4c96c1SJed Brown PetscValidIntPointer(fields,3); 11615d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 11620971522cSBarry Smith PetscFunctionReturn(0); 11630971522cSBarry Smith } 11640971522cSBarry Smith 11650971522cSBarry Smith #undef __FUNCT__ 1166704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1167704ba839SBarry Smith /*@ 1168704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1169704ba839SBarry Smith 1170ad4df100SBarry Smith Logically Collective on PC 1171704ba839SBarry Smith 1172704ba839SBarry Smith Input Parameters: 1173704ba839SBarry Smith + pc - the preconditioner context 1174a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1175db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1176704ba839SBarry Smith 1177d32f9abdSBarry Smith 1178a6ffb8dbSJed Brown Notes: 1179a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1180a6ffb8dbSJed Brown 1181db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1182db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1183d32f9abdSBarry Smith 1184704ba839SBarry Smith Level: intermediate 1185704ba839SBarry Smith 1186704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1187704ba839SBarry Smith 1188704ba839SBarry Smith @*/ 11897087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1190704ba839SBarry Smith { 11914ac538c5SBarry Smith PetscErrorCode ierr; 1192704ba839SBarry Smith 1193704ba839SBarry Smith PetscFunctionBegin; 11940700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11957b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1196db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11974ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1198704ba839SBarry Smith PetscFunctionReturn(0); 1199704ba839SBarry Smith } 1200704ba839SBarry Smith 1201704ba839SBarry Smith #undef __FUNCT__ 120257a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 120357a9adfeSMatthew G Knepley /*@ 120457a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 120557a9adfeSMatthew G Knepley 120657a9adfeSMatthew G Knepley Logically Collective on PC 120757a9adfeSMatthew G Knepley 120857a9adfeSMatthew G Knepley Input Parameters: 120957a9adfeSMatthew G Knepley + pc - the preconditioner context 121057a9adfeSMatthew G Knepley - splitname - name of this split 121157a9adfeSMatthew G Knepley 121257a9adfeSMatthew G Knepley Output Parameter: 121357a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 121457a9adfeSMatthew G Knepley 121557a9adfeSMatthew G Knepley Level: intermediate 121657a9adfeSMatthew G Knepley 121757a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 121857a9adfeSMatthew G Knepley 121957a9adfeSMatthew G Knepley @*/ 122057a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 122157a9adfeSMatthew G Knepley { 122257a9adfeSMatthew G Knepley PetscErrorCode ierr; 122357a9adfeSMatthew G Knepley 122457a9adfeSMatthew G Knepley PetscFunctionBegin; 122557a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 122657a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 122757a9adfeSMatthew G Knepley PetscValidPointer(is,3); 122857a9adfeSMatthew G Knepley { 122957a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 123057a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 123157a9adfeSMatthew G Knepley PetscBool found; 123257a9adfeSMatthew G Knepley 123357a9adfeSMatthew G Knepley *is = PETSC_NULL; 123457a9adfeSMatthew G Knepley while(ilink) { 123557a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 123657a9adfeSMatthew G Knepley if (found) { 123757a9adfeSMatthew G Knepley *is = ilink->is; 123857a9adfeSMatthew G Knepley break; 123957a9adfeSMatthew G Knepley } 124057a9adfeSMatthew G Knepley ilink = ilink->next; 124157a9adfeSMatthew G Knepley } 124257a9adfeSMatthew G Knepley } 124357a9adfeSMatthew G Knepley PetscFunctionReturn(0); 124457a9adfeSMatthew G Knepley } 124557a9adfeSMatthew G Knepley 124657a9adfeSMatthew G Knepley #undef __FUNCT__ 124751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 124851f519a2SBarry Smith /*@ 124951f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 125051f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 125151f519a2SBarry Smith 1252ad4df100SBarry Smith Logically Collective on PC 125351f519a2SBarry Smith 125451f519a2SBarry Smith Input Parameters: 125551f519a2SBarry Smith + pc - the preconditioner context 125651f519a2SBarry Smith - bs - the block size 125751f519a2SBarry Smith 125851f519a2SBarry Smith Level: intermediate 125951f519a2SBarry Smith 126051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 126151f519a2SBarry Smith 126251f519a2SBarry Smith @*/ 12637087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 126451f519a2SBarry Smith { 12654ac538c5SBarry Smith PetscErrorCode ierr; 126651f519a2SBarry Smith 126751f519a2SBarry Smith PetscFunctionBegin; 12680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1269c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 12704ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 127151f519a2SBarry Smith PetscFunctionReturn(0); 127251f519a2SBarry Smith } 127351f519a2SBarry Smith 127451f519a2SBarry Smith #undef __FUNCT__ 127569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 12760971522cSBarry Smith /*@C 127769a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 12780971522cSBarry Smith 127969a612a9SBarry Smith Collective on KSP 12800971522cSBarry Smith 12810971522cSBarry Smith Input Parameter: 12820971522cSBarry Smith . pc - the preconditioner context 12830971522cSBarry Smith 12840971522cSBarry Smith Output Parameters: 128513e0d083SBarry Smith + n - the number of splits 128669a612a9SBarry Smith - pc - the array of KSP contexts 12870971522cSBarry Smith 12880971522cSBarry Smith Note: 1289d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1290d32f9abdSBarry Smith (not the KSP just the array that contains them). 12910971522cSBarry Smith 129269a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12930971522cSBarry Smith 12940971522cSBarry Smith Level: advanced 12950971522cSBarry Smith 12960971522cSBarry Smith .seealso: PCFIELDSPLIT 12970971522cSBarry Smith @*/ 12987087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12990971522cSBarry Smith { 13004ac538c5SBarry Smith PetscErrorCode ierr; 13010971522cSBarry Smith 13020971522cSBarry Smith PetscFunctionBegin; 13030700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 130413e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 13054ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 13060971522cSBarry Smith PetscFunctionReturn(0); 13070971522cSBarry Smith } 13080971522cSBarry Smith 1309e69d4d44SBarry Smith #undef __FUNCT__ 1310e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1311e69d4d44SBarry Smith /*@ 1312e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1313a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1314e69d4d44SBarry Smith 1315e69d4d44SBarry Smith Collective on PC 1316e69d4d44SBarry Smith 1317e69d4d44SBarry Smith Input Parameters: 1318e69d4d44SBarry Smith + pc - the preconditioner context 1319fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1320084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1321084e4875SJed Brown 1322e69d4d44SBarry Smith Options Database: 1323084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1324e69d4d44SBarry Smith 1325fd1303e9SJungho Lee Notes: 1326fd1303e9SJungho Lee $ If ptype is 1327fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1328fd1303e9SJungho Lee $ to this function). 1329fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1330fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1331fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1332fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1333fd1303e9SJungho Lee $ preconditioner 1334fd1303e9SJungho Lee 1335fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1336fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1337fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1338fd1303e9SJungho Lee 1339fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1340fd1303e9SJungho Lee the name since it sets a proceedure to use. 1341fd1303e9SJungho Lee 1342e69d4d44SBarry Smith Level: intermediate 1343e69d4d44SBarry Smith 1344fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1345e69d4d44SBarry Smith 1346e69d4d44SBarry Smith @*/ 13477087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1348e69d4d44SBarry Smith { 13494ac538c5SBarry Smith PetscErrorCode ierr; 1350e69d4d44SBarry Smith 1351e69d4d44SBarry Smith PetscFunctionBegin; 13520700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13534ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1354e69d4d44SBarry Smith PetscFunctionReturn(0); 1355e69d4d44SBarry Smith } 1356e69d4d44SBarry Smith 1357e69d4d44SBarry Smith EXTERN_C_BEGIN 1358e69d4d44SBarry Smith #undef __FUNCT__ 1359e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 13607087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1361e69d4d44SBarry Smith { 1362e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1363084e4875SJed Brown PetscErrorCode ierr; 1364e69d4d44SBarry Smith 1365e69d4d44SBarry Smith PetscFunctionBegin; 1366084e4875SJed Brown jac->schurpre = ptype; 1367084e4875SJed Brown if (pre) { 13686bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1369084e4875SJed Brown jac->schur_user = pre; 1370084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1371084e4875SJed Brown } 1372e69d4d44SBarry Smith PetscFunctionReturn(0); 1373e69d4d44SBarry Smith } 1374e69d4d44SBarry Smith EXTERN_C_END 1375e69d4d44SBarry Smith 137630ad9308SMatthew Knepley #undef __FUNCT__ 1377c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1378ab1df9f5SJed Brown /*@ 1379c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1380ab1df9f5SJed Brown 1381ab1df9f5SJed Brown Collective on PC 1382ab1df9f5SJed Brown 1383ab1df9f5SJed Brown Input Parameters: 1384ab1df9f5SJed Brown + pc - the preconditioner context 1385c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1386ab1df9f5SJed Brown 1387ab1df9f5SJed Brown Options Database: 1388c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1389ab1df9f5SJed Brown 1390ab1df9f5SJed Brown 1391ab1df9f5SJed Brown Level: intermediate 1392ab1df9f5SJed Brown 1393ab1df9f5SJed Brown Notes: 1394ab1df9f5SJed Brown The FULL factorization is 1395ab1df9f5SJed Brown 1396ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1397ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1398ab1df9f5SJed Brown 13996be4592eSBarry 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, 14006be4592eSBarry 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). 1401ab1df9f5SJed Brown 14026be4592eSBarry 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 14036be4592eSBarry 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 14046be4592eSBarry 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, 14056be4592eSBarry 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. 1406ab1df9f5SJed Brown 14076be4592eSBarry 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 14086be4592eSBarry 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). 1409ab1df9f5SJed Brown 1410ab1df9f5SJed Brown References: 1411ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1412ab1df9f5SJed Brown 1413ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1414ab1df9f5SJed Brown 1415ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1416ab1df9f5SJed Brown @*/ 1417c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1418ab1df9f5SJed Brown { 1419ab1df9f5SJed Brown PetscErrorCode ierr; 1420ab1df9f5SJed Brown 1421ab1df9f5SJed Brown PetscFunctionBegin; 1422ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1423c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1424ab1df9f5SJed Brown PetscFunctionReturn(0); 1425ab1df9f5SJed Brown } 1426ab1df9f5SJed Brown 1427ab1df9f5SJed Brown #undef __FUNCT__ 1428c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1429c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1430ab1df9f5SJed Brown { 1431ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1432ab1df9f5SJed Brown 1433ab1df9f5SJed Brown PetscFunctionBegin; 1434ab1df9f5SJed Brown jac->schurfactorization = ftype; 1435ab1df9f5SJed Brown PetscFunctionReturn(0); 1436ab1df9f5SJed Brown } 1437ab1df9f5SJed Brown 1438ab1df9f5SJed Brown #undef __FUNCT__ 143930ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 144030ad9308SMatthew Knepley /*@C 14418c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 144230ad9308SMatthew Knepley 144330ad9308SMatthew Knepley Collective on KSP 144430ad9308SMatthew Knepley 144530ad9308SMatthew Knepley Input Parameter: 144630ad9308SMatthew Knepley . pc - the preconditioner context 144730ad9308SMatthew Knepley 144830ad9308SMatthew Knepley Output Parameters: 1449a04f6461SBarry Smith + A00 - the (0,0) block 1450a04f6461SBarry Smith . A01 - the (0,1) block 1451a04f6461SBarry Smith . A10 - the (1,0) block 1452a04f6461SBarry Smith - A11 - the (1,1) block 145330ad9308SMatthew Knepley 145430ad9308SMatthew Knepley Level: advanced 145530ad9308SMatthew Knepley 145630ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 145730ad9308SMatthew Knepley @*/ 1458a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 145930ad9308SMatthew Knepley { 146030ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 146130ad9308SMatthew Knepley 146230ad9308SMatthew Knepley PetscFunctionBegin; 14630700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1464c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1465a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1466a04f6461SBarry Smith if (A01) *A01 = jac->B; 1467a04f6461SBarry Smith if (A10) *A10 = jac->C; 1468a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 146930ad9308SMatthew Knepley PetscFunctionReturn(0); 147030ad9308SMatthew Knepley } 147130ad9308SMatthew Knepley 147279416396SBarry Smith EXTERN_C_BEGIN 147379416396SBarry Smith #undef __FUNCT__ 147479416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 14757087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 147679416396SBarry Smith { 147779416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1478e69d4d44SBarry Smith PetscErrorCode ierr; 147979416396SBarry Smith 148079416396SBarry Smith PetscFunctionBegin; 148179416396SBarry Smith jac->type = type; 14823b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 14833b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 14843b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1485e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1486e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1487c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1488e69d4d44SBarry Smith 14893b224e63SBarry Smith } else { 14903b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 14913b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1492e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14939e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1494c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 14953b224e63SBarry Smith } 149679416396SBarry Smith PetscFunctionReturn(0); 149779416396SBarry Smith } 149879416396SBarry Smith EXTERN_C_END 149979416396SBarry Smith 150051f519a2SBarry Smith EXTERN_C_BEGIN 150151f519a2SBarry Smith #undef __FUNCT__ 150251f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 15037087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 150451f519a2SBarry Smith { 150551f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 150651f519a2SBarry Smith 150751f519a2SBarry Smith PetscFunctionBegin; 150865e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 150965e19b50SBarry 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); 151051f519a2SBarry Smith jac->bs = bs; 151151f519a2SBarry Smith PetscFunctionReturn(0); 151251f519a2SBarry Smith } 151351f519a2SBarry Smith EXTERN_C_END 151451f519a2SBarry Smith 151579416396SBarry Smith #undef __FUNCT__ 151679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1517bc08b0f1SBarry Smith /*@ 151879416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 151979416396SBarry Smith 152079416396SBarry Smith Collective on PC 152179416396SBarry Smith 152279416396SBarry Smith Input Parameter: 152379416396SBarry Smith . pc - the preconditioner context 152481540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 152579416396SBarry Smith 152679416396SBarry Smith Options Database Key: 1527a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 152879416396SBarry Smith 1529b02e2d75SMatthew G Knepley Level: Intermediate 153079416396SBarry Smith 153179416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 153279416396SBarry Smith 153379416396SBarry Smith .seealso: PCCompositeSetType() 153479416396SBarry Smith 153579416396SBarry Smith @*/ 15367087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 153779416396SBarry Smith { 15384ac538c5SBarry Smith PetscErrorCode ierr; 153979416396SBarry Smith 154079416396SBarry Smith PetscFunctionBegin; 15410700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15424ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 154379416396SBarry Smith PetscFunctionReturn(0); 154479416396SBarry Smith } 154579416396SBarry Smith 1546b02e2d75SMatthew G Knepley #undef __FUNCT__ 1547b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1548b02e2d75SMatthew G Knepley /*@ 1549b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1550b02e2d75SMatthew G Knepley 1551b02e2d75SMatthew G Knepley Not collective 1552b02e2d75SMatthew G Knepley 1553b02e2d75SMatthew G Knepley Input Parameter: 1554b02e2d75SMatthew G Knepley . pc - the preconditioner context 1555b02e2d75SMatthew G Knepley 1556b02e2d75SMatthew G Knepley Output Parameter: 1557b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1558b02e2d75SMatthew G Knepley 1559b02e2d75SMatthew G Knepley Level: Intermediate 1560b02e2d75SMatthew G Knepley 1561b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1562b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1563b02e2d75SMatthew G Knepley @*/ 1564b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1565b02e2d75SMatthew G Knepley { 1566b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1567b02e2d75SMatthew G Knepley 1568b02e2d75SMatthew G Knepley PetscFunctionBegin; 1569b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1570b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1571b02e2d75SMatthew G Knepley *type = jac->type; 1572b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1573b02e2d75SMatthew G Knepley } 1574b02e2d75SMatthew G Knepley 15750971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 15760971522cSBarry Smith /*MC 1577a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1578a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 15790971522cSBarry Smith 1580edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1581edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 15820971522cSBarry Smith 1583a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 158469a612a9SBarry Smith and set the options directly on the resulting KSP object 15850971522cSBarry Smith 15860971522cSBarry Smith Level: intermediate 15870971522cSBarry Smith 158879416396SBarry Smith Options Database Keys: 158981540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 159081540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 159181540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 159281540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 15930f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 15940f188ba9SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag 1595435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 159679416396SBarry Smith 15975d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 15985d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 15995d4c12cdSJungho Lee 1600c8a0d604SMatthew G Knepley Notes: 1601c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1602d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1603d32f9abdSBarry Smith 1604d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1605d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1606d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1607d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1608d32f9abdSBarry Smith 1609c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1610c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1611c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1612c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1613c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1614a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1615c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1616c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1617c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1618a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1619c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1620c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 16215668aaf4SBarry Smith diag gives 1622c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1623c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 16245668aaf4SBarry 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 1625c8a0d604SMatthew G Knepley $ ( A00 0 ) 1626c8a0d604SMatthew G Knepley $ ( A10 S ) 1627c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1628c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1629c8a0d604SMatthew G Knepley $ ( 0 S ) 1630c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1631e69d4d44SBarry Smith 1632edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1633edf189efSBarry Smith is used automatically for a second block. 1634edf189efSBarry Smith 1635ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1636ff218e97SBarry Smith Generally it should be used with the AIJ format. 1637ff218e97SBarry Smith 1638ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1639ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1640ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 16410716a85fSBarry Smith 1642a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 16430971522cSBarry Smith 16447e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1645e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 16460971522cSBarry Smith M*/ 16470971522cSBarry Smith 16480971522cSBarry Smith EXTERN_C_BEGIN 16490971522cSBarry Smith #undef __FUNCT__ 16500971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 16517087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 16520971522cSBarry Smith { 16530971522cSBarry Smith PetscErrorCode ierr; 16540971522cSBarry Smith PC_FieldSplit *jac; 16550971522cSBarry Smith 16560971522cSBarry Smith PetscFunctionBegin; 165738f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 16580971522cSBarry Smith jac->bs = -1; 16590971522cSBarry Smith jac->nsplits = 0; 16603e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1661e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1662c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 166351f519a2SBarry Smith 16640971522cSBarry Smith pc->data = (void*)jac; 16650971522cSBarry Smith 16660971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1667421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 16680971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1669574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 16700971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 16710971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 16720971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 16730971522cSBarry Smith pc->ops->applyrichardson = 0; 16740971522cSBarry Smith 167569a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 167669a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 16770971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 16780971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1679704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1680704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 168179416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 168279416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 168351f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 168451f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 16850971522cSBarry Smith PetscFunctionReturn(0); 16860971522cSBarry Smith } 16870971522cSBarry Smith EXTERN_C_END 16880971522cSBarry Smith 16890971522cSBarry Smith 1690a541d17aSBarry Smith 1691