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 */ 49c1570756SJed 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; 2222fa5ba8aSJed Brown DM *dms; 2238b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2248b8307b2SJed Brown ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 2258b8307b2SJed Brown ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 2262fa5ba8aSJed Brown ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr); 2272fa5ba8aSJed Brown ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr); 2288b8307b2SJed Brown for (i=0; i<nDM; i++) { 2292fa5ba8aSJed Brown char buf[256]; 2302fa5ba8aSJed Brown const char *splitname; 2312fa5ba8aSJed Brown /* Split naming precedence: object name, prefix, number */ 2322fa5ba8aSJed Brown splitname = ((PetscObject)pc->dm)->name; 2332fa5ba8aSJed Brown if (!splitname) { 2342fa5ba8aSJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr); 2352fa5ba8aSJed Brown if (splitname) { 2362fa5ba8aSJed Brown size_t len; 2372fa5ba8aSJed Brown ierr = PetscStrncpy(buf,splitname,sizeof buf);CHKERRQ(ierr); 2382fa5ba8aSJed Brown buf[sizeof buf - 1] = 0; 2392fa5ba8aSJed Brown ierr = PetscStrlen(buf,&len);CHKERRQ(ierr); 2402fa5ba8aSJed Brown if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */ 2412fa5ba8aSJed Brown splitname = buf; 2422fa5ba8aSJed Brown } 2432fa5ba8aSJed Brown } 2442fa5ba8aSJed Brown if (!splitname) { 2452fa5ba8aSJed Brown ierr = PetscSNPrintf(buf,sizeof buf,"%D",i);CHKERRQ(ierr); 2462fa5ba8aSJed Brown splitname = buf; 2472fa5ba8aSJed Brown } 2488b8307b2SJed Brown ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr); 249fcfd50ebSBarry Smith ierr = ISDestroy(&fields[i]);CHKERRQ(ierr); 2508b8307b2SJed Brown } 2518b8307b2SJed Brown ierr = PetscFree(fields);CHKERRQ(ierr); 2522fa5ba8aSJed Brown for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) { 2532fa5ba8aSJed Brown ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr); 2542fa5ba8aSJed Brown ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr); 2552fa5ba8aSJed Brown } 2562fa5ba8aSJed Brown ierr = PetscFree(dms);CHKERRQ(ierr); 2578b8307b2SJed Brown } 25866ffff09SJed Brown } else { 259521d7252SBarry Smith if (jac->bs <= 0) { 260704ba839SBarry Smith if (pc->pmat) { 261521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 262704ba839SBarry Smith } else { 263704ba839SBarry Smith jac->bs = 1; 264704ba839SBarry Smith } 265521d7252SBarry Smith } 266d32f9abdSBarry Smith 267acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 2686ce1633cSBarry Smith if (stokes) { 2696ce1633cSBarry Smith IS zerodiags,rest; 2706ce1633cSBarry Smith PetscInt nmin,nmax; 2716ce1633cSBarry Smith 2726ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 2736ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 2746ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 2756ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 2766ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 277fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 278fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 2796ce1633cSBarry Smith } else { 280d32f9abdSBarry Smith if (!flg) { 281d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 282d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2836c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2846c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 285d32f9abdSBarry Smith } 2866c924f48SJed Brown if (flg || !jac->splitdefined) { 287d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 288db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2896c924f48SJed Brown char splitname[8]; 2906c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 291db4c96c1SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr); 29279416396SBarry Smith } 29397bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 294521d7252SBarry Smith } 29566ffff09SJed Brown } 2966ce1633cSBarry Smith } 297edf189efSBarry Smith } else if (jac->nsplits == 1) { 298edf189efSBarry Smith if (ilink->is) { 299edf189efSBarry Smith IS is2; 300edf189efSBarry Smith PetscInt nmin,nmax; 301edf189efSBarry Smith 302edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 303edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 304db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 305fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 306e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 30763ec74ffSBarry Smith } else if (jac->reset) { 30863ec74ffSBarry Smith /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt 309d0af7cd3SBarry Smith This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed 310d0af7cd3SBarry Smith since they already exist. This should be totally rewritten */ 311d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 312d0af7cd3SBarry Smith if (pc->dm && !stokes) { 313d0af7cd3SBarry Smith PetscBool dmcomposite; 314d0af7cd3SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 315d0af7cd3SBarry Smith if (dmcomposite) { 316d0af7cd3SBarry Smith PetscInt nDM; 317d0af7cd3SBarry Smith IS *fields; 3182fa5ba8aSJed Brown DM *dms; 319d0af7cd3SBarry Smith ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 320d0af7cd3SBarry Smith ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 321d0af7cd3SBarry Smith ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 3222fa5ba8aSJed Brown ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr); 3232fa5ba8aSJed Brown ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr); 324d0af7cd3SBarry Smith for (i=0; i<nDM; i++) { 3252fa5ba8aSJed Brown ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr); 3262fa5ba8aSJed Brown ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr); 327d0af7cd3SBarry Smith ilink->is = fields[i]; 328d0af7cd3SBarry Smith ilink = ilink->next; 329edf189efSBarry Smith } 330d0af7cd3SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 3312fa5ba8aSJed Brown ierr = PetscFree(dms);CHKERRQ(ierr); 332d0af7cd3SBarry Smith } 333d0af7cd3SBarry Smith } else { 334d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 335d0af7cd3SBarry Smith if (stokes) { 336d0af7cd3SBarry Smith IS zerodiags,rest; 337d0af7cd3SBarry Smith PetscInt nmin,nmax; 338d0af7cd3SBarry Smith 339d0af7cd3SBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 340d0af7cd3SBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 341d0af7cd3SBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 34220b26d62SBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 34320b26d62SBarry Smith ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr); 344d0af7cd3SBarry Smith ilink->is = rest; 345d0af7cd3SBarry Smith ilink->next->is = zerodiags; 34620b26d62SBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 347d0af7cd3SBarry Smith } 348d0af7cd3SBarry Smith } 349d0af7cd3SBarry Smith 350e7e72b3dSBarry Smith if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 35169a612a9SBarry Smith PetscFunctionReturn(0); 35269a612a9SBarry Smith } 35369a612a9SBarry Smith 35469a612a9SBarry Smith #undef __FUNCT__ 35569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 35669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 35769a612a9SBarry Smith { 35869a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 35969a612a9SBarry Smith PetscErrorCode ierr; 3605a9f2f41SSatish Balay PC_FieldSplitLink ilink; 361cf502942SBarry Smith PetscInt i,nsplit,ccsize; 36269a612a9SBarry Smith MatStructure flag = pc->flag; 363ace3abfcSBarry Smith PetscBool sorted; 36469a612a9SBarry Smith 36569a612a9SBarry Smith PetscFunctionBegin; 36669a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 36797bbdb24SBarry Smith nsplit = jac->nsplits; 3685a9f2f41SSatish Balay ilink = jac->head; 36997bbdb24SBarry Smith 37097bbdb24SBarry Smith /* get the matrices for each split */ 371704ba839SBarry Smith if (!jac->issetup) { 3721b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 37397bbdb24SBarry Smith 374704ba839SBarry Smith jac->issetup = PETSC_TRUE; 375704ba839SBarry Smith 376704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 37751f519a2SBarry Smith bs = jac->bs; 37897bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 379cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 3801b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3811b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3821b9fc7fcSBarry Smith if (jac->defaultsplit) { 383704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 384704ba839SBarry Smith } else if (!ilink->is) { 385ccb205f8SBarry Smith if (ilink->nfields > 1) { 3865a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 3875a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3881b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3891b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3901b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 39197bbdb24SBarry Smith } 39297bbdb24SBarry Smith } 393d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 394ccb205f8SBarry Smith } else { 395704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 396ccb205f8SBarry Smith } 3973e197d65SBarry Smith } 398704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 399e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 400704ba839SBarry Smith ilink = ilink->next; 4011b9fc7fcSBarry Smith } 4021b9fc7fcSBarry Smith } 4031b9fc7fcSBarry Smith 404704ba839SBarry Smith ilink = jac->head; 40597bbdb24SBarry Smith if (!jac->pmat) { 406cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 407cf502942SBarry Smith for (i=0; i<nsplit; i++) { 4084aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 409704ba839SBarry Smith ilink = ilink->next; 410cf502942SBarry Smith } 41197bbdb24SBarry Smith } else { 412cf502942SBarry Smith for (i=0; i<nsplit; i++) { 4134aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 414704ba839SBarry Smith ilink = ilink->next; 415cf502942SBarry Smith } 41697bbdb24SBarry Smith } 417519d70e2SJed Brown if (jac->realdiagonal) { 418519d70e2SJed Brown ilink = jac->head; 419519d70e2SJed Brown if (!jac->mat) { 420519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 421519d70e2SJed Brown for (i=0; i<nsplit; i++) { 422519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 423519d70e2SJed Brown ilink = ilink->next; 424519d70e2SJed Brown } 425519d70e2SJed Brown } else { 426519d70e2SJed Brown for (i=0; i<nsplit; i++) { 427966e74d7SJed Brown if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 428519d70e2SJed Brown ilink = ilink->next; 429519d70e2SJed Brown } 430519d70e2SJed Brown } 431519d70e2SJed Brown } else { 432519d70e2SJed Brown jac->mat = jac->pmat; 433519d70e2SJed Brown } 43497bbdb24SBarry Smith 4356c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 43668dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 43768dd23aaSBarry Smith ilink = jac->head; 43868dd23aaSBarry Smith if (!jac->Afield) { 43968dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 44068dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4414aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 44268dd23aaSBarry Smith ilink = ilink->next; 44368dd23aaSBarry Smith } 44468dd23aaSBarry Smith } else { 44568dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 4464aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 44768dd23aaSBarry Smith ilink = ilink->next; 44868dd23aaSBarry Smith } 44968dd23aaSBarry Smith } 45068dd23aaSBarry Smith } 45168dd23aaSBarry Smith 4523b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 4533b224e63SBarry Smith IS ccis; 4544aa3045dSJed Brown PetscInt rstart,rend; 455e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 45668dd23aaSBarry Smith 457e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 458e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 459e6cab6aaSJed Brown 4603b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 4613b224e63SBarry Smith if (jac->schur) { 4623b224e63SBarry Smith ilink = jac->head; 4634aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4644aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 465fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4663b224e63SBarry Smith ilink = ilink->next; 4674aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4684aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 469fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 470519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 471084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4723b224e63SBarry Smith 4733b224e63SBarry Smith } else { 4741cee3971SBarry Smith KSP ksp; 4756c924f48SJed Brown char schurprefix[256]; 4763b224e63SBarry Smith 477a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 4783b224e63SBarry Smith ilink = jac->head; 4794aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4804aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 481fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4823b224e63SBarry Smith ilink = ilink->next; 4834aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4844aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 485fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 4867d50072eSJed Brown /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */ 487519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 488a04f6461SBarry Smith /* set tabbing and options prefix of KSP inside the MatSchur */ 4891cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 490e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 491a04f6461SBarry Smith ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr); 492a04f6461SBarry Smith ierr = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr); 49320b26d62SBarry Smith /* Need to call this everytime because new matrix is being created */ 49420b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4953b224e63SBarry Smith 4963b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4979005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4981cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 499084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 500084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 501084e4875SJed Brown PC pc; 502cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 503084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 504084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 505e69d4d44SBarry Smith } 5066c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 5076c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 5083b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 50920b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 51020b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 5113b224e63SBarry Smith 5123b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 5133b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 5143b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 5153b224e63SBarry Smith ilink = jac->head; 5163b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 5173b224e63SBarry Smith ilink = ilink->next; 5183b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 5193b224e63SBarry Smith } 5203b224e63SBarry Smith } else { 52197bbdb24SBarry Smith /* set up the individual PCs */ 52297bbdb24SBarry Smith i = 0; 5235a9f2f41SSatish Balay ilink = jac->head; 5245a9f2f41SSatish Balay while (ilink) { 525519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 5263b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 527c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 5285a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 52997bbdb24SBarry Smith i++; 5305a9f2f41SSatish Balay ilink = ilink->next; 5310971522cSBarry Smith } 53297bbdb24SBarry Smith 53397bbdb24SBarry Smith /* create work vectors for each split */ 5341b9fc7fcSBarry Smith if (!jac->x) { 53597bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 5365a9f2f41SSatish Balay ilink = jac->head; 53797bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 538906ed7ccSBarry Smith Vec *vl,*vr; 539906ed7ccSBarry Smith 540906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 541906ed7ccSBarry Smith ilink->x = *vr; 542906ed7ccSBarry Smith ilink->y = *vl; 543906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 544906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 5455a9f2f41SSatish Balay jac->x[i] = ilink->x; 5465a9f2f41SSatish Balay jac->y[i] = ilink->y; 5475a9f2f41SSatish Balay ilink = ilink->next; 54897bbdb24SBarry Smith } 5493b224e63SBarry Smith } 5503b224e63SBarry Smith } 5513b224e63SBarry Smith 5523b224e63SBarry Smith 553d0f46423SBarry Smith if (!jac->head->sctx) { 5543b224e63SBarry Smith Vec xtmp; 5553b224e63SBarry Smith 55679416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5571b9fc7fcSBarry Smith 5585a9f2f41SSatish Balay ilink = jac->head; 5591b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 5601b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 561704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 5625a9f2f41SSatish Balay ilink = ilink->next; 5631b9fc7fcSBarry Smith } 5646bf464f9SBarry Smith ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 5651b9fc7fcSBarry Smith } 566c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 5670971522cSBarry Smith PetscFunctionReturn(0); 5680971522cSBarry Smith } 5690971522cSBarry Smith 5705a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 571ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 572ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 5735a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 574ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 575ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 57679416396SBarry Smith 5770971522cSBarry Smith #undef __FUNCT__ 5783b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5793b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5803b224e63SBarry Smith { 5813b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5823b224e63SBarry Smith PetscErrorCode ierr; 5833b224e63SBarry Smith KSP ksp; 5843b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5853b224e63SBarry Smith 5863b224e63SBarry Smith PetscFunctionBegin; 5873b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5883b224e63SBarry Smith 589c5d2311dSJed Brown switch (jac->schurfactorization) { 590c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 591a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 592c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 593c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 594c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 595c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 596c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 597c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 598c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 599c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 600c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 601c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 602c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 603c5d2311dSJed Brown break; 604c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 605a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 606c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 607c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 608c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 609c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 610c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 611c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 612c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 613c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 614c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 615c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 616c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 617c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 618c5d2311dSJed Brown break; 619c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 620a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 621c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 622c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 623c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 624c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 625c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 626c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 627c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 628c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 629c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 630c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 631c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 632c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 633c5d2311dSJed Brown break; 634c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 635a04f6461SBarry 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 */ 6363b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6373b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6383b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6393b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 6403b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 6413b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6423b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6433b224e63SBarry Smith 6443b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 6453b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6463b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6473b224e63SBarry Smith 6483b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 6493b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 6503b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 6513b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6523b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 653c5d2311dSJed Brown } 6543b224e63SBarry Smith PetscFunctionReturn(0); 6553b224e63SBarry Smith } 6563b224e63SBarry Smith 6573b224e63SBarry Smith #undef __FUNCT__ 6580971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 6590971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 6600971522cSBarry Smith { 6610971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6620971522cSBarry Smith PetscErrorCode ierr; 6635a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 664af93645bSJed Brown PetscInt cnt; 6650971522cSBarry Smith 6660971522cSBarry Smith PetscFunctionBegin; 66751f519a2SBarry Smith CHKMEMQ; 66851f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 66951f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 67051f519a2SBarry Smith 67179416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 6721b9fc7fcSBarry Smith if (jac->defaultsplit) { 6730971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6745a9f2f41SSatish Balay while (ilink) { 6755a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6765a9f2f41SSatish Balay ilink = ilink->next; 6770971522cSBarry Smith } 6780971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6791b9fc7fcSBarry Smith } else { 680efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6815a9f2f41SSatish Balay while (ilink) { 6825a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6835a9f2f41SSatish Balay ilink = ilink->next; 6841b9fc7fcSBarry Smith } 6851b9fc7fcSBarry Smith } 68616913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 68779416396SBarry Smith if (!jac->w1) { 68879416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 68979416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 69079416396SBarry Smith } 691efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6925a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6933e197d65SBarry Smith cnt = 1; 6945a9f2f41SSatish Balay while (ilink->next) { 6955a9f2f41SSatish Balay ilink = ilink->next; 6963e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6973e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6983e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6993e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7003e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7013e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 7023e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7033e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7043e197d65SBarry Smith } 70551f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 70611755939SBarry Smith cnt -= 2; 70751f519a2SBarry Smith while (ilink->previous) { 70851f519a2SBarry Smith ilink = ilink->previous; 70911755939SBarry Smith /* compute the residual only over the part of the vector needed */ 71011755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 71111755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 71211755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 71311755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 71411755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 71511755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 71611755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 71751f519a2SBarry Smith } 71811755939SBarry Smith } 71965e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 72051f519a2SBarry Smith CHKMEMQ; 7210971522cSBarry Smith PetscFunctionReturn(0); 7220971522cSBarry Smith } 7230971522cSBarry Smith 724421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 725ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 726ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 727421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 728ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 729ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 730421e10b8SBarry Smith 731421e10b8SBarry Smith #undef __FUNCT__ 7328c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 733421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 734421e10b8SBarry Smith { 735421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 736421e10b8SBarry Smith PetscErrorCode ierr; 737421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 738421e10b8SBarry Smith 739421e10b8SBarry Smith PetscFunctionBegin; 740421e10b8SBarry Smith CHKMEMQ; 741421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 742421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 743421e10b8SBarry Smith 744421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 745421e10b8SBarry Smith if (jac->defaultsplit) { 746421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 747421e10b8SBarry Smith while (ilink) { 748421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 749421e10b8SBarry Smith ilink = ilink->next; 750421e10b8SBarry Smith } 751421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 752421e10b8SBarry Smith } else { 753421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 754421e10b8SBarry Smith while (ilink) { 755421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 756421e10b8SBarry Smith ilink = ilink->next; 757421e10b8SBarry Smith } 758421e10b8SBarry Smith } 759421e10b8SBarry Smith } else { 760421e10b8SBarry Smith if (!jac->w1) { 761421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 762421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 763421e10b8SBarry Smith } 764421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 765421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 766421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 767421e10b8SBarry Smith while (ilink->next) { 768421e10b8SBarry Smith ilink = ilink->next; 7699989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 770421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 771421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 772421e10b8SBarry Smith } 773421e10b8SBarry Smith while (ilink->previous) { 774421e10b8SBarry Smith ilink = ilink->previous; 7759989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 776421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 777421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 778421e10b8SBarry Smith } 779421e10b8SBarry Smith } else { 780421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 781421e10b8SBarry Smith ilink = ilink->next; 782421e10b8SBarry Smith } 783421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 784421e10b8SBarry Smith while (ilink->previous) { 785421e10b8SBarry Smith ilink = ilink->previous; 7869989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 787421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 788421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 789421e10b8SBarry Smith } 790421e10b8SBarry Smith } 791421e10b8SBarry Smith } 792421e10b8SBarry Smith CHKMEMQ; 793421e10b8SBarry Smith PetscFunctionReturn(0); 794421e10b8SBarry Smith } 795421e10b8SBarry Smith 7960971522cSBarry Smith #undef __FUNCT__ 797574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 798574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 7990971522cSBarry Smith { 8000971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8010971522cSBarry Smith PetscErrorCode ierr; 8025a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 8030971522cSBarry Smith 8040971522cSBarry Smith PetscFunctionBegin; 8055a9f2f41SSatish Balay while (ilink) { 806574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 807fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 808fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 809fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 810fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 8115a9f2f41SSatish Balay next = ilink->next; 8125a9f2f41SSatish Balay ilink = next; 8130971522cSBarry Smith } 81405b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 815574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 816574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 817574deadeSBarry Smith } else if (jac->mat) { 818574deadeSBarry Smith jac->mat = PETSC_NULL; 819574deadeSBarry Smith } 82097bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 82168dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 8226bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 8236bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 8246bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 8256bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 8266bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 8276bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 8286bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 82963ec74ffSBarry Smith jac->reset = PETSC_TRUE; 830574deadeSBarry Smith PetscFunctionReturn(0); 831574deadeSBarry Smith } 832574deadeSBarry Smith 833574deadeSBarry Smith #undef __FUNCT__ 834574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 835574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 836574deadeSBarry Smith { 837574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 838574deadeSBarry Smith PetscErrorCode ierr; 839574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 840574deadeSBarry Smith 841574deadeSBarry Smith PetscFunctionBegin; 842574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 843574deadeSBarry Smith while (ilink) { 8446bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 845574deadeSBarry Smith next = ilink->next; 846574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 847574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 848574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 849574deadeSBarry Smith ilink = next; 850574deadeSBarry Smith } 851574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 852c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 8530971522cSBarry Smith PetscFunctionReturn(0); 8540971522cSBarry Smith } 8550971522cSBarry Smith 8560971522cSBarry Smith #undef __FUNCT__ 8570971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 8580971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 8590971522cSBarry Smith { 8601b9fc7fcSBarry Smith PetscErrorCode ierr; 8616c924f48SJed Brown PetscInt bs; 862bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 8639dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8643b224e63SBarry Smith PCCompositeType ctype; 8651b9fc7fcSBarry Smith 8660971522cSBarry Smith PetscFunctionBegin; 8671b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 868acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 86951f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 87051f519a2SBarry Smith if (flg) { 87151f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 87251f519a2SBarry Smith } 873704ba839SBarry Smith 874435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 875c0adfefeSBarry Smith if (stokes) { 876c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 877c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 878c0adfefeSBarry Smith } 879c0adfefeSBarry Smith 8803b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 8813b224e63SBarry Smith if (flg) { 8823b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 8833b224e63SBarry Smith } 884d32f9abdSBarry Smith 885c30613efSMatthew Knepley /* Only setup fields once */ 886c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 887d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 888d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 8896c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 8906c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 891d32f9abdSBarry Smith } 892c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 893c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 894084e4875SJed 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); 895c5d2311dSJed Brown } 8961b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8970971522cSBarry Smith PetscFunctionReturn(0); 8980971522cSBarry Smith } 8990971522cSBarry Smith 9000971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 9010971522cSBarry Smith 9020971522cSBarry Smith EXTERN_C_BEGIN 9030971522cSBarry Smith #undef __FUNCT__ 9040971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 9057087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 9060971522cSBarry Smith { 90797bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9080971522cSBarry Smith PetscErrorCode ierr; 9095a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 91069a612a9SBarry Smith char prefix[128]; 91151f519a2SBarry Smith PetscInt i; 9120971522cSBarry Smith 9130971522cSBarry Smith PetscFunctionBegin; 9146c924f48SJed Brown if (jac->splitdefined) { 9156c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 9166c924f48SJed Brown PetscFunctionReturn(0); 9176c924f48SJed Brown } 91851f519a2SBarry Smith for (i=0; i<n; i++) { 919e32f2f54SBarry 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); 920e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 92151f519a2SBarry Smith } 922704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 923a04f6461SBarry Smith if (splitname) { 924db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 925a04f6461SBarry Smith } else { 926a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 927a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 928a04f6461SBarry Smith } 929704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 9305a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 9315a9f2f41SSatish Balay ilink->nfields = n; 9325a9f2f41SSatish Balay ilink->next = PETSC_NULL; 9337adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9341cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 9355a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9369005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 93769a612a9SBarry Smith 938a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 9395a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 9400971522cSBarry Smith 9410971522cSBarry Smith if (!next) { 9425a9f2f41SSatish Balay jac->head = ilink; 94351f519a2SBarry Smith ilink->previous = PETSC_NULL; 9440971522cSBarry Smith } else { 9450971522cSBarry Smith while (next->next) { 9460971522cSBarry Smith next = next->next; 9470971522cSBarry Smith } 9485a9f2f41SSatish Balay next->next = ilink; 94951f519a2SBarry Smith ilink->previous = next; 9500971522cSBarry Smith } 9510971522cSBarry Smith jac->nsplits++; 9520971522cSBarry Smith PetscFunctionReturn(0); 9530971522cSBarry Smith } 9540971522cSBarry Smith EXTERN_C_END 9550971522cSBarry Smith 956e69d4d44SBarry Smith EXTERN_C_BEGIN 957e69d4d44SBarry Smith #undef __FUNCT__ 958e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 9597087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 960e69d4d44SBarry Smith { 961e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 962e69d4d44SBarry Smith PetscErrorCode ierr; 963e69d4d44SBarry Smith 964e69d4d44SBarry Smith PetscFunctionBegin; 965519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 966e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 967e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 96813e0d083SBarry Smith if (n) *n = jac->nsplits; 969e69d4d44SBarry Smith PetscFunctionReturn(0); 970e69d4d44SBarry Smith } 971e69d4d44SBarry Smith EXTERN_C_END 9720971522cSBarry Smith 9730971522cSBarry Smith EXTERN_C_BEGIN 9740971522cSBarry Smith #undef __FUNCT__ 97569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 9767087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 9770971522cSBarry Smith { 9780971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9790971522cSBarry Smith PetscErrorCode ierr; 9800971522cSBarry Smith PetscInt cnt = 0; 9815a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 9820971522cSBarry Smith 9830971522cSBarry Smith PetscFunctionBegin; 9845d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 9855a9f2f41SSatish Balay while (ilink) { 9865a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 9875a9f2f41SSatish Balay ilink = ilink->next; 9880971522cSBarry Smith } 9895d480477SMatthew 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); 99013e0d083SBarry Smith if (n) *n = jac->nsplits; 9910971522cSBarry Smith PetscFunctionReturn(0); 9920971522cSBarry Smith } 9930971522cSBarry Smith EXTERN_C_END 9940971522cSBarry Smith 995704ba839SBarry Smith EXTERN_C_BEGIN 996704ba839SBarry Smith #undef __FUNCT__ 997704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 9987087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 999704ba839SBarry Smith { 1000704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1001704ba839SBarry Smith PetscErrorCode ierr; 1002704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1003704ba839SBarry Smith char prefix[128]; 1004704ba839SBarry Smith 1005704ba839SBarry Smith PetscFunctionBegin; 10066c924f48SJed Brown if (jac->splitdefined) { 10076c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10086c924f48SJed Brown PetscFunctionReturn(0); 10096c924f48SJed Brown } 101016913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1011a04f6461SBarry Smith if (splitname) { 1012db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1013a04f6461SBarry Smith } else { 1014a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1015a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1016a04f6461SBarry Smith } 10171ab39975SBarry Smith ilink->is = is; 10181ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1019704ba839SBarry Smith ilink->next = PETSC_NULL; 1020704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 10211cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1022704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 10239005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1024704ba839SBarry Smith 1025a04f6461SBarry Smith ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1026704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1027704ba839SBarry Smith 1028704ba839SBarry Smith if (!next) { 1029704ba839SBarry Smith jac->head = ilink; 1030704ba839SBarry Smith ilink->previous = PETSC_NULL; 1031704ba839SBarry Smith } else { 1032704ba839SBarry Smith while (next->next) { 1033704ba839SBarry Smith next = next->next; 1034704ba839SBarry Smith } 1035704ba839SBarry Smith next->next = ilink; 1036704ba839SBarry Smith ilink->previous = next; 1037704ba839SBarry Smith } 1038704ba839SBarry Smith jac->nsplits++; 1039704ba839SBarry Smith 1040704ba839SBarry Smith PetscFunctionReturn(0); 1041704ba839SBarry Smith } 1042704ba839SBarry Smith EXTERN_C_END 1043704ba839SBarry Smith 10440971522cSBarry Smith #undef __FUNCT__ 10450971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 10460971522cSBarry Smith /*@ 10470971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 10480971522cSBarry Smith 1049ad4df100SBarry Smith Logically Collective on PC 10500971522cSBarry Smith 10510971522cSBarry Smith Input Parameters: 10520971522cSBarry Smith + pc - the preconditioner context 1053a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 10540971522cSBarry Smith . n - the number of fields in this split 1055db4c96c1SJed Brown - fields - the fields in this split 10560971522cSBarry Smith 10570971522cSBarry Smith Level: intermediate 10580971522cSBarry Smith 1059d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1060d32f9abdSBarry Smith 1061d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 1062d32f9abdSBarry 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 1063d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1064d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1065d32f9abdSBarry Smith 1066db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1067db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1068db4c96c1SJed Brown 1069d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 10700971522cSBarry Smith 10710971522cSBarry Smith @*/ 10727087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 10730971522cSBarry Smith { 10744ac538c5SBarry Smith PetscErrorCode ierr; 10750971522cSBarry Smith 10760971522cSBarry Smith PetscFunctionBegin; 10770700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1078db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1079db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1080db4c96c1SJed Brown PetscValidIntPointer(fields,3); 10814ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr); 10820971522cSBarry Smith PetscFunctionReturn(0); 10830971522cSBarry Smith } 10840971522cSBarry Smith 10850971522cSBarry Smith #undef __FUNCT__ 1086704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1087704ba839SBarry Smith /*@ 1088704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1089704ba839SBarry Smith 1090ad4df100SBarry Smith Logically Collective on PC 1091704ba839SBarry Smith 1092704ba839SBarry Smith Input Parameters: 1093704ba839SBarry Smith + pc - the preconditioner context 1094a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1095db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1096704ba839SBarry Smith 1097d32f9abdSBarry Smith 1098a6ffb8dbSJed Brown Notes: 1099a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1100a6ffb8dbSJed Brown 1101db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1102db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1103d32f9abdSBarry Smith 1104704ba839SBarry Smith Level: intermediate 1105704ba839SBarry Smith 1106704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1107704ba839SBarry Smith 1108704ba839SBarry Smith @*/ 11097087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1110704ba839SBarry Smith { 11114ac538c5SBarry Smith PetscErrorCode ierr; 1112704ba839SBarry Smith 1113704ba839SBarry Smith PetscFunctionBegin; 11140700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1115db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1116db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 11174ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1118704ba839SBarry Smith PetscFunctionReturn(0); 1119704ba839SBarry Smith } 1120704ba839SBarry Smith 1121704ba839SBarry Smith #undef __FUNCT__ 112257a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 112357a9adfeSMatthew G Knepley /*@ 112457a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 112557a9adfeSMatthew G Knepley 112657a9adfeSMatthew G Knepley Logically Collective on PC 112757a9adfeSMatthew G Knepley 112857a9adfeSMatthew G Knepley Input Parameters: 112957a9adfeSMatthew G Knepley + pc - the preconditioner context 113057a9adfeSMatthew G Knepley - splitname - name of this split 113157a9adfeSMatthew G Knepley 113257a9adfeSMatthew G Knepley Output Parameter: 113357a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 113457a9adfeSMatthew G Knepley 113557a9adfeSMatthew G Knepley Level: intermediate 113657a9adfeSMatthew G Knepley 113757a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 113857a9adfeSMatthew G Knepley 113957a9adfeSMatthew G Knepley @*/ 114057a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 114157a9adfeSMatthew G Knepley { 114257a9adfeSMatthew G Knepley PetscErrorCode ierr; 114357a9adfeSMatthew G Knepley 114457a9adfeSMatthew G Knepley PetscFunctionBegin; 114557a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 114657a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 114757a9adfeSMatthew G Knepley PetscValidPointer(is,3); 114857a9adfeSMatthew G Knepley { 114957a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 115057a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 115157a9adfeSMatthew G Knepley PetscBool found; 115257a9adfeSMatthew G Knepley 115357a9adfeSMatthew G Knepley *is = PETSC_NULL; 115457a9adfeSMatthew G Knepley while(ilink) { 115557a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 115657a9adfeSMatthew G Knepley if (found) { 115757a9adfeSMatthew G Knepley *is = ilink->is; 115857a9adfeSMatthew G Knepley break; 115957a9adfeSMatthew G Knepley } 116057a9adfeSMatthew G Knepley ilink = ilink->next; 116157a9adfeSMatthew G Knepley } 116257a9adfeSMatthew G Knepley } 116357a9adfeSMatthew G Knepley PetscFunctionReturn(0); 116457a9adfeSMatthew G Knepley } 116557a9adfeSMatthew G Knepley 116657a9adfeSMatthew G Knepley #undef __FUNCT__ 116751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 116851f519a2SBarry Smith /*@ 116951f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 117051f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 117151f519a2SBarry Smith 1172ad4df100SBarry Smith Logically Collective on PC 117351f519a2SBarry Smith 117451f519a2SBarry Smith Input Parameters: 117551f519a2SBarry Smith + pc - the preconditioner context 117651f519a2SBarry Smith - bs - the block size 117751f519a2SBarry Smith 117851f519a2SBarry Smith Level: intermediate 117951f519a2SBarry Smith 118051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 118151f519a2SBarry Smith 118251f519a2SBarry Smith @*/ 11837087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 118451f519a2SBarry Smith { 11854ac538c5SBarry Smith PetscErrorCode ierr; 118651f519a2SBarry Smith 118751f519a2SBarry Smith PetscFunctionBegin; 11880700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1189c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 11904ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 119151f519a2SBarry Smith PetscFunctionReturn(0); 119251f519a2SBarry Smith } 119351f519a2SBarry Smith 119451f519a2SBarry Smith #undef __FUNCT__ 119569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 11960971522cSBarry Smith /*@C 119769a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 11980971522cSBarry Smith 119969a612a9SBarry Smith Collective on KSP 12000971522cSBarry Smith 12010971522cSBarry Smith Input Parameter: 12020971522cSBarry Smith . pc - the preconditioner context 12030971522cSBarry Smith 12040971522cSBarry Smith Output Parameters: 120513e0d083SBarry Smith + n - the number of splits 120669a612a9SBarry Smith - pc - the array of KSP contexts 12070971522cSBarry Smith 12080971522cSBarry Smith Note: 1209d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1210d32f9abdSBarry Smith (not the KSP just the array that contains them). 12110971522cSBarry Smith 121269a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 12130971522cSBarry Smith 12140971522cSBarry Smith Level: advanced 12150971522cSBarry Smith 12160971522cSBarry Smith .seealso: PCFIELDSPLIT 12170971522cSBarry Smith @*/ 12187087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 12190971522cSBarry Smith { 12204ac538c5SBarry Smith PetscErrorCode ierr; 12210971522cSBarry Smith 12220971522cSBarry Smith PetscFunctionBegin; 12230700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 122413e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 12254ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 12260971522cSBarry Smith PetscFunctionReturn(0); 12270971522cSBarry Smith } 12280971522cSBarry Smith 1229e69d4d44SBarry Smith #undef __FUNCT__ 1230e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1231e69d4d44SBarry Smith /*@ 1232e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1233a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1234e69d4d44SBarry Smith 1235e69d4d44SBarry Smith Collective on PC 1236e69d4d44SBarry Smith 1237e69d4d44SBarry Smith Input Parameters: 1238e69d4d44SBarry Smith + pc - the preconditioner context 1239fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1240084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1241084e4875SJed Brown 1242e69d4d44SBarry Smith Options Database: 1243084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1244e69d4d44SBarry Smith 1245fd1303e9SJungho Lee Notes: 1246fd1303e9SJungho Lee $ If ptype is 1247fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1248fd1303e9SJungho Lee $ to this function). 1249fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1250fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1251fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1252fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1253fd1303e9SJungho Lee $ preconditioner 1254fd1303e9SJungho Lee 1255fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1256fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1257fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1258fd1303e9SJungho Lee 1259fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1260fd1303e9SJungho Lee the name since it sets a proceedure to use. 1261fd1303e9SJungho Lee 1262e69d4d44SBarry Smith Level: intermediate 1263e69d4d44SBarry Smith 1264fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1265e69d4d44SBarry Smith 1266e69d4d44SBarry Smith @*/ 12677087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1268e69d4d44SBarry Smith { 12694ac538c5SBarry Smith PetscErrorCode ierr; 1270e69d4d44SBarry Smith 1271e69d4d44SBarry Smith PetscFunctionBegin; 12720700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12734ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1274e69d4d44SBarry Smith PetscFunctionReturn(0); 1275e69d4d44SBarry Smith } 1276e69d4d44SBarry Smith 1277e69d4d44SBarry Smith EXTERN_C_BEGIN 1278e69d4d44SBarry Smith #undef __FUNCT__ 1279e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 12807087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1281e69d4d44SBarry Smith { 1282e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1283084e4875SJed Brown PetscErrorCode ierr; 1284e69d4d44SBarry Smith 1285e69d4d44SBarry Smith PetscFunctionBegin; 1286084e4875SJed Brown jac->schurpre = ptype; 1287084e4875SJed Brown if (pre) { 12886bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1289084e4875SJed Brown jac->schur_user = pre; 1290084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1291084e4875SJed Brown } 1292e69d4d44SBarry Smith PetscFunctionReturn(0); 1293e69d4d44SBarry Smith } 1294e69d4d44SBarry Smith EXTERN_C_END 1295e69d4d44SBarry Smith 129630ad9308SMatthew Knepley #undef __FUNCT__ 129730ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 129830ad9308SMatthew Knepley /*@C 12998c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 130030ad9308SMatthew Knepley 130130ad9308SMatthew Knepley Collective on KSP 130230ad9308SMatthew Knepley 130330ad9308SMatthew Knepley Input Parameter: 130430ad9308SMatthew Knepley . pc - the preconditioner context 130530ad9308SMatthew Knepley 130630ad9308SMatthew Knepley Output Parameters: 1307a04f6461SBarry Smith + A00 - the (0,0) block 1308a04f6461SBarry Smith . A01 - the (0,1) block 1309a04f6461SBarry Smith . A10 - the (1,0) block 1310a04f6461SBarry Smith - A11 - the (1,1) block 131130ad9308SMatthew Knepley 131230ad9308SMatthew Knepley Level: advanced 131330ad9308SMatthew Knepley 131430ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 131530ad9308SMatthew Knepley @*/ 1316a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 131730ad9308SMatthew Knepley { 131830ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 131930ad9308SMatthew Knepley 132030ad9308SMatthew Knepley PetscFunctionBegin; 13210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1322c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1323a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1324a04f6461SBarry Smith if (A01) *A01 = jac->B; 1325a04f6461SBarry Smith if (A10) *A10 = jac->C; 1326a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 132730ad9308SMatthew Knepley PetscFunctionReturn(0); 132830ad9308SMatthew Knepley } 132930ad9308SMatthew Knepley 133079416396SBarry Smith EXTERN_C_BEGIN 133179416396SBarry Smith #undef __FUNCT__ 133279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 13337087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 133479416396SBarry Smith { 133579416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1336e69d4d44SBarry Smith PetscErrorCode ierr; 133779416396SBarry Smith 133879416396SBarry Smith PetscFunctionBegin; 133979416396SBarry Smith jac->type = type; 13403b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 13413b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 13423b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1343e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1344e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1345e69d4d44SBarry Smith 13463b224e63SBarry Smith } else { 13473b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 13483b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1349e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13509e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 13513b224e63SBarry Smith } 135279416396SBarry Smith PetscFunctionReturn(0); 135379416396SBarry Smith } 135479416396SBarry Smith EXTERN_C_END 135579416396SBarry Smith 135651f519a2SBarry Smith EXTERN_C_BEGIN 135751f519a2SBarry Smith #undef __FUNCT__ 135851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 13597087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 136051f519a2SBarry Smith { 136151f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 136251f519a2SBarry Smith 136351f519a2SBarry Smith PetscFunctionBegin; 136465e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 136565e19b50SBarry 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); 136651f519a2SBarry Smith jac->bs = bs; 136751f519a2SBarry Smith PetscFunctionReturn(0); 136851f519a2SBarry Smith } 136951f519a2SBarry Smith EXTERN_C_END 137051f519a2SBarry Smith 137179416396SBarry Smith #undef __FUNCT__ 137279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1373bc08b0f1SBarry Smith /*@ 137479416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 137579416396SBarry Smith 137679416396SBarry Smith Collective on PC 137779416396SBarry Smith 137879416396SBarry Smith Input Parameter: 137979416396SBarry Smith . pc - the preconditioner context 138081540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 138179416396SBarry Smith 138279416396SBarry Smith Options Database Key: 1383a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 138479416396SBarry Smith 138579416396SBarry Smith Level: Developer 138679416396SBarry Smith 138779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 138879416396SBarry Smith 138979416396SBarry Smith .seealso: PCCompositeSetType() 139079416396SBarry Smith 139179416396SBarry Smith @*/ 13927087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 139379416396SBarry Smith { 13944ac538c5SBarry Smith PetscErrorCode ierr; 139579416396SBarry Smith 139679416396SBarry Smith PetscFunctionBegin; 13970700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13984ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 139979416396SBarry Smith PetscFunctionReturn(0); 140079416396SBarry Smith } 140179416396SBarry Smith 14020971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 14030971522cSBarry Smith /*MC 1404a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1405a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 14060971522cSBarry Smith 1407edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1408edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 14090971522cSBarry Smith 1410a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 141169a612a9SBarry Smith and set the options directly on the resulting KSP object 14120971522cSBarry Smith 14130971522cSBarry Smith Level: intermediate 14140971522cSBarry Smith 141579416396SBarry Smith Options Database Keys: 141681540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 141781540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 141881540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 141981540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 142081540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1421e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 1422435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 142379416396SBarry Smith 1424edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 14253b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 14263b224e63SBarry Smith 1427d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1428d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1429d32f9abdSBarry Smith 1430d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1431d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1432d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1433d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1434d32f9abdSBarry Smith 1435a04f6461SBarry Smith For the Schur complement preconditioner if J = ( A00 A01 ) 1436a04f6461SBarry Smith ( A10 A11 ) 1437e69d4d44SBarry Smith the preconditioner is 1438a04f6461SBarry Smith (I -B ksp(A00)) ( inv(A00) 0 (I 0 ) 1439a04f6461SBarry Smith (0 I ) ( 0 ksp(S) ) (-A10 ksp(A00) I ) 1440a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1441a04f6461SBarry 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). 1442a04f6461SBarry Smith For PCFieldSplitGetKSP() when field number is 1443edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1444a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 14457e8cb189SBarry Smith option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc 1446e69d4d44SBarry Smith 1447edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1448edf189efSBarry Smith is used automatically for a second block. 1449edf189efSBarry Smith 1450*ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1451*ff218e97SBarry Smith Generally it should be used with the AIJ format. 1452*ff218e97SBarry Smith 1453*ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1454*ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1455*ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 14560716a85fSBarry Smith 1457a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 14580971522cSBarry Smith 14597e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1460e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 14610971522cSBarry Smith M*/ 14620971522cSBarry Smith 14630971522cSBarry Smith EXTERN_C_BEGIN 14640971522cSBarry Smith #undef __FUNCT__ 14650971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 14667087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 14670971522cSBarry Smith { 14680971522cSBarry Smith PetscErrorCode ierr; 14690971522cSBarry Smith PC_FieldSplit *jac; 14700971522cSBarry Smith 14710971522cSBarry Smith PetscFunctionBegin; 147238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 14730971522cSBarry Smith jac->bs = -1; 14740971522cSBarry Smith jac->nsplits = 0; 14753e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1476e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1477c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 147851f519a2SBarry Smith 14790971522cSBarry Smith pc->data = (void*)jac; 14800971522cSBarry Smith 14810971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1482421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 14830971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1484574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 14850971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 14860971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 14870971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 14880971522cSBarry Smith pc->ops->applyrichardson = 0; 14890971522cSBarry Smith 149069a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 149169a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 14920971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 14930971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1494704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1495704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 149679416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 149779416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 149851f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 149951f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 15000971522cSBarry Smith PetscFunctionReturn(0); 15010971522cSBarry Smith } 15020971522cSBarry Smith EXTERN_C_END 15030971522cSBarry Smith 15040971522cSBarry Smith 1505a541d17aSBarry Smith 1506