1dba47a55SKris Buschelman 2c6db04a5SJed Brown #include <private/pcimpl.h> /*I "petscpc.h" I*/ 33c48a1e8SJed Brown #include <petscdmcomposite.h> /*I "petscdmcomposite.h" I*/ 40971522cSBarry Smith 5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0}; 7c5d2311dSJed Brown 8c5d2311dSJed Brown typedef enum { 9c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG, 10c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER, 11c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER, 12c5d2311dSJed Brown PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL 13c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType; 14084e4875SJed Brown 150971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 160971522cSBarry Smith struct _PC_FieldSplitLink { 1769a612a9SBarry Smith KSP ksp; 180971522cSBarry Smith Vec x,y; 19db4c96c1SJed Brown char *splitname; 200971522cSBarry Smith PetscInt nfields; 210971522cSBarry Smith PetscInt *fields; 221b9fc7fcSBarry Smith VecScatter sctx; 234aa3045dSJed Brown IS is; 2451f519a2SBarry Smith PC_FieldSplitLink next,previous; 250971522cSBarry Smith }; 260971522cSBarry Smith 270971522cSBarry Smith typedef struct { 2868dd23aaSBarry Smith PCCompositeType type; 29ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 30ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 31ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 3230ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 3330ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 3479416396SBarry Smith Vec *x,*y,w1,w2; 35519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 36519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3730ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 38ace3abfcSBarry Smith PetscBool issetup; 3930ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 4030ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 4130ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 42a04f6461SBarry Smith Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01 */ 43084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 44084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 45c5d2311dSJed Brown PCFieldSplitSchurFactorizationType schurfactorization; 4630ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 4797bbdb24SBarry Smith PC_FieldSplitLink head; 4863ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 490971522cSBarry Smith } PC_FieldSplit; 500971522cSBarry Smith 5116913363SBarry Smith /* 5216913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5316913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5416913363SBarry Smith PC you could change this. 5516913363SBarry Smith */ 56084e4875SJed Brown 57e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 58084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 59084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 60084e4875SJed Brown { 61084e4875SJed Brown switch (jac->schurpre) { 62084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 63084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 64084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 65084e4875SJed Brown default: 66084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 67084e4875SJed Brown } 68084e4875SJed Brown } 69084e4875SJed Brown 70084e4875SJed Brown 710971522cSBarry Smith #undef __FUNCT__ 720971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 730971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 740971522cSBarry Smith { 750971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 760971522cSBarry Smith PetscErrorCode ierr; 77ace3abfcSBarry Smith PetscBool iascii; 780971522cSBarry Smith PetscInt i,j; 795a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 800971522cSBarry Smith 810971522cSBarry Smith PetscFunctionBegin; 822692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 830971522cSBarry Smith if (iascii) { 8451f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 8569a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 860971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 870971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 881ab39975SBarry Smith if (ilink->fields) { 890971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 9079416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 915a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9279416396SBarry Smith if (j > 0) { 9379416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9479416396SBarry Smith } 955a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 960971522cSBarry Smith } 970971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 991ab39975SBarry Smith } else { 1001ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1011ab39975SBarry Smith } 1025a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1035a9f2f41SSatish Balay ilink = ilink->next; 1040971522cSBarry Smith } 1050971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1060971522cSBarry Smith } else { 10765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1080971522cSBarry Smith } 1090971522cSBarry Smith PetscFunctionReturn(0); 1100971522cSBarry Smith } 1110971522cSBarry Smith 1120971522cSBarry Smith #undef __FUNCT__ 1133b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1143b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1153b224e63SBarry Smith { 1163b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1173b224e63SBarry Smith PetscErrorCode ierr; 118ace3abfcSBarry Smith PetscBool iascii; 1193b224e63SBarry Smith PetscInt i,j; 1203b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1213b224e63SBarry Smith KSP ksp; 1223b224e63SBarry Smith 1233b224e63SBarry Smith PetscFunctionBegin; 1242692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1253b224e63SBarry Smith if (iascii) { 126c5d2311dSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr); 1273b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1283b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1293b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1303b224e63SBarry Smith if (ilink->fields) { 1313b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1323b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1333b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1343b224e63SBarry Smith if (j > 0) { 1353b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1363b224e63SBarry Smith } 1373b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1383b224e63SBarry Smith } 1393b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1403b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1413b224e63SBarry Smith } else { 1423b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1433b224e63SBarry Smith } 1443b224e63SBarry Smith ilink = ilink->next; 1453b224e63SBarry Smith } 146435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1473b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14812cae6f2SJed Brown if (jac->schur) { 14912cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1503b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 15112cae6f2SJed Brown } else { 15212cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 15312cae6f2SJed Brown } 1543b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 155435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1563b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15712cae6f2SJed Brown if (jac->kspschur) { 1583b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 15912cae6f2SJed Brown } else { 16012cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 16112cae6f2SJed Brown } 1623b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1633b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1643b224e63SBarry Smith } else { 16565e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1663b224e63SBarry Smith } 1673b224e63SBarry Smith PetscFunctionReturn(0); 1683b224e63SBarry Smith } 1693b224e63SBarry Smith 1703b224e63SBarry Smith #undef __FUNCT__ 1716c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1726c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1736c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1746c924f48SJed Brown { 1756c924f48SJed Brown PetscErrorCode ierr; 1766c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1776c924f48SJed Brown PetscInt i,nfields,*ifields; 178ace3abfcSBarry Smith PetscBool flg; 1796c924f48SJed Brown char optionname[128],splitname[8]; 1806c924f48SJed Brown 1816c924f48SJed Brown PetscFunctionBegin; 1826c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 1836c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 1846c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 1856c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 1866c924f48SJed Brown nfields = jac->bs; 1876c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 1886c924f48SJed Brown if (!flg) break; 1896c924f48SJed Brown if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1906c924f48SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr); 1916c924f48SJed Brown } 1926c924f48SJed Brown if (i > 0) { 1936c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 1946c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 1956c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 1966c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 1976c924f48SJed Brown } 1986c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 1996c924f48SJed Brown PetscFunctionReturn(0); 2006c924f48SJed Brown } 2016c924f48SJed Brown 2026c924f48SJed Brown #undef __FUNCT__ 20369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 20469a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2050971522cSBarry Smith { 2060971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2070971522cSBarry Smith PetscErrorCode ierr; 2085a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2096ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2106c924f48SJed Brown PetscInt i; 2110971522cSBarry Smith 2120971522cSBarry Smith PetscFunctionBegin; 213d32f9abdSBarry Smith if (!ilink) { 214d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 215d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2168b8307b2SJed Brown PetscBool dmcomposite; 2178b8307b2SJed Brown ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 2188b8307b2SJed Brown if (dmcomposite) { 2198b8307b2SJed Brown PetscInt nDM; 2208b8307b2SJed Brown IS *fields; 2218b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2228b8307b2SJed Brown ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 2238b8307b2SJed Brown ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 2248b8307b2SJed Brown for (i=0; i<nDM; i++) { 2258b8307b2SJed Brown char splitname[8]; 2268b8307b2SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2278b8307b2SJed Brown ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr); 228fcfd50ebSBarry Smith ierr = ISDestroy(&fields[i]);CHKERRQ(ierr); 2298b8307b2SJed Brown } 2308b8307b2SJed Brown ierr = PetscFree(fields);CHKERRQ(ierr); 2318b8307b2SJed Brown } 23266ffff09SJed Brown } else { 233521d7252SBarry Smith if (jac->bs <= 0) { 234704ba839SBarry Smith if (pc->pmat) { 235521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 236704ba839SBarry Smith } else { 237704ba839SBarry Smith jac->bs = 1; 238704ba839SBarry Smith } 239521d7252SBarry Smith } 240d32f9abdSBarry Smith 241acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2426ce1633cSBarry Smith if (stokes) { 2436ce1633cSBarry Smith IS zerodiags,rest; 2446ce1633cSBarry Smith PetscInt nmin,nmax; 2456ce1633cSBarry Smith 2466ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2476ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2486ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2496ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2506ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 251fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 252fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2536ce1633cSBarry Smith } else { 254d32f9abdSBarry Smith if (!flg) { 255d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 256d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2576c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2586c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 259d32f9abdSBarry Smith } 2606c924f48SJed Brown if (flg || !jac->splitdefined) { 261d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 262db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2636c924f48SJed Brown char splitname[8]; 2646c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 265db4c96c1SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr); 26679416396SBarry Smith } 26797bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 268521d7252SBarry Smith } 26966ffff09SJed Brown } 2706ce1633cSBarry Smith } 271edf189efSBarry Smith } else if (jac->nsplits == 1) { 272edf189efSBarry Smith if (ilink->is) { 273edf189efSBarry Smith IS is2; 274edf189efSBarry Smith PetscInt nmin,nmax; 275edf189efSBarry Smith 276edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 277edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 278db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 279fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 280e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 28163ec74ffSBarry Smith } else if (jac->reset) { 28263ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 283d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 284d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 285d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 286d0af7cd3SBarry Smith if (pc->dm && !stokes) { 287d0af7cd3SBarry Smith PetscBool dmcomposite; 288d0af7cd3SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 289d0af7cd3SBarry Smith if (dmcomposite) { 290d0af7cd3SBarry Smith PetscInt nDM; 291d0af7cd3SBarry Smith IS *fields; 292d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 293d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 294d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 295d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 296d0af7cd3SBarry Smith ilink->is = fields[i]; 297d0af7cd3SBarry Smith ilink = ilink->next; 298edf189efSBarry Smith } 299d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 300d0af7cd3SBarry Smith } 301d0af7cd3SBarry Smith } else { 302d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 303d0af7cd3SBarry Smith if (stokes) { 304d0af7cd3SBarry Smith IS zerodiags,rest; 305d0af7cd3SBarry Smith PetscInt nmin,nmax; 306d0af7cd3SBarry Smith 307d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 308d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 309d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 310d0af7cd3SBarry Smith ilink->is = rest; 311d0af7cd3SBarry Smith ilink->next->is = zerodiags; 312d0af7cd3SBarry Smith } else { 313d0af7cd3SBarry Smith SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 314d0af7cd3SBarry Smith } 315d0af7cd3SBarry Smith } 316d0af7cd3SBarry Smith } 317d0af7cd3SBarry Smith 318e7e72b3dSBarry Smith if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 31969a612a9SBarry Smith PetscFunctionReturn(0); 32069a612a9SBarry Smith } 32169a612a9SBarry Smith 32269a612a9SBarry Smith #undef __FUNCT__ 32369a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 32469a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 32569a612a9SBarry Smith { 32669a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 32769a612a9SBarry Smith PetscErrorCode ierr; 3285a9f2f41SSatish Balay PC_FieldSplitLink ilink; 329cf502942SBarry Smith PetscInt i,nsplit,ccsize; 33069a612a9SBarry Smith MatStructure flag = pc->flag; 331ace3abfcSBarry Smith PetscBool sorted; 33269a612a9SBarry Smith 33369a612a9SBarry Smith PetscFunctionBegin; 33469a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 33597bbdb24SBarry Smith nsplit = jac->nsplits; 3365a9f2f41SSatish Balay ilink = jac->head; 33797bbdb24SBarry Smith 33897bbdb24SBarry Smith /* get the matrices for each split */ 339704ba839SBarry Smith if (!jac->issetup) { 3401b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 34197bbdb24SBarry Smith 342704ba839SBarry Smith jac->issetup = PETSC_TRUE; 343704ba839SBarry Smith 344704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 34551f519a2SBarry Smith bs = jac->bs; 34697bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 347cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 3481b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3491b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3501b9fc7fcSBarry Smith if (jac->defaultsplit) { 351704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 352704ba839SBarry Smith } else if (!ilink->is) { 353ccb205f8SBarry Smith if (ilink->nfields > 1) { 3545a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 3555a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3561b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3571b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3581b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 35997bbdb24SBarry Smith } 36097bbdb24SBarry Smith } 361d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 362ccb205f8SBarry Smith } else { 363704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 364ccb205f8SBarry Smith } 3653e197d65SBarry Smith } 366704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 367e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 368704ba839SBarry Smith ilink = ilink->next; 3691b9fc7fcSBarry Smith } 3701b9fc7fcSBarry Smith } 3711b9fc7fcSBarry Smith 372704ba839SBarry Smith ilink = jac->head; 37397bbdb24SBarry Smith if (!jac->pmat) { 374cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 375cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3764aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 377704ba839SBarry Smith ilink = ilink->next; 378cf502942SBarry Smith } 37997bbdb24SBarry Smith } else { 380cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3814aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 382704ba839SBarry Smith ilink = ilink->next; 383cf502942SBarry Smith } 38497bbdb24SBarry Smith } 385519d70e2SJed Brown if (jac->realdiagonal) { 386519d70e2SJed Brown ilink = jac->head; 387519d70e2SJed Brown if (!jac->mat) { 388519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 389519d70e2SJed Brown for (i=0; i<nsplit; i++) { 390519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 391519d70e2SJed Brown ilink = ilink->next; 392519d70e2SJed Brown } 393519d70e2SJed Brown } else { 394519d70e2SJed Brown for (i=0; i<nsplit; i++) { 395966e74d7SJed Brown if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 396519d70e2SJed Brown ilink = ilink->next; 397519d70e2SJed Brown } 398519d70e2SJed Brown } 399519d70e2SJed Brown } else { 400519d70e2SJed Brown jac->mat = jac->pmat; 401519d70e2SJed Brown } 40297bbdb24SBarry Smith 4036c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 40468dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 40568dd23aaSBarry Smith ilink = jac->head; 40668dd23aaSBarry Smith if (!jac->Afield) { 40768dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 40868dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4094aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 41068dd23aaSBarry Smith ilink = ilink->next; 41168dd23aaSBarry Smith } 41268dd23aaSBarry Smith } else { 41368dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4144aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 41568dd23aaSBarry Smith ilink = ilink->next; 41668dd23aaSBarry Smith } 41768dd23aaSBarry Smith } 41868dd23aaSBarry Smith } 41968dd23aaSBarry Smith 4203b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 4213b224e63SBarry Smith IS ccis; 4224aa3045dSJed Brown PetscInt rstart,rend; 423e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 42468dd23aaSBarry Smith 425e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 426e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 427e6cab6aaSJed Brown 4283b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4293b224e63SBarry Smith if (jac->schur) { 4303b224e63SBarry Smith ilink = jac->head; 4314aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4324aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 433fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4343b224e63SBarry Smith ilink = ilink->next; 4354aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4364aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 437fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 438519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 439084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4403b224e63SBarry Smith 4413b224e63SBarry Smith } else { 4421cee3971SBarry Smith KSP ksp; 4436c924f48SJed Brown char schurprefix[256]; 4443b224e63SBarry Smith 445a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4463b224e63SBarry Smith ilink = jac->head; 4474aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4484aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 449fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4503b224e63SBarry Smith ilink = ilink->next; 4514aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4524aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 453fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4547d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 455519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 456a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4571cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 458e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 459a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 460a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 4613b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4623b224e63SBarry Smith 4633b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4649005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4651cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 466084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 467084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 468084e4875SJed Brown PC pc; 469cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 470084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 471084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 472e69d4d44SBarry Smith } 4736c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 4746c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 4753b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 4763b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 4773b224e63SBarry Smith 4783b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 4793b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 4803b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 4813b224e63SBarry Smith ilink = jac->head; 4823b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 4833b224e63SBarry Smith ilink = ilink->next; 4843b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 4853b224e63SBarry Smith } 4863b224e63SBarry Smith } else { 48797bbdb24SBarry Smith /* set up the individual PCs */ 48897bbdb24SBarry Smith i = 0; 4895a9f2f41SSatish Balay ilink = jac->head; 4905a9f2f41SSatish Balay while (ilink) { 491519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 4923b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 4935a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 4945a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 49597bbdb24SBarry Smith i++; 4965a9f2f41SSatish Balay ilink = ilink->next; 4970971522cSBarry Smith } 49897bbdb24SBarry Smith 49997bbdb24SBarry Smith /* create work vectors for each split */ 5001b9fc7fcSBarry Smith if (!jac->x) { 50197bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5025a9f2f41SSatish Balay ilink = jac->head; 50397bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 504906ed7ccSBarry Smith Vec *vl,*vr; 505906ed7ccSBarry Smith 506906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 507906ed7ccSBarry Smith ilink->x = *vr; 508906ed7ccSBarry Smith ilink->y = *vl; 509906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 510906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5115a9f2f41SSatish Balay jac->x[i] = ilink->x; 5125a9f2f41SSatish Balay jac->y[i] = ilink->y; 5135a9f2f41SSatish Balay ilink = ilink->next; 51497bbdb24SBarry Smith } 5153b224e63SBarry Smith } 5163b224e63SBarry Smith } 5173b224e63SBarry Smith 5183b224e63SBarry Smith 519d0f46423SBarry Smith if (!jac->head->sctx) { 5203b224e63SBarry Smith Vec xtmp; 5213b224e63SBarry Smith 52279416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5231b9fc7fcSBarry Smith 5245a9f2f41SSatish Balay ilink = jac->head; 5251b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5261b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 527704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5285a9f2f41SSatish Balay ilink = ilink->next; 5291b9fc7fcSBarry Smith } 5306bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5311b9fc7fcSBarry Smith } 5320971522cSBarry Smith PetscFunctionReturn(0); 5330971522cSBarry Smith } 5340971522cSBarry Smith 5355a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 536ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 537ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5385a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 539ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 540ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 54179416396SBarry Smith 5420971522cSBarry Smith #undef __FUNCT__ 5433b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5443b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5453b224e63SBarry Smith { 5463b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5473b224e63SBarry Smith PetscErrorCode ierr; 5483b224e63SBarry Smith KSP ksp; 5493b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5503b224e63SBarry Smith 5513b224e63SBarry Smith PetscFunctionBegin; 5523b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5533b224e63SBarry Smith 554c5d2311dSJed Brown switch (jac->schurfactorization) { 555c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 556a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 557c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 558c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 559c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 560c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 561c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 562c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 563c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 564c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 565c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 566c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 567c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 568c5d2311dSJed Brown break; 569c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 570a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 571c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 572c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 573c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 574c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 575c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 576c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 577c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 578c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 580c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 581c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 582c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 583c5d2311dSJed Brown break; 584c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 585a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 586c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 587c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 588c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 589c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 590c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 591c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 592c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 593c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 594c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 595c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 596c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 598c5d2311dSJed Brown break; 599c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 600a04f6461SBarry 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 */ 6013b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6023b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6033b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6043b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6053b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6063b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6073b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6083b224e63SBarry Smith 6093b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6103b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6113b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6123b224e63SBarry Smith 6133b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6143b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6153b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6163b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6173b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 618c5d2311dSJed Brown } 6193b224e63SBarry Smith PetscFunctionReturn(0); 6203b224e63SBarry Smith } 6213b224e63SBarry Smith 6223b224e63SBarry Smith #undef __FUNCT__ 6230971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6240971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6250971522cSBarry Smith { 6260971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6270971522cSBarry Smith PetscErrorCode ierr; 6285a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 629af93645bSJed Brown PetscInt cnt; 6300971522cSBarry Smith 6310971522cSBarry Smith PetscFunctionBegin; 63251f519a2SBarry Smith CHKMEMQ; 63351f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 63451f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 63551f519a2SBarry Smith 63679416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6371b9fc7fcSBarry Smith if (jac->defaultsplit) { 6380971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6395a9f2f41SSatish Balay while (ilink) { 6405a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6415a9f2f41SSatish Balay ilink = ilink->next; 6420971522cSBarry Smith } 6430971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6441b9fc7fcSBarry Smith } else { 645efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6465a9f2f41SSatish Balay while (ilink) { 6475a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6485a9f2f41SSatish Balay ilink = ilink->next; 6491b9fc7fcSBarry Smith } 6501b9fc7fcSBarry Smith } 65116913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 65279416396SBarry Smith if (!jac->w1) { 65379416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 65479416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 65579416396SBarry Smith } 656efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6575a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6583e197d65SBarry Smith cnt = 1; 6595a9f2f41SSatish Balay while (ilink->next) { 6605a9f2f41SSatish Balay ilink = ilink->next; 6613e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6623e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6633e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6643e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6653e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6663e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6673e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6683e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6693e197d65SBarry Smith } 67051f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 67111755939SBarry Smith cnt -= 2; 67251f519a2SBarry Smith while (ilink->previous) { 67351f519a2SBarry Smith ilink = ilink->previous; 67411755939SBarry Smith /* compute the residual only over the part of the vector needed */ 67511755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 67611755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 67711755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 67811755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 67911755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 68011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 68111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 68251f519a2SBarry Smith } 68311755939SBarry Smith } 68465e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 68551f519a2SBarry Smith CHKMEMQ; 6860971522cSBarry Smith PetscFunctionReturn(0); 6870971522cSBarry Smith } 6880971522cSBarry Smith 689421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 690ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 691ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 692421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 693ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 694ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 695421e10b8SBarry Smith 696421e10b8SBarry Smith #undef __FUNCT__ 6978c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 698421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 699421e10b8SBarry Smith { 700421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 701421e10b8SBarry Smith PetscErrorCode ierr; 702421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 703421e10b8SBarry Smith 704421e10b8SBarry Smith PetscFunctionBegin; 705421e10b8SBarry Smith CHKMEMQ; 706421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 707421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 708421e10b8SBarry Smith 709421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 710421e10b8SBarry Smith if (jac->defaultsplit) { 711421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 712421e10b8SBarry Smith while (ilink) { 713421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 714421e10b8SBarry Smith ilink = ilink->next; 715421e10b8SBarry Smith } 716421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 717421e10b8SBarry Smith } else { 718421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 719421e10b8SBarry Smith while (ilink) { 720421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 721421e10b8SBarry Smith ilink = ilink->next; 722421e10b8SBarry Smith } 723421e10b8SBarry Smith } 724421e10b8SBarry Smith } else { 725421e10b8SBarry Smith if (!jac->w1) { 726421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 727421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 728421e10b8SBarry Smith } 729421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 730421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 731421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 732421e10b8SBarry Smith while (ilink->next) { 733421e10b8SBarry Smith ilink = ilink->next; 7349989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 735421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 736421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 737421e10b8SBarry Smith } 738421e10b8SBarry Smith while (ilink->previous) { 739421e10b8SBarry Smith ilink = ilink->previous; 7409989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 741421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 742421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 743421e10b8SBarry Smith } 744421e10b8SBarry Smith } else { 745421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 746421e10b8SBarry Smith ilink = ilink->next; 747421e10b8SBarry Smith } 748421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 749421e10b8SBarry Smith while (ilink->previous) { 750421e10b8SBarry Smith ilink = ilink->previous; 7519989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 752421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 753421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 754421e10b8SBarry Smith } 755421e10b8SBarry Smith } 756421e10b8SBarry Smith } 757421e10b8SBarry Smith CHKMEMQ; 758421e10b8SBarry Smith PetscFunctionReturn(0); 759421e10b8SBarry Smith } 760421e10b8SBarry Smith 7610971522cSBarry Smith #undef __FUNCT__ 762574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 763574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 7640971522cSBarry Smith { 7650971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7660971522cSBarry Smith PetscErrorCode ierr; 7675a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 7680971522cSBarry Smith 7690971522cSBarry Smith PetscFunctionBegin; 7705a9f2f41SSatish Balay while (ilink) { 771574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 772fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 773fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 774fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 775fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 7765a9f2f41SSatish Balay next = ilink->next; 7775a9f2f41SSatish Balay ilink = next; 7780971522cSBarry Smith } 77905b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 780574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 781574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 782574deadeSBarry Smith } else if (jac->mat) { 783574deadeSBarry Smith jac->mat = PETSC_NULL; 784574deadeSBarry Smith } 78597bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 78668dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 7876bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 7886bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 7896bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 7906bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 7916bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 7926bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 7936bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 79463ec74ffSBarry Smith jac->reset = PETSC_TRUE; 795574deadeSBarry Smith PetscFunctionReturn(0); 796574deadeSBarry Smith } 797574deadeSBarry Smith 798574deadeSBarry Smith #undef __FUNCT__ 799574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 800574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 801574deadeSBarry Smith { 802574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 803574deadeSBarry Smith PetscErrorCode ierr; 804574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 805574deadeSBarry Smith 806574deadeSBarry Smith PetscFunctionBegin; 807574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 808574deadeSBarry Smith while (ilink) { 8096bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 810574deadeSBarry Smith next = ilink->next; 811574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 812574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 813574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 814574deadeSBarry Smith ilink = next; 815574deadeSBarry Smith } 816574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 817c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8180971522cSBarry Smith PetscFunctionReturn(0); 8190971522cSBarry Smith } 8200971522cSBarry Smith 8210971522cSBarry Smith #undef __FUNCT__ 8220971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8230971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8240971522cSBarry Smith { 8251b9fc7fcSBarry Smith PetscErrorCode ierr; 8266c924f48SJed Brown PetscInt bs; 827bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8289dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8293b224e63SBarry Smith PCCompositeType ctype; 8301b9fc7fcSBarry Smith 8310971522cSBarry Smith PetscFunctionBegin; 8321b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 833acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 83451f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 83551f519a2SBarry Smith if (flg) { 83651f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 83751f519a2SBarry Smith } 838704ba839SBarry Smith 839435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 840c0adfefeSBarry Smith if (stokes) { 841c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 842c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 843c0adfefeSBarry Smith } 844c0adfefeSBarry Smith 8453b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 8463b224e63SBarry Smith if (flg) { 8473b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 8483b224e63SBarry Smith } 849d32f9abdSBarry Smith 850c30613efSMatthew Knepley /* Only setup fields once */ 851c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 852d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 853d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 8546c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 8556c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 856d32f9abdSBarry Smith } 857c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 858c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 859084e4875SJed 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); 860c5d2311dSJed Brown } 8611b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8620971522cSBarry Smith PetscFunctionReturn(0); 8630971522cSBarry Smith } 8640971522cSBarry Smith 8650971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 8660971522cSBarry Smith 8670971522cSBarry Smith EXTERN_C_BEGIN 8680971522cSBarry Smith #undef __FUNCT__ 8690971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 8707087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 8710971522cSBarry Smith { 87297bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8730971522cSBarry Smith PetscErrorCode ierr; 8745a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 87569a612a9SBarry Smith char prefix[128]; 87651f519a2SBarry Smith PetscInt i; 8770971522cSBarry Smith 8780971522cSBarry Smith PetscFunctionBegin; 8796c924f48SJed Brown if (jac->splitdefined) { 8806c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 8816c924f48SJed Brown PetscFunctionReturn(0); 8826c924f48SJed Brown } 88351f519a2SBarry Smith for (i=0; i<n; i++) { 884e32f2f54SBarry 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); 885e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 88651f519a2SBarry Smith } 887704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 888a04f6461SBarry Smith if (splitname) { 889db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 890a04f6461SBarry Smith } else { 891a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 892a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 893a04f6461SBarry Smith } 894704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 8955a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 8965a9f2f41SSatish Balay ilink->nfields = n; 8975a9f2f41SSatish Balay ilink->next = PETSC_NULL; 8987adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 8991cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9005a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9019005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 90269a612a9SBarry Smith 903a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9045a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9050971522cSBarry Smith 9060971522cSBarry Smith if (!next) { 9075a9f2f41SSatish Balay jac->head = ilink; 90851f519a2SBarry Smith ilink->previous = PETSC_NULL; 9090971522cSBarry Smith } else { 9100971522cSBarry Smith while (next->next) { 9110971522cSBarry Smith next = next->next; 9120971522cSBarry Smith } 9135a9f2f41SSatish Balay next->next = ilink; 91451f519a2SBarry Smith ilink->previous = next; 9150971522cSBarry Smith } 9160971522cSBarry Smith jac->nsplits++; 9170971522cSBarry Smith PetscFunctionReturn(0); 9180971522cSBarry Smith } 9190971522cSBarry Smith EXTERN_C_END 9200971522cSBarry Smith 921e69d4d44SBarry Smith EXTERN_C_BEGIN 922e69d4d44SBarry Smith #undef __FUNCT__ 923e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9247087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 925e69d4d44SBarry Smith { 926e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 927e69d4d44SBarry Smith PetscErrorCode ierr; 928e69d4d44SBarry Smith 929e69d4d44SBarry Smith PetscFunctionBegin; 930519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 931e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 932e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 93313e0d083SBarry Smith if (n) *n = jac->nsplits; 934e69d4d44SBarry Smith PetscFunctionReturn(0); 935e69d4d44SBarry Smith } 936e69d4d44SBarry Smith EXTERN_C_END 9370971522cSBarry Smith 9380971522cSBarry Smith EXTERN_C_BEGIN 9390971522cSBarry Smith #undef __FUNCT__ 94069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 9417087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 9420971522cSBarry Smith { 9430971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9440971522cSBarry Smith PetscErrorCode ierr; 9450971522cSBarry Smith PetscInt cnt = 0; 9465a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 9470971522cSBarry Smith 9480971522cSBarry Smith PetscFunctionBegin; 9495d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 9505a9f2f41SSatish Balay while (ilink) { 9515a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 9525a9f2f41SSatish Balay ilink = ilink->next; 9530971522cSBarry Smith } 9545d480477SMatthew 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); 95513e0d083SBarry Smith if (n) *n = jac->nsplits; 9560971522cSBarry Smith PetscFunctionReturn(0); 9570971522cSBarry Smith } 9580971522cSBarry Smith EXTERN_C_END 9590971522cSBarry Smith 960704ba839SBarry Smith EXTERN_C_BEGIN 961704ba839SBarry Smith #undef __FUNCT__ 962704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 9637087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 964704ba839SBarry Smith { 965704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 966704ba839SBarry Smith PetscErrorCode ierr; 967704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 968704ba839SBarry Smith char prefix[128]; 969704ba839SBarry Smith 970704ba839SBarry Smith PetscFunctionBegin; 9716c924f48SJed Brown if (jac->splitdefined) { 9726c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9736c924f48SJed Brown PetscFunctionReturn(0); 9746c924f48SJed Brown } 97516913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 976a04f6461SBarry Smith if (splitname) { 977db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 978a04f6461SBarry Smith } else { 979a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 980a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 981a04f6461SBarry Smith } 9821ab39975SBarry Smith ilink->is = is; 9831ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 984704ba839SBarry Smith ilink->next = PETSC_NULL; 985704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9861cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 987704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9889005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 989704ba839SBarry Smith 990a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 991704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 992704ba839SBarry Smith 993704ba839SBarry Smith if (!next) { 994704ba839SBarry Smith jac->head = ilink; 995704ba839SBarry Smith ilink->previous = PETSC_NULL; 996704ba839SBarry Smith } else { 997704ba839SBarry Smith while (next->next) { 998704ba839SBarry Smith next = next->next; 999704ba839SBarry Smith } 1000704ba839SBarry Smith next->next = ilink; 1001704ba839SBarry Smith ilink->previous = next; 1002704ba839SBarry Smith } 1003704ba839SBarry Smith jac->nsplits++; 1004704ba839SBarry Smith 1005704ba839SBarry Smith PetscFunctionReturn(0); 1006704ba839SBarry Smith } 1007704ba839SBarry Smith EXTERN_C_END 1008704ba839SBarry Smith 10090971522cSBarry Smith #undef __FUNCT__ 10100971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10110971522cSBarry Smith /*@ 10120971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10130971522cSBarry Smith 1014ad4df100SBarry Smith Logically Collective on PC 10150971522cSBarry Smith 10160971522cSBarry Smith Input Parameters: 10170971522cSBarry Smith + pc - the preconditioner context 1018a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10190971522cSBarry Smith . n - the number of fields in this split 1020db4c96c1SJed Brown - fields - the fields in this split 10210971522cSBarry Smith 10220971522cSBarry Smith Level: intermediate 10230971522cSBarry Smith 1024d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1025d32f9abdSBarry Smith 1026d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1027d32f9abdSBarry 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 1028d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1029d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1030d32f9abdSBarry Smith 1031db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1032db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1033db4c96c1SJed Brown 1034d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 10350971522cSBarry Smith 10360971522cSBarry Smith @*/ 10377087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 10380971522cSBarry Smith { 10394ac538c5SBarry Smith PetscErrorCode ierr; 10400971522cSBarry Smith 10410971522cSBarry Smith PetscFunctionBegin; 10420700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1043db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1044db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1045db4c96c1SJed Brown PetscValidIntPointer(fields,3); 10464ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr); 10470971522cSBarry Smith PetscFunctionReturn(0); 10480971522cSBarry Smith } 10490971522cSBarry Smith 10500971522cSBarry Smith #undef __FUNCT__ 1051704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1052704ba839SBarry Smith /*@ 1053704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1054704ba839SBarry Smith 1055ad4df100SBarry Smith Logically Collective on PC 1056704ba839SBarry Smith 1057704ba839SBarry Smith Input Parameters: 1058704ba839SBarry Smith + pc - the preconditioner context 1059a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1060db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1061704ba839SBarry Smith 1062d32f9abdSBarry Smith 1063a6ffb8dbSJed Brown Notes: 1064a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1065a6ffb8dbSJed Brown 1066db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1067db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1068d32f9abdSBarry Smith 1069704ba839SBarry Smith Level: intermediate 1070704ba839SBarry Smith 1071704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1072704ba839SBarry Smith 1073704ba839SBarry Smith @*/ 10747087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1075704ba839SBarry Smith { 10764ac538c5SBarry Smith PetscErrorCode ierr; 1077704ba839SBarry Smith 1078704ba839SBarry Smith PetscFunctionBegin; 10790700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1080db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1081db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 10824ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1083704ba839SBarry Smith PetscFunctionReturn(0); 1084704ba839SBarry Smith } 1085704ba839SBarry Smith 1086704ba839SBarry Smith #undef __FUNCT__ 108757a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 108857a9adfeSMatthew G Knepley /*@ 108957a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 109057a9adfeSMatthew G Knepley 109157a9adfeSMatthew G Knepley Logically Collective on PC 109257a9adfeSMatthew G Knepley 109357a9adfeSMatthew G Knepley Input Parameters: 109457a9adfeSMatthew G Knepley + pc - the preconditioner context 109557a9adfeSMatthew G Knepley - splitname - name of this split 109657a9adfeSMatthew G Knepley 109757a9adfeSMatthew G Knepley Output Parameter: 109857a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 109957a9adfeSMatthew G Knepley 110057a9adfeSMatthew G Knepley Level: intermediate 110157a9adfeSMatthew G Knepley 110257a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 110357a9adfeSMatthew G Knepley 110457a9adfeSMatthew G Knepley @*/ 110557a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 110657a9adfeSMatthew G Knepley { 110757a9adfeSMatthew G Knepley PetscErrorCode ierr; 110857a9adfeSMatthew G Knepley 110957a9adfeSMatthew G Knepley PetscFunctionBegin; 111057a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 111157a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 111257a9adfeSMatthew G Knepley PetscValidPointer(is,3); 111357a9adfeSMatthew G Knepley { 111457a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 111557a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 111657a9adfeSMatthew G Knepley PetscBool found; 111757a9adfeSMatthew G Knepley 111857a9adfeSMatthew G Knepley *is = PETSC_NULL; 111957a9adfeSMatthew G Knepley while(ilink) { 112057a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 112157a9adfeSMatthew G Knepley if (found) { 112257a9adfeSMatthew G Knepley *is = ilink->is; 112357a9adfeSMatthew G Knepley break; 112457a9adfeSMatthew G Knepley } 112557a9adfeSMatthew G Knepley ilink = ilink->next; 112657a9adfeSMatthew G Knepley } 112757a9adfeSMatthew G Knepley } 112857a9adfeSMatthew G Knepley PetscFunctionReturn(0); 112957a9adfeSMatthew G Knepley } 113057a9adfeSMatthew G Knepley 113157a9adfeSMatthew G Knepley #undef __FUNCT__ 113251f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 113351f519a2SBarry Smith /*@ 113451f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 113551f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 113651f519a2SBarry Smith 1137ad4df100SBarry Smith Logically Collective on PC 113851f519a2SBarry Smith 113951f519a2SBarry Smith Input Parameters: 114051f519a2SBarry Smith + pc - the preconditioner context 114151f519a2SBarry Smith - bs - the block size 114251f519a2SBarry Smith 114351f519a2SBarry Smith Level: intermediate 114451f519a2SBarry Smith 114551f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 114651f519a2SBarry Smith 114751f519a2SBarry Smith @*/ 11487087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 114951f519a2SBarry Smith { 11504ac538c5SBarry Smith PetscErrorCode ierr; 115151f519a2SBarry Smith 115251f519a2SBarry Smith PetscFunctionBegin; 11530700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1154c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 11554ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 115651f519a2SBarry Smith PetscFunctionReturn(0); 115751f519a2SBarry Smith } 115851f519a2SBarry Smith 115951f519a2SBarry Smith #undef __FUNCT__ 116069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 11610971522cSBarry Smith /*@C 116269a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 11630971522cSBarry Smith 116469a612a9SBarry Smith Collective on KSP 11650971522cSBarry Smith 11660971522cSBarry Smith Input Parameter: 11670971522cSBarry Smith . pc - the preconditioner context 11680971522cSBarry Smith 11690971522cSBarry Smith Output Parameters: 117013e0d083SBarry Smith + n - the number of splits 117169a612a9SBarry Smith - pc - the array of KSP contexts 11720971522cSBarry Smith 11730971522cSBarry Smith Note: 1174d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1175d32f9abdSBarry Smith (not the KSP just the array that contains them). 11760971522cSBarry Smith 117769a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 11780971522cSBarry Smith 11790971522cSBarry Smith Level: advanced 11800971522cSBarry Smith 11810971522cSBarry Smith .seealso: PCFIELDSPLIT 11820971522cSBarry Smith @*/ 11837087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 11840971522cSBarry Smith { 11854ac538c5SBarry Smith PetscErrorCode ierr; 11860971522cSBarry Smith 11870971522cSBarry Smith PetscFunctionBegin; 11880700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 118913e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 11904ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 11910971522cSBarry Smith PetscFunctionReturn(0); 11920971522cSBarry Smith } 11930971522cSBarry Smith 1194e69d4d44SBarry Smith #undef __FUNCT__ 1195e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1196e69d4d44SBarry Smith /*@ 1197e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1198a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1199e69d4d44SBarry Smith 1200e69d4d44SBarry Smith Collective on PC 1201e69d4d44SBarry Smith 1202e69d4d44SBarry Smith Input Parameters: 1203e69d4d44SBarry Smith + pc - the preconditioner context 1204*fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1205084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1206084e4875SJed Brown 1207e69d4d44SBarry Smith Options Database: 1208084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1209e69d4d44SBarry Smith 1210*fd1303e9SJungho Lee Notes: 1211*fd1303e9SJungho Lee $ If ptype is 1212*fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1213*fd1303e9SJungho Lee $ to this function). 1214*fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1215*fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1216*fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1217*fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1218*fd1303e9SJungho Lee $ preconditioner 1219*fd1303e9SJungho Lee 1220*fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1221*fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1222*fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1223*fd1303e9SJungho Lee 1224*fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1225*fd1303e9SJungho Lee the name since it sets a proceedure to use. 1226*fd1303e9SJungho Lee 1227e69d4d44SBarry Smith Level: intermediate 1228e69d4d44SBarry Smith 1229*fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1230e69d4d44SBarry Smith 1231e69d4d44SBarry Smith @*/ 12327087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1233e69d4d44SBarry Smith { 12344ac538c5SBarry Smith PetscErrorCode ierr; 1235e69d4d44SBarry Smith 1236e69d4d44SBarry Smith PetscFunctionBegin; 12370700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12384ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1239e69d4d44SBarry Smith PetscFunctionReturn(0); 1240e69d4d44SBarry Smith } 1241e69d4d44SBarry Smith 1242e69d4d44SBarry Smith EXTERN_C_BEGIN 1243e69d4d44SBarry Smith #undef __FUNCT__ 1244e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 12457087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1246e69d4d44SBarry Smith { 1247e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1248084e4875SJed Brown PetscErrorCode ierr; 1249e69d4d44SBarry Smith 1250e69d4d44SBarry Smith PetscFunctionBegin; 1251084e4875SJed Brown jac->schurpre = ptype; 1252084e4875SJed Brown if (pre) { 12536bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1254084e4875SJed Brown jac->schur_user = pre; 1255084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1256084e4875SJed Brown } 1257e69d4d44SBarry Smith PetscFunctionReturn(0); 1258e69d4d44SBarry Smith } 1259e69d4d44SBarry Smith EXTERN_C_END 1260e69d4d44SBarry Smith 126130ad9308SMatthew Knepley #undef __FUNCT__ 126230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 126330ad9308SMatthew Knepley /*@C 12648c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 126530ad9308SMatthew Knepley 126630ad9308SMatthew Knepley Collective on KSP 126730ad9308SMatthew Knepley 126830ad9308SMatthew Knepley Input Parameter: 126930ad9308SMatthew Knepley . pc - the preconditioner context 127030ad9308SMatthew Knepley 127130ad9308SMatthew Knepley Output Parameters: 1272a04f6461SBarry Smith + A00 - the (0,0) block 1273a04f6461SBarry Smith . A01 - the (0,1) block 1274a04f6461SBarry Smith . A10 - the (1,0) block 1275a04f6461SBarry Smith - A11 - the (1,1) block 127630ad9308SMatthew Knepley 127730ad9308SMatthew Knepley Level: advanced 127830ad9308SMatthew Knepley 127930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 128030ad9308SMatthew Knepley @*/ 1281a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 128230ad9308SMatthew Knepley { 128330ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 128430ad9308SMatthew Knepley 128530ad9308SMatthew Knepley PetscFunctionBegin; 12860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1287c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1288a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1289a04f6461SBarry Smith if (A01) *A01 = jac->B; 1290a04f6461SBarry Smith if (A10) *A10 = jac->C; 1291a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 129230ad9308SMatthew Knepley PetscFunctionReturn(0); 129330ad9308SMatthew Knepley } 129430ad9308SMatthew Knepley 129579416396SBarry Smith EXTERN_C_BEGIN 129679416396SBarry Smith #undef __FUNCT__ 129779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 12987087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 129979416396SBarry Smith { 130079416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1301e69d4d44SBarry Smith PetscErrorCode ierr; 130279416396SBarry Smith 130379416396SBarry Smith PetscFunctionBegin; 130479416396SBarry Smith jac->type = type; 13053b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 13063b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 13073b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1308e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1309e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1310e69d4d44SBarry Smith 13113b224e63SBarry Smith } else { 13123b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 13133b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1314e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13159e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 13163b224e63SBarry Smith } 131779416396SBarry Smith PetscFunctionReturn(0); 131879416396SBarry Smith } 131979416396SBarry Smith EXTERN_C_END 132079416396SBarry Smith 132151f519a2SBarry Smith EXTERN_C_BEGIN 132251f519a2SBarry Smith #undef __FUNCT__ 132351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 13247087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 132551f519a2SBarry Smith { 132651f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 132751f519a2SBarry Smith 132851f519a2SBarry Smith PetscFunctionBegin; 132965e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 133065e19b50SBarry 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); 133151f519a2SBarry Smith jac->bs = bs; 133251f519a2SBarry Smith PetscFunctionReturn(0); 133351f519a2SBarry Smith } 133451f519a2SBarry Smith EXTERN_C_END 133551f519a2SBarry Smith 133679416396SBarry Smith #undef __FUNCT__ 133779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1338bc08b0f1SBarry Smith /*@ 133979416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 134079416396SBarry Smith 134179416396SBarry Smith Collective on PC 134279416396SBarry Smith 134379416396SBarry Smith Input Parameter: 134479416396SBarry Smith . pc - the preconditioner context 134581540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 134679416396SBarry Smith 134779416396SBarry Smith Options Database Key: 1348a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 134979416396SBarry Smith 135079416396SBarry Smith Level: Developer 135179416396SBarry Smith 135279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 135379416396SBarry Smith 135479416396SBarry Smith .seealso: PCCompositeSetType() 135579416396SBarry Smith 135679416396SBarry Smith @*/ 13577087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 135879416396SBarry Smith { 13594ac538c5SBarry Smith PetscErrorCode ierr; 136079416396SBarry Smith 136179416396SBarry Smith PetscFunctionBegin; 13620700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13634ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 136479416396SBarry Smith PetscFunctionReturn(0); 136579416396SBarry Smith } 136679416396SBarry Smith 13670971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 13680971522cSBarry Smith /*MC 1369a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1370a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 13710971522cSBarry Smith 1372edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1373edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 13740971522cSBarry Smith 1375a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 137669a612a9SBarry Smith and set the options directly on the resulting KSP object 13770971522cSBarry Smith 13780971522cSBarry Smith Level: intermediate 13790971522cSBarry Smith 138079416396SBarry Smith Options Database Keys: 138181540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 138281540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 138381540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 138481540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 138581540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1386e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1387435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 138879416396SBarry Smith 1389edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 13903b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 13913b224e63SBarry Smith 1392d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1393d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1394d32f9abdSBarry Smith 1395d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1396d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1397d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1398d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1399d32f9abdSBarry Smith 1400a04f6461SBarry Smith For the Schur complement preconditioner if J = ( A00 A01 ) 1401a04f6461SBarry Smith ( A10 A11 ) 1402e69d4d44SBarry Smith the preconditioner is 1403a04f6461SBarry Smith (I -B ksp(A00)) ( inv(A00) 0 (I 0 ) 1404a04f6461SBarry Smith (0 I ) ( 0 ksp(S) ) (-A10 ksp(A00) I ) 1405a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1406a04f6461SBarry Smith ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given in providing the SECOND split or 1 if not give). 1407a04f6461SBarry Smith For PCFieldSplitGetKSP() when field number is 1408edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1409a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 14107e8cb189SBarry Smith option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc 1411e69d4d44SBarry Smith 1412edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1413edf189efSBarry Smith is used automatically for a second block. 1414edf189efSBarry Smith 14150716a85fSBarry Smith The fieldsplit preconditioner cannot be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. Generally it should be used with the AIJ format. 14160716a85fSBarry Smith 1417a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 14180971522cSBarry Smith 14197e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1420e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 14210971522cSBarry Smith M*/ 14220971522cSBarry Smith 14230971522cSBarry Smith EXTERN_C_BEGIN 14240971522cSBarry Smith #undef __FUNCT__ 14250971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 14267087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 14270971522cSBarry Smith { 14280971522cSBarry Smith PetscErrorCode ierr; 14290971522cSBarry Smith PC_FieldSplit *jac; 14300971522cSBarry Smith 14310971522cSBarry Smith PetscFunctionBegin; 143238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 14330971522cSBarry Smith jac->bs = -1; 14340971522cSBarry Smith jac->nsplits = 0; 14353e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1436e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1437c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 143851f519a2SBarry Smith 14390971522cSBarry Smith pc->data = (void*)jac; 14400971522cSBarry Smith 14410971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1442421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 14430971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1444574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 14450971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 14460971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 14470971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 14480971522cSBarry Smith pc->ops->applyrichardson = 0; 14490971522cSBarry Smith 145069a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 145169a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14520971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 14530971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1454704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1455704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 145679416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 145779416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 145851f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 145951f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 14600971522cSBarry Smith PetscFunctionReturn(0); 14610971522cSBarry Smith } 14620971522cSBarry Smith EXTERN_C_END 14630971522cSBarry Smith 14640971522cSBarry Smith 1465a541d17aSBarry Smith 1466