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 */ 49*c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 500971522cSBarry Smith } PC_FieldSplit; 510971522cSBarry Smith 5216913363SBarry Smith /* 5316913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5416913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5516913363SBarry Smith PC you could change this. 5616913363SBarry Smith */ 57084e4875SJed Brown 58e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 59084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 60084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 61084e4875SJed Brown { 62084e4875SJed Brown switch (jac->schurpre) { 63084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 64084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 65084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 66084e4875SJed Brown default: 67084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 68084e4875SJed Brown } 69084e4875SJed Brown } 70084e4875SJed Brown 71084e4875SJed Brown 720971522cSBarry Smith #undef __FUNCT__ 730971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 740971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 750971522cSBarry Smith { 760971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 770971522cSBarry Smith PetscErrorCode ierr; 78ace3abfcSBarry Smith PetscBool iascii; 790971522cSBarry Smith PetscInt i,j; 805a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 810971522cSBarry Smith 820971522cSBarry Smith PetscFunctionBegin; 832692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 840971522cSBarry Smith if (iascii) { 8551f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 8669a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 870971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 880971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 891ab39975SBarry Smith if (ilink->fields) { 900971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 9179416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 925a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9379416396SBarry Smith if (j > 0) { 9479416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9579416396SBarry Smith } 965a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 970971522cSBarry Smith } 980971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1001ab39975SBarry Smith } else { 1011ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1021ab39975SBarry Smith } 1035a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1045a9f2f41SSatish Balay ilink = ilink->next; 1050971522cSBarry Smith } 1060971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1070971522cSBarry Smith } else { 10865e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1090971522cSBarry Smith } 1100971522cSBarry Smith PetscFunctionReturn(0); 1110971522cSBarry Smith } 1120971522cSBarry Smith 1130971522cSBarry Smith #undef __FUNCT__ 1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1163b224e63SBarry Smith { 1173b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1183b224e63SBarry Smith PetscErrorCode ierr; 119ace3abfcSBarry Smith PetscBool iascii; 1203b224e63SBarry Smith PetscInt i,j; 1213b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1223b224e63SBarry Smith KSP ksp; 1233b224e63SBarry Smith 1243b224e63SBarry Smith PetscFunctionBegin; 1252692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1263b224e63SBarry Smith if (iascii) { 127c5d2311dSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr); 1283b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1293b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1303b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1313b224e63SBarry Smith if (ilink->fields) { 1323b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1333b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1343b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1353b224e63SBarry Smith if (j > 0) { 1363b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1373b224e63SBarry Smith } 1383b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1393b224e63SBarry Smith } 1403b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1413b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1423b224e63SBarry Smith } else { 1433b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1443b224e63SBarry Smith } 1453b224e63SBarry Smith ilink = ilink->next; 1463b224e63SBarry Smith } 147435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1483b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14912cae6f2SJed Brown if (jac->schur) { 15012cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1513b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 15212cae6f2SJed Brown } else { 15312cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 15412cae6f2SJed Brown } 1553b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 156435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 1573b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15812cae6f2SJed Brown if (jac->kspschur) { 1593b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 16012cae6f2SJed Brown } else { 16112cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 16212cae6f2SJed Brown } 1633b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1643b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1653b224e63SBarry Smith } else { 16665e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1673b224e63SBarry Smith } 1683b224e63SBarry Smith PetscFunctionReturn(0); 1693b224e63SBarry Smith } 1703b224e63SBarry Smith 1713b224e63SBarry Smith #undef __FUNCT__ 1726c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1736c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1746c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1756c924f48SJed Brown { 1766c924f48SJed Brown PetscErrorCode ierr; 1776c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1786c924f48SJed Brown PetscInt i,nfields,*ifields; 179ace3abfcSBarry Smith PetscBool flg; 1806c924f48SJed Brown char optionname[128],splitname[8]; 1816c924f48SJed Brown 1826c924f48SJed Brown PetscFunctionBegin; 1836c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 1846c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 1856c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 1866c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 1876c924f48SJed Brown nfields = jac->bs; 1886c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 1896c924f48SJed Brown if (!flg) break; 1906c924f48SJed Brown if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1916c924f48SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr); 1926c924f48SJed Brown } 1936c924f48SJed Brown if (i > 0) { 1946c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 1956c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 1966c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 1976c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 1986c924f48SJed Brown } 1996c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2006c924f48SJed Brown PetscFunctionReturn(0); 2016c924f48SJed Brown } 2026c924f48SJed Brown 2036c924f48SJed Brown #undef __FUNCT__ 20469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 20569a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2060971522cSBarry Smith { 2070971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2080971522cSBarry Smith PetscErrorCode ierr; 2095a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2106ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2116c924f48SJed Brown PetscInt i; 2120971522cSBarry Smith 2130971522cSBarry Smith PetscFunctionBegin; 214d32f9abdSBarry Smith if (!ilink) { 215d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 216d0af7cd3SBarry Smith if (pc->dm && !stokes) { 2178b8307b2SJed Brown PetscBool dmcomposite; 2188b8307b2SJed Brown ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 2198b8307b2SJed Brown if (dmcomposite) { 2208b8307b2SJed Brown PetscInt nDM; 2218b8307b2SJed Brown IS *fields; 2228b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2238b8307b2SJed Brown ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 2248b8307b2SJed Brown ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 2258b8307b2SJed Brown for (i=0; i<nDM; i++) { 2268b8307b2SJed Brown char splitname[8]; 2278b8307b2SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2288b8307b2SJed Brown ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr); 229fcfd50ebSBarry Smith ierr = ISDestroy(&fields[i]);CHKERRQ(ierr); 2308b8307b2SJed Brown } 2318b8307b2SJed Brown ierr = PetscFree(fields);CHKERRQ(ierr); 2328b8307b2SJed Brown } 23366ffff09SJed Brown } else { 234521d7252SBarry Smith if (jac->bs <= 0) { 235704ba839SBarry Smith if (pc->pmat) { 236521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 237704ba839SBarry Smith } else { 238704ba839SBarry Smith jac->bs = 1; 239704ba839SBarry Smith } 240521d7252SBarry Smith } 241d32f9abdSBarry Smith 242acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2436ce1633cSBarry Smith if (stokes) { 2446ce1633cSBarry Smith IS zerodiags,rest; 2456ce1633cSBarry Smith PetscInt nmin,nmax; 2466ce1633cSBarry Smith 2476ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2486ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2496ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2506ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2516ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 252fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 253fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2546ce1633cSBarry Smith } else { 255d32f9abdSBarry Smith if (!flg) { 256d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 257d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2586c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2596c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 260d32f9abdSBarry Smith } 2616c924f48SJed Brown if (flg || !jac->splitdefined) { 262d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 263db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2646c924f48SJed Brown char splitname[8]; 2656c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 266db4c96c1SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr); 26779416396SBarry Smith } 26897bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 269521d7252SBarry Smith } 27066ffff09SJed Brown } 2716ce1633cSBarry Smith } 272edf189efSBarry Smith } else if (jac->nsplits == 1) { 273edf189efSBarry Smith if (ilink->is) { 274edf189efSBarry Smith IS is2; 275edf189efSBarry Smith PetscInt nmin,nmax; 276edf189efSBarry Smith 277edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 278edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 279db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 280fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 281e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 28263ec74ffSBarry Smith } else if (jac->reset) { 28363ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 284d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 285d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 286d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 287d0af7cd3SBarry Smith if (pc->dm && !stokes) { 288d0af7cd3SBarry Smith PetscBool dmcomposite; 289d0af7cd3SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 290d0af7cd3SBarry Smith if (dmcomposite) { 291d0af7cd3SBarry Smith PetscInt nDM; 292d0af7cd3SBarry Smith IS *fields; 293d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 294d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 295d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 296d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 297d0af7cd3SBarry Smith ilink->is = fields[i]; 298d0af7cd3SBarry Smith ilink = ilink->next; 299edf189efSBarry Smith } 300d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 301d0af7cd3SBarry Smith } 302d0af7cd3SBarry Smith } else { 303d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 304d0af7cd3SBarry Smith if (stokes) { 305d0af7cd3SBarry Smith IS zerodiags,rest; 306d0af7cd3SBarry Smith PetscInt nmin,nmax; 307d0af7cd3SBarry Smith 308d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 309d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 310d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 311d0af7cd3SBarry Smith ilink->is = rest; 312d0af7cd3SBarry Smith ilink->next->is = zerodiags; 313d0af7cd3SBarry Smith } else { 314d0af7cd3SBarry Smith SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 315d0af7cd3SBarry Smith } 316d0af7cd3SBarry Smith } 317d0af7cd3SBarry Smith } 318d0af7cd3SBarry Smith 319e7e72b3dSBarry Smith if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 32069a612a9SBarry Smith PetscFunctionReturn(0); 32169a612a9SBarry Smith } 32269a612a9SBarry Smith 32369a612a9SBarry Smith #undef __FUNCT__ 32469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 32569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 32669a612a9SBarry Smith { 32769a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 32869a612a9SBarry Smith PetscErrorCode ierr; 3295a9f2f41SSatish Balay PC_FieldSplitLink ilink; 330cf502942SBarry Smith PetscInt i,nsplit,ccsize; 33169a612a9SBarry Smith MatStructure flag = pc->flag; 332ace3abfcSBarry Smith PetscBool sorted; 33369a612a9SBarry Smith 33469a612a9SBarry Smith PetscFunctionBegin; 33569a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 33697bbdb24SBarry Smith nsplit = jac->nsplits; 3375a9f2f41SSatish Balay ilink = jac->head; 33897bbdb24SBarry Smith 33997bbdb24SBarry Smith /* get the matrices for each split */ 340704ba839SBarry Smith if (!jac->issetup) { 3411b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 34297bbdb24SBarry Smith 343704ba839SBarry Smith jac->issetup = PETSC_TRUE; 344704ba839SBarry Smith 345704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 34651f519a2SBarry Smith bs = jac->bs; 34797bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 348cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 3491b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3501b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3511b9fc7fcSBarry Smith if (jac->defaultsplit) { 352704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 353704ba839SBarry Smith } else if (!ilink->is) { 354ccb205f8SBarry Smith if (ilink->nfields > 1) { 3555a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 3565a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3571b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3581b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3591b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 36097bbdb24SBarry Smith } 36197bbdb24SBarry Smith } 362d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 363ccb205f8SBarry Smith } else { 364704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 365ccb205f8SBarry Smith } 3663e197d65SBarry Smith } 367704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 368e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 369704ba839SBarry Smith ilink = ilink->next; 3701b9fc7fcSBarry Smith } 3711b9fc7fcSBarry Smith } 3721b9fc7fcSBarry Smith 373704ba839SBarry Smith ilink = jac->head; 37497bbdb24SBarry Smith if (!jac->pmat) { 375cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 376cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3774aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 378704ba839SBarry Smith ilink = ilink->next; 379cf502942SBarry Smith } 38097bbdb24SBarry Smith } else { 381cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3824aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 383704ba839SBarry Smith ilink = ilink->next; 384cf502942SBarry Smith } 38597bbdb24SBarry Smith } 386519d70e2SJed Brown if (jac->realdiagonal) { 387519d70e2SJed Brown ilink = jac->head; 388519d70e2SJed Brown if (!jac->mat) { 389519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 390519d70e2SJed Brown for (i=0; i<nsplit; i++) { 391519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 392519d70e2SJed Brown ilink = ilink->next; 393519d70e2SJed Brown } 394519d70e2SJed Brown } else { 395519d70e2SJed Brown for (i=0; i<nsplit; i++) { 396966e74d7SJed Brown if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 397519d70e2SJed Brown ilink = ilink->next; 398519d70e2SJed Brown } 399519d70e2SJed Brown } 400519d70e2SJed Brown } else { 401519d70e2SJed Brown jac->mat = jac->pmat; 402519d70e2SJed Brown } 40397bbdb24SBarry Smith 4046c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 40568dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 40668dd23aaSBarry Smith ilink = jac->head; 40768dd23aaSBarry Smith if (!jac->Afield) { 40868dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 40968dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4104aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 41168dd23aaSBarry Smith ilink = ilink->next; 41268dd23aaSBarry Smith } 41368dd23aaSBarry Smith } else { 41468dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4154aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 41668dd23aaSBarry Smith ilink = ilink->next; 41768dd23aaSBarry Smith } 41868dd23aaSBarry Smith } 41968dd23aaSBarry Smith } 42068dd23aaSBarry Smith 4213b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 4223b224e63SBarry Smith IS ccis; 4234aa3045dSJed Brown PetscInt rstart,rend; 424e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 42568dd23aaSBarry Smith 426e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 427e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 428e6cab6aaSJed Brown 4293b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4303b224e63SBarry Smith if (jac->schur) { 4313b224e63SBarry Smith ilink = jac->head; 4324aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4334aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 434fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4353b224e63SBarry Smith ilink = ilink->next; 4364aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4374aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 438fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 439519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 440084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4413b224e63SBarry Smith 4423b224e63SBarry Smith } else { 4431cee3971SBarry Smith KSP ksp; 4446c924f48SJed Brown char schurprefix[256]; 4453b224e63SBarry Smith 446a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4473b224e63SBarry Smith ilink = jac->head; 4484aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4494aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 450fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4513b224e63SBarry Smith ilink = ilink->next; 4524aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4534aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 454fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4557d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 456519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 457a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4581cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 459e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 460a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 461a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 462*c1570756SJed Brown if (!jac->suboptionsset) {ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);} 4633b224e63SBarry Smith 4643b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4659005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4661cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 467084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 468084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 469084e4875SJed Brown PC pc; 470cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 471084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 472084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 473e69d4d44SBarry Smith } 4746c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 4756c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 4763b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 477*c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);} 4783b224e63SBarry Smith 4793b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 4803b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 4813b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 4823b224e63SBarry Smith ilink = jac->head; 4833b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 4843b224e63SBarry Smith ilink = ilink->next; 4853b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 4863b224e63SBarry Smith } 4873b224e63SBarry Smith } else { 48897bbdb24SBarry Smith /* set up the individual PCs */ 48997bbdb24SBarry Smith i = 0; 4905a9f2f41SSatish Balay ilink = jac->head; 4915a9f2f41SSatish Balay while (ilink) { 492519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 4933b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 494*c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 4955a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 49697bbdb24SBarry Smith i++; 4975a9f2f41SSatish Balay ilink = ilink->next; 4980971522cSBarry Smith } 49997bbdb24SBarry Smith 50097bbdb24SBarry Smith /* create work vectors for each split */ 5011b9fc7fcSBarry Smith if (!jac->x) { 50297bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5035a9f2f41SSatish Balay ilink = jac->head; 50497bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 505906ed7ccSBarry Smith Vec *vl,*vr; 506906ed7ccSBarry Smith 507906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 508906ed7ccSBarry Smith ilink->x = *vr; 509906ed7ccSBarry Smith ilink->y = *vl; 510906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 511906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5125a9f2f41SSatish Balay jac->x[i] = ilink->x; 5135a9f2f41SSatish Balay jac->y[i] = ilink->y; 5145a9f2f41SSatish Balay ilink = ilink->next; 51597bbdb24SBarry Smith } 5163b224e63SBarry Smith } 5173b224e63SBarry Smith } 5183b224e63SBarry Smith 5193b224e63SBarry Smith 520d0f46423SBarry Smith if (!jac->head->sctx) { 5213b224e63SBarry Smith Vec xtmp; 5223b224e63SBarry Smith 52379416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5241b9fc7fcSBarry Smith 5255a9f2f41SSatish Balay ilink = jac->head; 5261b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5271b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 528704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5295a9f2f41SSatish Balay ilink = ilink->next; 5301b9fc7fcSBarry Smith } 5316bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5321b9fc7fcSBarry Smith } 533*c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 5340971522cSBarry Smith PetscFunctionReturn(0); 5350971522cSBarry Smith } 5360971522cSBarry Smith 5375a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 538ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 539ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5405a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 541ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 542ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 54379416396SBarry Smith 5440971522cSBarry Smith #undef __FUNCT__ 5453b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5463b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5473b224e63SBarry Smith { 5483b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5493b224e63SBarry Smith PetscErrorCode ierr; 5503b224e63SBarry Smith KSP ksp; 5513b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5523b224e63SBarry Smith 5533b224e63SBarry Smith PetscFunctionBegin; 5543b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5553b224e63SBarry Smith 556c5d2311dSJed Brown switch (jac->schurfactorization) { 557c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 558a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 559c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 560c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 561c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 562c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 563c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 564c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 565c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 566c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 567c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 568c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 569c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 570c5d2311dSJed Brown break; 571c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 572a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 573c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 574c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 575c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 576c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 577c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 578c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 580c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 581c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 582c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 583c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 584c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 585c5d2311dSJed Brown break; 586c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 587a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 588c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 589c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 590c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 591c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 592c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 593c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 594c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 595c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 596c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 598c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 599c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 600c5d2311dSJed Brown break; 601c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 602a04f6461SBarry 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 */ 6033b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6043b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6053b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6063b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6073b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6083b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6093b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6103b224e63SBarry Smith 6113b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6123b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6133b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6143b224e63SBarry Smith 6153b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6163b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6173b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6183b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6193b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 620c5d2311dSJed Brown } 6213b224e63SBarry Smith PetscFunctionReturn(0); 6223b224e63SBarry Smith } 6233b224e63SBarry Smith 6243b224e63SBarry Smith #undef __FUNCT__ 6250971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6260971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6270971522cSBarry Smith { 6280971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6290971522cSBarry Smith PetscErrorCode ierr; 6305a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 631af93645bSJed Brown PetscInt cnt; 6320971522cSBarry Smith 6330971522cSBarry Smith PetscFunctionBegin; 63451f519a2SBarry Smith CHKMEMQ; 63551f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 63651f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 63751f519a2SBarry Smith 63879416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6391b9fc7fcSBarry Smith if (jac->defaultsplit) { 6400971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6415a9f2f41SSatish Balay while (ilink) { 6425a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6435a9f2f41SSatish Balay ilink = ilink->next; 6440971522cSBarry Smith } 6450971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6461b9fc7fcSBarry Smith } else { 647efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6485a9f2f41SSatish Balay while (ilink) { 6495a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6505a9f2f41SSatish Balay ilink = ilink->next; 6511b9fc7fcSBarry Smith } 6521b9fc7fcSBarry Smith } 65316913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 65479416396SBarry Smith if (!jac->w1) { 65579416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 65679416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 65779416396SBarry Smith } 658efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6595a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6603e197d65SBarry Smith cnt = 1; 6615a9f2f41SSatish Balay while (ilink->next) { 6625a9f2f41SSatish Balay ilink = ilink->next; 6633e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6643e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6653e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6663e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6673e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6683e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6693e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6703e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6713e197d65SBarry Smith } 67251f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 67311755939SBarry Smith cnt -= 2; 67451f519a2SBarry Smith while (ilink->previous) { 67551f519a2SBarry Smith ilink = ilink->previous; 67611755939SBarry Smith /* compute the residual only over the part of the vector needed */ 67711755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 67811755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 67911755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 68011755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 68111755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 68211755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 68311755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 68451f519a2SBarry Smith } 68511755939SBarry Smith } 68665e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 68751f519a2SBarry Smith CHKMEMQ; 6880971522cSBarry Smith PetscFunctionReturn(0); 6890971522cSBarry Smith } 6900971522cSBarry Smith 691421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 692ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 693ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 694421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 695ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 696ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 697421e10b8SBarry Smith 698421e10b8SBarry Smith #undef __FUNCT__ 6998c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 700421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 701421e10b8SBarry Smith { 702421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 703421e10b8SBarry Smith PetscErrorCode ierr; 704421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 705421e10b8SBarry Smith 706421e10b8SBarry Smith PetscFunctionBegin; 707421e10b8SBarry Smith CHKMEMQ; 708421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 709421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 710421e10b8SBarry Smith 711421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 712421e10b8SBarry Smith if (jac->defaultsplit) { 713421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 714421e10b8SBarry Smith while (ilink) { 715421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 716421e10b8SBarry Smith ilink = ilink->next; 717421e10b8SBarry Smith } 718421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 719421e10b8SBarry Smith } else { 720421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 721421e10b8SBarry Smith while (ilink) { 722421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 723421e10b8SBarry Smith ilink = ilink->next; 724421e10b8SBarry Smith } 725421e10b8SBarry Smith } 726421e10b8SBarry Smith } else { 727421e10b8SBarry Smith if (!jac->w1) { 728421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 729421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 730421e10b8SBarry Smith } 731421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 732421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 733421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 734421e10b8SBarry Smith while (ilink->next) { 735421e10b8SBarry Smith ilink = ilink->next; 7369989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 737421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 738421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 739421e10b8SBarry Smith } 740421e10b8SBarry Smith while (ilink->previous) { 741421e10b8SBarry Smith ilink = ilink->previous; 7429989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 743421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 744421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 745421e10b8SBarry Smith } 746421e10b8SBarry Smith } else { 747421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 748421e10b8SBarry Smith ilink = ilink->next; 749421e10b8SBarry Smith } 750421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 751421e10b8SBarry Smith while (ilink->previous) { 752421e10b8SBarry Smith ilink = ilink->previous; 7539989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 754421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 755421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 756421e10b8SBarry Smith } 757421e10b8SBarry Smith } 758421e10b8SBarry Smith } 759421e10b8SBarry Smith CHKMEMQ; 760421e10b8SBarry Smith PetscFunctionReturn(0); 761421e10b8SBarry Smith } 762421e10b8SBarry Smith 7630971522cSBarry Smith #undef __FUNCT__ 764574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 765574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 7660971522cSBarry Smith { 7670971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7680971522cSBarry Smith PetscErrorCode ierr; 7695a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 7700971522cSBarry Smith 7710971522cSBarry Smith PetscFunctionBegin; 7725a9f2f41SSatish Balay while (ilink) { 773574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 774fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 775fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 776fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 777fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 7785a9f2f41SSatish Balay next = ilink->next; 7795a9f2f41SSatish Balay ilink = next; 7800971522cSBarry Smith } 78105b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 782574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 783574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 784574deadeSBarry Smith } else if (jac->mat) { 785574deadeSBarry Smith jac->mat = PETSC_NULL; 786574deadeSBarry Smith } 78797bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 78868dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 7896bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 7906bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 7916bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 7926bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 7936bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 7946bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 7956bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 79663ec74ffSBarry Smith jac->reset = PETSC_TRUE; 797574deadeSBarry Smith PetscFunctionReturn(0); 798574deadeSBarry Smith } 799574deadeSBarry Smith 800574deadeSBarry Smith #undef __FUNCT__ 801574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 802574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 803574deadeSBarry Smith { 804574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 805574deadeSBarry Smith PetscErrorCode ierr; 806574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 807574deadeSBarry Smith 808574deadeSBarry Smith PetscFunctionBegin; 809574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 810574deadeSBarry Smith while (ilink) { 8116bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 812574deadeSBarry Smith next = ilink->next; 813574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 814574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 815574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 816574deadeSBarry Smith ilink = next; 817574deadeSBarry Smith } 818574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 819c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8200971522cSBarry Smith PetscFunctionReturn(0); 8210971522cSBarry Smith } 8220971522cSBarry Smith 8230971522cSBarry Smith #undef __FUNCT__ 8240971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8250971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8260971522cSBarry Smith { 8271b9fc7fcSBarry Smith PetscErrorCode ierr; 8286c924f48SJed Brown PetscInt bs; 829bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8309dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8313b224e63SBarry Smith PCCompositeType ctype; 8321b9fc7fcSBarry Smith 8330971522cSBarry Smith PetscFunctionBegin; 8341b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 835acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 83651f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 83751f519a2SBarry Smith if (flg) { 83851f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 83951f519a2SBarry Smith } 840704ba839SBarry Smith 841435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 842c0adfefeSBarry Smith if (stokes) { 843c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 844c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 845c0adfefeSBarry Smith } 846c0adfefeSBarry Smith 8473b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 8483b224e63SBarry Smith if (flg) { 8493b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 8503b224e63SBarry Smith } 851d32f9abdSBarry Smith 852c30613efSMatthew Knepley /* Only setup fields once */ 853c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 854d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 855d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 8566c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 8576c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 858d32f9abdSBarry Smith } 859c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 860c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 861084e4875SJed 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); 862c5d2311dSJed Brown } 8631b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8640971522cSBarry Smith PetscFunctionReturn(0); 8650971522cSBarry Smith } 8660971522cSBarry Smith 8670971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 8680971522cSBarry Smith 8690971522cSBarry Smith EXTERN_C_BEGIN 8700971522cSBarry Smith #undef __FUNCT__ 8710971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 8727087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 8730971522cSBarry Smith { 87497bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8750971522cSBarry Smith PetscErrorCode ierr; 8765a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 87769a612a9SBarry Smith char prefix[128]; 87851f519a2SBarry Smith PetscInt i; 8790971522cSBarry Smith 8800971522cSBarry Smith PetscFunctionBegin; 8816c924f48SJed Brown if (jac->splitdefined) { 8826c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 8836c924f48SJed Brown PetscFunctionReturn(0); 8846c924f48SJed Brown } 88551f519a2SBarry Smith for (i=0; i<n; i++) { 886e32f2f54SBarry 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); 887e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 88851f519a2SBarry Smith } 889704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 890a04f6461SBarry Smith if (splitname) { 891db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 892a04f6461SBarry Smith } else { 893a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 894a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 895a04f6461SBarry Smith } 896704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 8975a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 8985a9f2f41SSatish Balay ilink->nfields = n; 8995a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9007adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9011cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9025a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9039005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 90469a612a9SBarry Smith 905a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9065a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9070971522cSBarry Smith 9080971522cSBarry Smith if (!next) { 9095a9f2f41SSatish Balay jac->head = ilink; 91051f519a2SBarry Smith ilink->previous = PETSC_NULL; 9110971522cSBarry Smith } else { 9120971522cSBarry Smith while (next->next) { 9130971522cSBarry Smith next = next->next; 9140971522cSBarry Smith } 9155a9f2f41SSatish Balay next->next = ilink; 91651f519a2SBarry Smith ilink->previous = next; 9170971522cSBarry Smith } 9180971522cSBarry Smith jac->nsplits++; 9190971522cSBarry Smith PetscFunctionReturn(0); 9200971522cSBarry Smith } 9210971522cSBarry Smith EXTERN_C_END 9220971522cSBarry Smith 923e69d4d44SBarry Smith EXTERN_C_BEGIN 924e69d4d44SBarry Smith #undef __FUNCT__ 925e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9267087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 927e69d4d44SBarry Smith { 928e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 929e69d4d44SBarry Smith PetscErrorCode ierr; 930e69d4d44SBarry Smith 931e69d4d44SBarry Smith PetscFunctionBegin; 932519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 933e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 934e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 93513e0d083SBarry Smith if (n) *n = jac->nsplits; 936e69d4d44SBarry Smith PetscFunctionReturn(0); 937e69d4d44SBarry Smith } 938e69d4d44SBarry Smith EXTERN_C_END 9390971522cSBarry Smith 9400971522cSBarry Smith EXTERN_C_BEGIN 9410971522cSBarry Smith #undef __FUNCT__ 94269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 9437087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 9440971522cSBarry Smith { 9450971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9460971522cSBarry Smith PetscErrorCode ierr; 9470971522cSBarry Smith PetscInt cnt = 0; 9485a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 9490971522cSBarry Smith 9500971522cSBarry Smith PetscFunctionBegin; 9515d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 9525a9f2f41SSatish Balay while (ilink) { 9535a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 9545a9f2f41SSatish Balay ilink = ilink->next; 9550971522cSBarry Smith } 9565d480477SMatthew 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); 95713e0d083SBarry Smith if (n) *n = jac->nsplits; 9580971522cSBarry Smith PetscFunctionReturn(0); 9590971522cSBarry Smith } 9600971522cSBarry Smith EXTERN_C_END 9610971522cSBarry Smith 962704ba839SBarry Smith EXTERN_C_BEGIN 963704ba839SBarry Smith #undef __FUNCT__ 964704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 9657087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 966704ba839SBarry Smith { 967704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 968704ba839SBarry Smith PetscErrorCode ierr; 969704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 970704ba839SBarry Smith char prefix[128]; 971704ba839SBarry Smith 972704ba839SBarry Smith PetscFunctionBegin; 9736c924f48SJed Brown if (jac->splitdefined) { 9746c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9756c924f48SJed Brown PetscFunctionReturn(0); 9766c924f48SJed Brown } 97716913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 978a04f6461SBarry Smith if (splitname) { 979db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 980a04f6461SBarry Smith } else { 981a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 982a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 983a04f6461SBarry Smith } 9841ab39975SBarry Smith ilink->is = is; 9851ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 986704ba839SBarry Smith ilink->next = PETSC_NULL; 987704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9881cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 989704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9909005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 991704ba839SBarry Smith 992a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 993704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 994704ba839SBarry Smith 995704ba839SBarry Smith if (!next) { 996704ba839SBarry Smith jac->head = ilink; 997704ba839SBarry Smith ilink->previous = PETSC_NULL; 998704ba839SBarry Smith } else { 999704ba839SBarry Smith while (next->next) { 1000704ba839SBarry Smith next = next->next; 1001704ba839SBarry Smith } 1002704ba839SBarry Smith next->next = ilink; 1003704ba839SBarry Smith ilink->previous = next; 1004704ba839SBarry Smith } 1005704ba839SBarry Smith jac->nsplits++; 1006704ba839SBarry Smith 1007704ba839SBarry Smith PetscFunctionReturn(0); 1008704ba839SBarry Smith } 1009704ba839SBarry Smith EXTERN_C_END 1010704ba839SBarry Smith 10110971522cSBarry Smith #undef __FUNCT__ 10120971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10130971522cSBarry Smith /*@ 10140971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10150971522cSBarry Smith 1016ad4df100SBarry Smith Logically Collective on PC 10170971522cSBarry Smith 10180971522cSBarry Smith Input Parameters: 10190971522cSBarry Smith + pc - the preconditioner context 1020a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10210971522cSBarry Smith . n - the number of fields in this split 1022db4c96c1SJed Brown - fields - the fields in this split 10230971522cSBarry Smith 10240971522cSBarry Smith Level: intermediate 10250971522cSBarry Smith 1026d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1027d32f9abdSBarry Smith 1028d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1029d32f9abdSBarry 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 1030d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1031d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1032d32f9abdSBarry Smith 1033db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1034db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1035db4c96c1SJed Brown 1036d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 10370971522cSBarry Smith 10380971522cSBarry Smith @*/ 10397087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 10400971522cSBarry Smith { 10414ac538c5SBarry Smith PetscErrorCode ierr; 10420971522cSBarry Smith 10430971522cSBarry Smith PetscFunctionBegin; 10440700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1045db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1046db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1047db4c96c1SJed Brown PetscValidIntPointer(fields,3); 10484ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr); 10490971522cSBarry Smith PetscFunctionReturn(0); 10500971522cSBarry Smith } 10510971522cSBarry Smith 10520971522cSBarry Smith #undef __FUNCT__ 1053704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1054704ba839SBarry Smith /*@ 1055704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1056704ba839SBarry Smith 1057ad4df100SBarry Smith Logically Collective on PC 1058704ba839SBarry Smith 1059704ba839SBarry Smith Input Parameters: 1060704ba839SBarry Smith + pc - the preconditioner context 1061a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1062db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1063704ba839SBarry Smith 1064d32f9abdSBarry Smith 1065a6ffb8dbSJed Brown Notes: 1066a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1067a6ffb8dbSJed Brown 1068db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1069db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1070d32f9abdSBarry Smith 1071704ba839SBarry Smith Level: intermediate 1072704ba839SBarry Smith 1073704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1074704ba839SBarry Smith 1075704ba839SBarry Smith @*/ 10767087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1077704ba839SBarry Smith { 10784ac538c5SBarry Smith PetscErrorCode ierr; 1079704ba839SBarry Smith 1080704ba839SBarry Smith PetscFunctionBegin; 10810700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1082db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1083db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 10844ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1085704ba839SBarry Smith PetscFunctionReturn(0); 1086704ba839SBarry Smith } 1087704ba839SBarry Smith 1088704ba839SBarry Smith #undef __FUNCT__ 108957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 109057a9adfeSMatthew G Knepley /*@ 109157a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 109257a9adfeSMatthew G Knepley 109357a9adfeSMatthew G Knepley Logically Collective on PC 109457a9adfeSMatthew G Knepley 109557a9adfeSMatthew G Knepley Input Parameters: 109657a9adfeSMatthew G Knepley + pc - the preconditioner context 109757a9adfeSMatthew G Knepley - splitname - name of this split 109857a9adfeSMatthew G Knepley 109957a9adfeSMatthew G Knepley Output Parameter: 110057a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 110157a9adfeSMatthew G Knepley 110257a9adfeSMatthew G Knepley Level: intermediate 110357a9adfeSMatthew G Knepley 110457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 110557a9adfeSMatthew G Knepley 110657a9adfeSMatthew G Knepley @*/ 110757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 110857a9adfeSMatthew G Knepley { 110957a9adfeSMatthew G Knepley PetscErrorCode ierr; 111057a9adfeSMatthew G Knepley 111157a9adfeSMatthew G Knepley PetscFunctionBegin; 111257a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 111357a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 111457a9adfeSMatthew G Knepley PetscValidPointer(is,3); 111557a9adfeSMatthew G Knepley { 111657a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 111757a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 111857a9adfeSMatthew G Knepley PetscBool found; 111957a9adfeSMatthew G Knepley 112057a9adfeSMatthew G Knepley *is = PETSC_NULL; 112157a9adfeSMatthew G Knepley while(ilink) { 112257a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 112357a9adfeSMatthew G Knepley if (found) { 112457a9adfeSMatthew G Knepley *is = ilink->is; 112557a9adfeSMatthew G Knepley break; 112657a9adfeSMatthew G Knepley } 112757a9adfeSMatthew G Knepley ilink = ilink->next; 112857a9adfeSMatthew G Knepley } 112957a9adfeSMatthew G Knepley } 113057a9adfeSMatthew G Knepley PetscFunctionReturn(0); 113157a9adfeSMatthew G Knepley } 113257a9adfeSMatthew G Knepley 113357a9adfeSMatthew G Knepley #undef __FUNCT__ 113451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 113551f519a2SBarry Smith /*@ 113651f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 113751f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 113851f519a2SBarry Smith 1139ad4df100SBarry Smith Logically Collective on PC 114051f519a2SBarry Smith 114151f519a2SBarry Smith Input Parameters: 114251f519a2SBarry Smith + pc - the preconditioner context 114351f519a2SBarry Smith - bs - the block size 114451f519a2SBarry Smith 114551f519a2SBarry Smith Level: intermediate 114651f519a2SBarry Smith 114751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 114851f519a2SBarry Smith 114951f519a2SBarry Smith @*/ 11507087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 115151f519a2SBarry Smith { 11524ac538c5SBarry Smith PetscErrorCode ierr; 115351f519a2SBarry Smith 115451f519a2SBarry Smith PetscFunctionBegin; 11550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1156c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 11574ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 115851f519a2SBarry Smith PetscFunctionReturn(0); 115951f519a2SBarry Smith } 116051f519a2SBarry Smith 116151f519a2SBarry Smith #undef __FUNCT__ 116269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 11630971522cSBarry Smith /*@C 116469a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 11650971522cSBarry Smith 116669a612a9SBarry Smith Collective on KSP 11670971522cSBarry Smith 11680971522cSBarry Smith Input Parameter: 11690971522cSBarry Smith . pc - the preconditioner context 11700971522cSBarry Smith 11710971522cSBarry Smith Output Parameters: 117213e0d083SBarry Smith + n - the number of splits 117369a612a9SBarry Smith - pc - the array of KSP contexts 11740971522cSBarry Smith 11750971522cSBarry Smith Note: 1176d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1177d32f9abdSBarry Smith (not the KSP just the array that contains them). 11780971522cSBarry Smith 117969a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 11800971522cSBarry Smith 11810971522cSBarry Smith Level: advanced 11820971522cSBarry Smith 11830971522cSBarry Smith .seealso: PCFIELDSPLIT 11840971522cSBarry Smith @*/ 11857087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 11860971522cSBarry Smith { 11874ac538c5SBarry Smith PetscErrorCode ierr; 11880971522cSBarry Smith 11890971522cSBarry Smith PetscFunctionBegin; 11900700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 119113e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 11924ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 11930971522cSBarry Smith PetscFunctionReturn(0); 11940971522cSBarry Smith } 11950971522cSBarry Smith 1196e69d4d44SBarry Smith #undef __FUNCT__ 1197e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1198e69d4d44SBarry Smith /*@ 1199e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1200a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1201e69d4d44SBarry Smith 1202e69d4d44SBarry Smith Collective on PC 1203e69d4d44SBarry Smith 1204e69d4d44SBarry Smith Input Parameters: 1205e69d4d44SBarry Smith + pc - the preconditioner context 1206fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1207084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1208084e4875SJed Brown 1209e69d4d44SBarry Smith Options Database: 1210084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1211e69d4d44SBarry Smith 1212fd1303e9SJungho Lee Notes: 1213fd1303e9SJungho Lee $ If ptype is 1214fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1215fd1303e9SJungho Lee $ to this function). 1216fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1217fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1218fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1219fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1220fd1303e9SJungho Lee $ preconditioner 1221fd1303e9SJungho Lee 1222fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1223fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1224fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1225fd1303e9SJungho Lee 1226fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1227fd1303e9SJungho Lee the name since it sets a proceedure to use. 1228fd1303e9SJungho Lee 1229e69d4d44SBarry Smith Level: intermediate 1230e69d4d44SBarry Smith 1231fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1232e69d4d44SBarry Smith 1233e69d4d44SBarry Smith @*/ 12347087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1235e69d4d44SBarry Smith { 12364ac538c5SBarry Smith PetscErrorCode ierr; 1237e69d4d44SBarry Smith 1238e69d4d44SBarry Smith PetscFunctionBegin; 12390700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12404ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1241e69d4d44SBarry Smith PetscFunctionReturn(0); 1242e69d4d44SBarry Smith } 1243e69d4d44SBarry Smith 1244e69d4d44SBarry Smith EXTERN_C_BEGIN 1245e69d4d44SBarry Smith #undef __FUNCT__ 1246e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 12477087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1248e69d4d44SBarry Smith { 1249e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1250084e4875SJed Brown PetscErrorCode ierr; 1251e69d4d44SBarry Smith 1252e69d4d44SBarry Smith PetscFunctionBegin; 1253084e4875SJed Brown jac->schurpre = ptype; 1254084e4875SJed Brown if (pre) { 12556bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1256084e4875SJed Brown jac->schur_user = pre; 1257084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1258084e4875SJed Brown } 1259e69d4d44SBarry Smith PetscFunctionReturn(0); 1260e69d4d44SBarry Smith } 1261e69d4d44SBarry Smith EXTERN_C_END 1262e69d4d44SBarry Smith 126330ad9308SMatthew Knepley #undef __FUNCT__ 126430ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 126530ad9308SMatthew Knepley /*@C 12668c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 126730ad9308SMatthew Knepley 126830ad9308SMatthew Knepley Collective on KSP 126930ad9308SMatthew Knepley 127030ad9308SMatthew Knepley Input Parameter: 127130ad9308SMatthew Knepley . pc - the preconditioner context 127230ad9308SMatthew Knepley 127330ad9308SMatthew Knepley Output Parameters: 1274a04f6461SBarry Smith + A00 - the (0,0) block 1275a04f6461SBarry Smith . A01 - the (0,1) block 1276a04f6461SBarry Smith . A10 - the (1,0) block 1277a04f6461SBarry Smith - A11 - the (1,1) block 127830ad9308SMatthew Knepley 127930ad9308SMatthew Knepley Level: advanced 128030ad9308SMatthew Knepley 128130ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 128230ad9308SMatthew Knepley @*/ 1283a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 128430ad9308SMatthew Knepley { 128530ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 128630ad9308SMatthew Knepley 128730ad9308SMatthew Knepley PetscFunctionBegin; 12880700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1289c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1290a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1291a04f6461SBarry Smith if (A01) *A01 = jac->B; 1292a04f6461SBarry Smith if (A10) *A10 = jac->C; 1293a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 129430ad9308SMatthew Knepley PetscFunctionReturn(0); 129530ad9308SMatthew Knepley } 129630ad9308SMatthew Knepley 129779416396SBarry Smith EXTERN_C_BEGIN 129879416396SBarry Smith #undef __FUNCT__ 129979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 13007087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 130179416396SBarry Smith { 130279416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1303e69d4d44SBarry Smith PetscErrorCode ierr; 130479416396SBarry Smith 130579416396SBarry Smith PetscFunctionBegin; 130679416396SBarry Smith jac->type = type; 13073b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 13083b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 13093b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1310e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1311e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1312e69d4d44SBarry Smith 13133b224e63SBarry Smith } else { 13143b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 13153b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1316e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13179e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 13183b224e63SBarry Smith } 131979416396SBarry Smith PetscFunctionReturn(0); 132079416396SBarry Smith } 132179416396SBarry Smith EXTERN_C_END 132279416396SBarry Smith 132351f519a2SBarry Smith EXTERN_C_BEGIN 132451f519a2SBarry Smith #undef __FUNCT__ 132551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 13267087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 132751f519a2SBarry Smith { 132851f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 132951f519a2SBarry Smith 133051f519a2SBarry Smith PetscFunctionBegin; 133165e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 133265e19b50SBarry 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); 133351f519a2SBarry Smith jac->bs = bs; 133451f519a2SBarry Smith PetscFunctionReturn(0); 133551f519a2SBarry Smith } 133651f519a2SBarry Smith EXTERN_C_END 133751f519a2SBarry Smith 133879416396SBarry Smith #undef __FUNCT__ 133979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1340bc08b0f1SBarry Smith /*@ 134179416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 134279416396SBarry Smith 134379416396SBarry Smith Collective on PC 134479416396SBarry Smith 134579416396SBarry Smith Input Parameter: 134679416396SBarry Smith . pc - the preconditioner context 134781540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 134879416396SBarry Smith 134979416396SBarry Smith Options Database Key: 1350a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 135179416396SBarry Smith 135279416396SBarry Smith Level: Developer 135379416396SBarry Smith 135479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 135579416396SBarry Smith 135679416396SBarry Smith .seealso: PCCompositeSetType() 135779416396SBarry Smith 135879416396SBarry Smith @*/ 13597087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 136079416396SBarry Smith { 13614ac538c5SBarry Smith PetscErrorCode ierr; 136279416396SBarry Smith 136379416396SBarry Smith PetscFunctionBegin; 13640700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13654ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 136679416396SBarry Smith PetscFunctionReturn(0); 136779416396SBarry Smith } 136879416396SBarry Smith 13690971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 13700971522cSBarry Smith /*MC 1371a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1372a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 13730971522cSBarry Smith 1374edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1375edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 13760971522cSBarry Smith 1377a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 137869a612a9SBarry Smith and set the options directly on the resulting KSP object 13790971522cSBarry Smith 13800971522cSBarry Smith Level: intermediate 13810971522cSBarry Smith 138279416396SBarry Smith Options Database Keys: 138381540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 138481540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 138581540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 138681540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 138781540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1388e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1389435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 139079416396SBarry Smith 1391edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 13923b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 13933b224e63SBarry Smith 1394d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1395d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1396d32f9abdSBarry Smith 1397d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1398d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1399d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1400d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1401d32f9abdSBarry Smith 1402a04f6461SBarry Smith For the Schur complement preconditioner if J = ( A00 A01 ) 1403a04f6461SBarry Smith ( A10 A11 ) 1404e69d4d44SBarry Smith the preconditioner is 1405a04f6461SBarry Smith (I -B ksp(A00)) ( inv(A00) 0 (I 0 ) 1406a04f6461SBarry Smith (0 I ) ( 0 ksp(S) ) (-A10 ksp(A00) I ) 1407a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1408a04f6461SBarry 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). 1409a04f6461SBarry Smith For PCFieldSplitGetKSP() when field number is 1410edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1411a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 14127e8cb189SBarry Smith option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc 1413e69d4d44SBarry Smith 1414edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1415edf189efSBarry Smith is used automatically for a second block. 1416edf189efSBarry Smith 14170716a85fSBarry 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. 14180716a85fSBarry Smith 1419a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 14200971522cSBarry Smith 14217e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1422e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 14230971522cSBarry Smith M*/ 14240971522cSBarry Smith 14250971522cSBarry Smith EXTERN_C_BEGIN 14260971522cSBarry Smith #undef __FUNCT__ 14270971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 14287087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 14290971522cSBarry Smith { 14300971522cSBarry Smith PetscErrorCode ierr; 14310971522cSBarry Smith PC_FieldSplit *jac; 14320971522cSBarry Smith 14330971522cSBarry Smith PetscFunctionBegin; 143438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 14350971522cSBarry Smith jac->bs = -1; 14360971522cSBarry Smith jac->nsplits = 0; 14373e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1438e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1439c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 144051f519a2SBarry Smith 14410971522cSBarry Smith pc->data = (void*)jac; 14420971522cSBarry Smith 14430971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1444421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 14450971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1446574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 14470971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 14480971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 14490971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 14500971522cSBarry Smith pc->ops->applyrichardson = 0; 14510971522cSBarry Smith 145269a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 145369a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14540971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 14550971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1456704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1457704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 145879416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 145979416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 146051f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 146151f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 14620971522cSBarry Smith PetscFunctionReturn(0); 14630971522cSBarry Smith } 14640971522cSBarry Smith EXTERN_C_END 14650971522cSBarry Smith 14660971522cSBarry Smith 1467a541d17aSBarry Smith 1468