1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 36356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.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 */ 4230ad9308SMatthew Knepley Mat schur; /* The Schur complement S = D - C A^{-1} B */ 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; 480971522cSBarry Smith } PC_FieldSplit; 490971522cSBarry Smith 5016913363SBarry Smith /* 5116913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5216913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5316913363SBarry Smith PC you could change this. 5416913363SBarry Smith */ 55084e4875SJed Brown 56e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 57084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 58084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 59084e4875SJed Brown { 60084e4875SJed Brown switch (jac->schurpre) { 61084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 62084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 63084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 64084e4875SJed Brown default: 65084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 66084e4875SJed Brown } 67084e4875SJed Brown } 68084e4875SJed Brown 69084e4875SJed Brown 700971522cSBarry Smith #undef __FUNCT__ 710971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 720971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 730971522cSBarry Smith { 740971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 750971522cSBarry Smith PetscErrorCode ierr; 76ace3abfcSBarry Smith PetscBool iascii; 770971522cSBarry Smith PetscInt i,j; 785a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 790971522cSBarry Smith 800971522cSBarry Smith PetscFunctionBegin; 812692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 820971522cSBarry Smith if (iascii) { 8351f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 8469a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 850971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 860971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 871ab39975SBarry Smith if (ilink->fields) { 880971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 8979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 905a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 9179416396SBarry Smith if (j > 0) { 9279416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 9379416396SBarry Smith } 945a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 950971522cSBarry Smith } 960971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9779416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 981ab39975SBarry Smith } else { 991ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1001ab39975SBarry Smith } 1015a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1025a9f2f41SSatish Balay ilink = ilink->next; 1030971522cSBarry Smith } 1040971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1050971522cSBarry Smith } else { 10665e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1070971522cSBarry Smith } 1080971522cSBarry Smith PetscFunctionReturn(0); 1090971522cSBarry Smith } 1100971522cSBarry Smith 1110971522cSBarry Smith #undef __FUNCT__ 1123b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1133b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1143b224e63SBarry Smith { 1153b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1163b224e63SBarry Smith PetscErrorCode ierr; 117ace3abfcSBarry Smith PetscBool iascii; 1183b224e63SBarry Smith PetscInt i,j; 1193b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1203b224e63SBarry Smith KSP ksp; 1213b224e63SBarry Smith 1223b224e63SBarry Smith PetscFunctionBegin; 1232692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1243b224e63SBarry Smith if (iascii) { 125c5d2311dSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr); 1263b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1273b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1283b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1293b224e63SBarry Smith if (ilink->fields) { 1303b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1313b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1323b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1333b224e63SBarry Smith if (j > 0) { 1343b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1353b224e63SBarry Smith } 1363b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1373b224e63SBarry Smith } 1383b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1393b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1403b224e63SBarry Smith } else { 1413b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1423b224e63SBarry Smith } 1433b224e63SBarry Smith ilink = ilink->next; 1443b224e63SBarry Smith } 1453b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1463b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14712cae6f2SJed Brown if (jac->schur) { 14812cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1493b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 15012cae6f2SJed Brown } else { 15112cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 15212cae6f2SJed Brown } 1533b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1543b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1553b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 15612cae6f2SJed Brown if (jac->kspschur) { 1573b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 15812cae6f2SJed Brown } else { 15912cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 16012cae6f2SJed Brown } 1613b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1623b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1633b224e63SBarry Smith } else { 16465e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1653b224e63SBarry Smith } 1663b224e63SBarry Smith PetscFunctionReturn(0); 1673b224e63SBarry Smith } 1683b224e63SBarry Smith 1693b224e63SBarry Smith #undef __FUNCT__ 1706c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1716c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1726c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1736c924f48SJed Brown { 1746c924f48SJed Brown PetscErrorCode ierr; 1756c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1766c924f48SJed Brown PetscInt i,nfields,*ifields; 177ace3abfcSBarry Smith PetscBool flg; 1786c924f48SJed Brown char optionname[128],splitname[8]; 1796c924f48SJed Brown 1806c924f48SJed Brown PetscFunctionBegin; 1816c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 1826c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 1836c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 1846c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 1856c924f48SJed Brown nfields = jac->bs; 1866c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 1876c924f48SJed Brown if (!flg) break; 1886c924f48SJed Brown if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1896c924f48SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr); 1906c924f48SJed Brown } 1916c924f48SJed Brown if (i > 0) { 1926c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 1936c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 1946c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 1956c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 1966c924f48SJed Brown } 1976c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 1986c924f48SJed Brown PetscFunctionReturn(0); 1996c924f48SJed Brown } 2006c924f48SJed Brown 2016c924f48SJed Brown #undef __FUNCT__ 20269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 20369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2040971522cSBarry Smith { 2050971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2060971522cSBarry Smith PetscErrorCode ierr; 2075a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 208*6ce1633cSBarry Smith PetscBool flg = PETSC_FALSE,stokes = PETSC_FALSE; 2096c924f48SJed Brown PetscInt i; 2100971522cSBarry Smith 2110971522cSBarry Smith PetscFunctionBegin; 212d32f9abdSBarry Smith if (!ilink) { 2138b8307b2SJed Brown if (pc->dm) { 2148b8307b2SJed Brown PetscBool dmcomposite; 2158b8307b2SJed Brown ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr); 2168b8307b2SJed Brown if (dmcomposite) { 2178b8307b2SJed Brown PetscInt nDM; 2188b8307b2SJed Brown IS *fields; 2198b8307b2SJed Brown ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 2208b8307b2SJed Brown ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr); 2218b8307b2SJed Brown ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr); 2228b8307b2SJed Brown for (i=0; i<nDM; i++) { 2238b8307b2SJed Brown char splitname[8]; 2248b8307b2SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 2258b8307b2SJed Brown ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr); 2268b8307b2SJed Brown ierr = ISDestroy(fields[i]);CHKERRQ(ierr); 2278b8307b2SJed Brown } 2288b8307b2SJed Brown ierr = PetscFree(fields);CHKERRQ(ierr); 2298b8307b2SJed Brown } 23066ffff09SJed Brown } else { 231521d7252SBarry Smith if (jac->bs <= 0) { 232704ba839SBarry Smith if (pc->pmat) { 233521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 234704ba839SBarry Smith } else { 235704ba839SBarry Smith jac->bs = 1; 236704ba839SBarry Smith } 237521d7252SBarry Smith } 238d32f9abdSBarry Smith 239acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 240*6ce1633cSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_stokes",&stokes,PETSC_NULL);CHKERRQ(ierr); 241*6ce1633cSBarry Smith if (stokes) { 242*6ce1633cSBarry Smith IS zerodiags,rest; 243*6ce1633cSBarry Smith PetscInt nmin,nmax; 244*6ce1633cSBarry Smith 245*6ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 246*6ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 247*6ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 248*6ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 249*6ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 250*6ce1633cSBarry Smith ierr = ISDestroy(zerodiags);CHKERRQ(ierr); 251*6ce1633cSBarry Smith ierr = ISDestroy(rest);CHKERRQ(ierr); 252*6ce1633cSBarry Smith } else { 253d32f9abdSBarry Smith if (!flg) { 254d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 255d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2566c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2576c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 258d32f9abdSBarry Smith } 2596c924f48SJed Brown if (flg || !jac->splitdefined) { 260d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 261db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2626c924f48SJed Brown char splitname[8]; 2636c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 264db4c96c1SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr); 26579416396SBarry Smith } 26697bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 267521d7252SBarry Smith } 26866ffff09SJed Brown } 269*6ce1633cSBarry Smith } 270edf189efSBarry Smith } else if (jac->nsplits == 1) { 271edf189efSBarry Smith if (ilink->is) { 272edf189efSBarry Smith IS is2; 273edf189efSBarry Smith PetscInt nmin,nmax; 274edf189efSBarry Smith 275edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 276edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 277db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 278edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 279e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 280edf189efSBarry Smith } 281e7e72b3dSBarry Smith if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 28269a612a9SBarry Smith PetscFunctionReturn(0); 28369a612a9SBarry Smith } 28469a612a9SBarry Smith 28569a612a9SBarry Smith 28669a612a9SBarry Smith #undef __FUNCT__ 28769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 28869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 28969a612a9SBarry Smith { 29069a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 29169a612a9SBarry Smith PetscErrorCode ierr; 2925a9f2f41SSatish Balay PC_FieldSplitLink ilink; 293cf502942SBarry Smith PetscInt i,nsplit,ccsize; 29469a612a9SBarry Smith MatStructure flag = pc->flag; 295ace3abfcSBarry Smith PetscBool sorted; 29669a612a9SBarry Smith 29769a612a9SBarry Smith PetscFunctionBegin; 29869a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 29997bbdb24SBarry Smith nsplit = jac->nsplits; 3005a9f2f41SSatish Balay ilink = jac->head; 30197bbdb24SBarry Smith 30297bbdb24SBarry Smith /* get the matrices for each split */ 303704ba839SBarry Smith if (!jac->issetup) { 3041b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 30597bbdb24SBarry Smith 306704ba839SBarry Smith jac->issetup = PETSC_TRUE; 307704ba839SBarry Smith 308704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 30951f519a2SBarry Smith bs = jac->bs; 31097bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 311cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 3121b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 3131b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 3141b9fc7fcSBarry Smith if (jac->defaultsplit) { 315704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 316704ba839SBarry Smith } else if (!ilink->is) { 317ccb205f8SBarry Smith if (ilink->nfields > 1) { 3185a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 3195a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 3201b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 3211b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 3221b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 32397bbdb24SBarry Smith } 32497bbdb24SBarry Smith } 325d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 326ccb205f8SBarry Smith } else { 327704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 328ccb205f8SBarry Smith } 3293e197d65SBarry Smith } 330704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 331e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 332704ba839SBarry Smith ilink = ilink->next; 3331b9fc7fcSBarry Smith } 3341b9fc7fcSBarry Smith } 3351b9fc7fcSBarry Smith 336704ba839SBarry Smith ilink = jac->head; 33797bbdb24SBarry Smith if (!jac->pmat) { 338cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 339cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3404aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 341704ba839SBarry Smith ilink = ilink->next; 342cf502942SBarry Smith } 34397bbdb24SBarry Smith } else { 344cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3454aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 346704ba839SBarry Smith ilink = ilink->next; 347cf502942SBarry Smith } 34897bbdb24SBarry Smith } 349519d70e2SJed Brown if (jac->realdiagonal) { 350519d70e2SJed Brown ilink = jac->head; 351519d70e2SJed Brown if (!jac->mat) { 352519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 353519d70e2SJed Brown for (i=0; i<nsplit; i++) { 354519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 355519d70e2SJed Brown ilink = ilink->next; 356519d70e2SJed Brown } 357519d70e2SJed Brown } else { 358519d70e2SJed Brown for (i=0; i<nsplit; i++) { 359519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 360519d70e2SJed Brown ilink = ilink->next; 361519d70e2SJed Brown } 362519d70e2SJed Brown } 363519d70e2SJed Brown } else { 364519d70e2SJed Brown jac->mat = jac->pmat; 365519d70e2SJed Brown } 36697bbdb24SBarry Smith 3676c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 36868dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 36968dd23aaSBarry Smith ilink = jac->head; 37068dd23aaSBarry Smith if (!jac->Afield) { 37168dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 37268dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3734aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 37468dd23aaSBarry Smith ilink = ilink->next; 37568dd23aaSBarry Smith } 37668dd23aaSBarry Smith } else { 37768dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3784aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 37968dd23aaSBarry Smith ilink = ilink->next; 38068dd23aaSBarry Smith } 38168dd23aaSBarry Smith } 38268dd23aaSBarry Smith } 38368dd23aaSBarry Smith 3843b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 3853b224e63SBarry Smith IS ccis; 3864aa3045dSJed Brown PetscInt rstart,rend; 387e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 38868dd23aaSBarry Smith 389e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 390e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 391e6cab6aaSJed Brown 3923b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 3933b224e63SBarry Smith if (jac->schur) { 3943b224e63SBarry Smith ilink = jac->head; 3954aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3964aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3973b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3983b224e63SBarry Smith ilink = ilink->next; 3994aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4004aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 4013b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 402519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 403084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 4043b224e63SBarry Smith 4053b224e63SBarry Smith } else { 4061cee3971SBarry Smith KSP ksp; 4076c924f48SJed Brown char schurprefix[256]; 4083b224e63SBarry Smith 4093b224e63SBarry Smith /* extract the B and C matrices */ 4103b224e63SBarry Smith ilink = jac->head; 4114aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4124aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 4133b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 4143b224e63SBarry Smith ilink = ilink->next; 4154aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 4164aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 4173b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 418084e4875SJed Brown /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 419519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 4201cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 421e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 4223b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 4233b224e63SBarry Smith 4243b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 4259005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 4261cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 427084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 428084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 429084e4875SJed Brown PC pc; 430cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 431084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 432084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 433e69d4d44SBarry Smith } 4346c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 4356c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 4363b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 4373b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 4383b224e63SBarry Smith 4393b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 4403b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 4413b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 4423b224e63SBarry Smith ilink = jac->head; 4433b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 4443b224e63SBarry Smith ilink = ilink->next; 4453b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 4463b224e63SBarry Smith } 4473b224e63SBarry Smith } else { 44897bbdb24SBarry Smith /* set up the individual PCs */ 44997bbdb24SBarry Smith i = 0; 4505a9f2f41SSatish Balay ilink = jac->head; 4515a9f2f41SSatish Balay while (ilink) { 452519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 4533b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 4545a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 4555a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 45697bbdb24SBarry Smith i++; 4575a9f2f41SSatish Balay ilink = ilink->next; 4580971522cSBarry Smith } 45997bbdb24SBarry Smith 46097bbdb24SBarry Smith /* create work vectors for each split */ 4611b9fc7fcSBarry Smith if (!jac->x) { 46297bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 4635a9f2f41SSatish Balay ilink = jac->head; 46497bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 465906ed7ccSBarry Smith Vec *vl,*vr; 466906ed7ccSBarry Smith 467906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 468906ed7ccSBarry Smith ilink->x = *vr; 469906ed7ccSBarry Smith ilink->y = *vl; 470906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 471906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 4725a9f2f41SSatish Balay jac->x[i] = ilink->x; 4735a9f2f41SSatish Balay jac->y[i] = ilink->y; 4745a9f2f41SSatish Balay ilink = ilink->next; 47597bbdb24SBarry Smith } 4763b224e63SBarry Smith } 4773b224e63SBarry Smith } 4783b224e63SBarry Smith 4793b224e63SBarry Smith 480d0f46423SBarry Smith if (!jac->head->sctx) { 4813b224e63SBarry Smith Vec xtmp; 4823b224e63SBarry Smith 48379416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4841b9fc7fcSBarry Smith 4855a9f2f41SSatish Balay ilink = jac->head; 4861b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4871b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 488704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4895a9f2f41SSatish Balay ilink = ilink->next; 4901b9fc7fcSBarry Smith } 4911b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4921b9fc7fcSBarry Smith } 4930971522cSBarry Smith PetscFunctionReturn(0); 4940971522cSBarry Smith } 4950971522cSBarry Smith 4965a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 497ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 498ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4995a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 500ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 501ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 50279416396SBarry Smith 5030971522cSBarry Smith #undef __FUNCT__ 5043b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 5053b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 5063b224e63SBarry Smith { 5073b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5083b224e63SBarry Smith PetscErrorCode ierr; 5093b224e63SBarry Smith KSP ksp; 5103b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 5113b224e63SBarry Smith 5123b224e63SBarry Smith PetscFunctionBegin; 5133b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 5143b224e63SBarry Smith 515c5d2311dSJed Brown switch (jac->schurfactorization) { 516c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG: 517c5d2311dSJed Brown /* [A 0; 0 -S], positive definite, suitable for MINRES */ 518c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 519c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 520c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 521c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 522c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 523c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 524c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 525c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 526c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 527c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 528c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 529c5d2311dSJed Brown break; 530c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER: 531c5d2311dSJed Brown /* [A 0; C S], suitable for left preconditioning */ 532c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 533c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 534c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 535c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 536c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 537c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 538c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 539c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 540c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 541c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 542c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 543c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 544c5d2311dSJed Brown break; 545c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER: 546c5d2311dSJed Brown /* [A B; 0 S], suitable for right preconditioning */ 547c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 548c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 549c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 550c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 551c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 552c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 553c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 554c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 555c5d2311dSJed Brown ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 556c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 557c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 558c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 559c5d2311dSJed Brown break; 560c5d2311dSJed Brown case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL: 561c5d2311dSJed Brown /* [1 0; CA^{-1} 1] [A 0; 0 S] [1 A^{-1}B; 0 1], an exact solve if applied exactly, needs one extra solve with A */ 5623b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5633b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5643b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 5653b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 5663b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 5673b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5683b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5693b224e63SBarry Smith 5703b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 5713b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5723b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5733b224e63SBarry Smith 5743b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 5753b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 5763b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 5773b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5783b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 579c5d2311dSJed Brown } 5803b224e63SBarry Smith PetscFunctionReturn(0); 5813b224e63SBarry Smith } 5823b224e63SBarry Smith 5833b224e63SBarry Smith #undef __FUNCT__ 5840971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 5850971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 5860971522cSBarry Smith { 5870971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5880971522cSBarry Smith PetscErrorCode ierr; 5895a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 590af93645bSJed Brown PetscInt cnt; 5910971522cSBarry Smith 5920971522cSBarry Smith PetscFunctionBegin; 59351f519a2SBarry Smith CHKMEMQ; 59451f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 59551f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 59651f519a2SBarry Smith 59779416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 5981b9fc7fcSBarry Smith if (jac->defaultsplit) { 5990971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 6005a9f2f41SSatish Balay while (ilink) { 6015a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6025a9f2f41SSatish Balay ilink = ilink->next; 6030971522cSBarry Smith } 6040971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 6051b9fc7fcSBarry Smith } else { 606efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6075a9f2f41SSatish Balay while (ilink) { 6085a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6095a9f2f41SSatish Balay ilink = ilink->next; 6101b9fc7fcSBarry Smith } 6111b9fc7fcSBarry Smith } 61216913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 61379416396SBarry Smith if (!jac->w1) { 61479416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 61579416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 61679416396SBarry Smith } 617efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 6185a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 6193e197d65SBarry Smith cnt = 1; 6205a9f2f41SSatish Balay while (ilink->next) { 6215a9f2f41SSatish Balay ilink = ilink->next; 6223e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 6233e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 6243e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 6253e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6263e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6273e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 6283e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6293e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6303e197d65SBarry Smith } 63151f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 63211755939SBarry Smith cnt -= 2; 63351f519a2SBarry Smith while (ilink->previous) { 63451f519a2SBarry Smith ilink = ilink->previous; 63511755939SBarry Smith /* compute the residual only over the part of the vector needed */ 63611755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 63711755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 63811755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 63911755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 64011755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 64111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 64211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 64351f519a2SBarry Smith } 64411755939SBarry Smith } 64565e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 64651f519a2SBarry Smith CHKMEMQ; 6470971522cSBarry Smith PetscFunctionReturn(0); 6480971522cSBarry Smith } 6490971522cSBarry Smith 650421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 651ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 652ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 653421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 654ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 655ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 656421e10b8SBarry Smith 657421e10b8SBarry Smith #undef __FUNCT__ 6588c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 659421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 660421e10b8SBarry Smith { 661421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 662421e10b8SBarry Smith PetscErrorCode ierr; 663421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 664421e10b8SBarry Smith 665421e10b8SBarry Smith PetscFunctionBegin; 666421e10b8SBarry Smith CHKMEMQ; 667421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 668421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 669421e10b8SBarry Smith 670421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 671421e10b8SBarry Smith if (jac->defaultsplit) { 672421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 673421e10b8SBarry Smith while (ilink) { 674421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 675421e10b8SBarry Smith ilink = ilink->next; 676421e10b8SBarry Smith } 677421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 678421e10b8SBarry Smith } else { 679421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 680421e10b8SBarry Smith while (ilink) { 681421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 682421e10b8SBarry Smith ilink = ilink->next; 683421e10b8SBarry Smith } 684421e10b8SBarry Smith } 685421e10b8SBarry Smith } else { 686421e10b8SBarry Smith if (!jac->w1) { 687421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 688421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 689421e10b8SBarry Smith } 690421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 691421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 692421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 693421e10b8SBarry Smith while (ilink->next) { 694421e10b8SBarry Smith ilink = ilink->next; 6959989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 696421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 697421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 698421e10b8SBarry Smith } 699421e10b8SBarry Smith while (ilink->previous) { 700421e10b8SBarry Smith ilink = ilink->previous; 7019989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 702421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 703421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 704421e10b8SBarry Smith } 705421e10b8SBarry Smith } else { 706421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 707421e10b8SBarry Smith ilink = ilink->next; 708421e10b8SBarry Smith } 709421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 710421e10b8SBarry Smith while (ilink->previous) { 711421e10b8SBarry Smith ilink = ilink->previous; 7129989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 713421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 714421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 715421e10b8SBarry Smith } 716421e10b8SBarry Smith } 717421e10b8SBarry Smith } 718421e10b8SBarry Smith CHKMEMQ; 719421e10b8SBarry Smith PetscFunctionReturn(0); 720421e10b8SBarry Smith } 721421e10b8SBarry Smith 7220971522cSBarry Smith #undef __FUNCT__ 7230971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 7240971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 7250971522cSBarry Smith { 7260971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7270971522cSBarry Smith PetscErrorCode ierr; 7285a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 7290971522cSBarry Smith 7300971522cSBarry Smith PetscFunctionBegin; 7315a9f2f41SSatish Balay while (ilink) { 7325a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 7335a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 7345a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 7355a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 736704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 7375a9f2f41SSatish Balay next = ilink->next; 7387e8c30b6SJed Brown ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 739704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 740704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 7415a9f2f41SSatish Balay ilink = next; 7420971522cSBarry Smith } 74305b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 744519d70e2SJed Brown if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 74597bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 74668dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 74779416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 74879416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 7493b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 750084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 751d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 7523b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 7533b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 7540971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 7550971522cSBarry Smith PetscFunctionReturn(0); 7560971522cSBarry Smith } 7570971522cSBarry Smith 7580971522cSBarry Smith #undef __FUNCT__ 7590971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 7600971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 7610971522cSBarry Smith { 7621b9fc7fcSBarry Smith PetscErrorCode ierr; 7636c924f48SJed Brown PetscInt bs; 764ace3abfcSBarry Smith PetscBool flg; 7659dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7663b224e63SBarry Smith PCCompositeType ctype; 7671b9fc7fcSBarry Smith 7680971522cSBarry Smith PetscFunctionBegin; 7691b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 770acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 77151f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 77251f519a2SBarry Smith if (flg) { 77351f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 77451f519a2SBarry Smith } 775704ba839SBarry Smith 7763b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 7773b224e63SBarry Smith if (flg) { 7783b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 7793b224e63SBarry Smith } 780d32f9abdSBarry Smith 781c30613efSMatthew Knepley /* Only setup fields once */ 782c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 783d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 784d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 7856c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 7866c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 787d32f9abdSBarry Smith } 788c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 789c5d2311dSJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr); 790084e4875SJed 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); 791c5d2311dSJed Brown } 7921b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7930971522cSBarry Smith PetscFunctionReturn(0); 7940971522cSBarry Smith } 7950971522cSBarry Smith 7960971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 7970971522cSBarry Smith 7980971522cSBarry Smith EXTERN_C_BEGIN 7990971522cSBarry Smith #undef __FUNCT__ 8000971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 8017087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 8020971522cSBarry Smith { 80397bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8040971522cSBarry Smith PetscErrorCode ierr; 8055a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 80669a612a9SBarry Smith char prefix[128]; 80751f519a2SBarry Smith PetscInt i; 8080971522cSBarry Smith 8090971522cSBarry Smith PetscFunctionBegin; 8106c924f48SJed Brown if (jac->splitdefined) { 8116c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 8126c924f48SJed Brown PetscFunctionReturn(0); 8136c924f48SJed Brown } 81451f519a2SBarry Smith for (i=0; i<n; i++) { 815e32f2f54SBarry 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); 816e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 81751f519a2SBarry Smith } 818704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 819db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 820704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 8215a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 8225a9f2f41SSatish Balay ilink->nfields = n; 8235a9f2f41SSatish Balay ilink->next = PETSC_NULL; 8247adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 8251cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 8265a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 8279005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 82869a612a9SBarry Smith 8296c924f48SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr); 8305a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 8310971522cSBarry Smith 8320971522cSBarry Smith if (!next) { 8335a9f2f41SSatish Balay jac->head = ilink; 83451f519a2SBarry Smith ilink->previous = PETSC_NULL; 8350971522cSBarry Smith } else { 8360971522cSBarry Smith while (next->next) { 8370971522cSBarry Smith next = next->next; 8380971522cSBarry Smith } 8395a9f2f41SSatish Balay next->next = ilink; 84051f519a2SBarry Smith ilink->previous = next; 8410971522cSBarry Smith } 8420971522cSBarry Smith jac->nsplits++; 8430971522cSBarry Smith PetscFunctionReturn(0); 8440971522cSBarry Smith } 8450971522cSBarry Smith EXTERN_C_END 8460971522cSBarry Smith 847e69d4d44SBarry Smith EXTERN_C_BEGIN 848e69d4d44SBarry Smith #undef __FUNCT__ 849e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 8507087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 851e69d4d44SBarry Smith { 852e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 853e69d4d44SBarry Smith PetscErrorCode ierr; 854e69d4d44SBarry Smith 855e69d4d44SBarry Smith PetscFunctionBegin; 856519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 857e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 858e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 859084e4875SJed Brown *n = jac->nsplits; 860e69d4d44SBarry Smith PetscFunctionReturn(0); 861e69d4d44SBarry Smith } 862e69d4d44SBarry Smith EXTERN_C_END 8630971522cSBarry Smith 8640971522cSBarry Smith EXTERN_C_BEGIN 8650971522cSBarry Smith #undef __FUNCT__ 86669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 8677087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 8680971522cSBarry Smith { 8690971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8700971522cSBarry Smith PetscErrorCode ierr; 8710971522cSBarry Smith PetscInt cnt = 0; 8725a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 8730971522cSBarry Smith 8740971522cSBarry Smith PetscFunctionBegin; 87569a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 8765a9f2f41SSatish Balay while (ilink) { 8775a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 8785a9f2f41SSatish Balay ilink = ilink->next; 8790971522cSBarry Smith } 880e32f2f54SBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 8810971522cSBarry Smith *n = jac->nsplits; 8820971522cSBarry Smith PetscFunctionReturn(0); 8830971522cSBarry Smith } 8840971522cSBarry Smith EXTERN_C_END 8850971522cSBarry Smith 886704ba839SBarry Smith EXTERN_C_BEGIN 887704ba839SBarry Smith #undef __FUNCT__ 888704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 8897087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 890704ba839SBarry Smith { 891704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 892704ba839SBarry Smith PetscErrorCode ierr; 893704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 894704ba839SBarry Smith char prefix[128]; 895704ba839SBarry Smith 896704ba839SBarry Smith PetscFunctionBegin; 8976c924f48SJed Brown if (jac->splitdefined) { 8986c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 8996c924f48SJed Brown PetscFunctionReturn(0); 9006c924f48SJed Brown } 90116913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 902db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 9031ab39975SBarry Smith ilink->is = is; 9041ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 905704ba839SBarry Smith ilink->next = PETSC_NULL; 906704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 9071cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 908704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 9099005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 910704ba839SBarry Smith 9116c924f48SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr); 912704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 913704ba839SBarry Smith 914704ba839SBarry Smith if (!next) { 915704ba839SBarry Smith jac->head = ilink; 916704ba839SBarry Smith ilink->previous = PETSC_NULL; 917704ba839SBarry Smith } else { 918704ba839SBarry Smith while (next->next) { 919704ba839SBarry Smith next = next->next; 920704ba839SBarry Smith } 921704ba839SBarry Smith next->next = ilink; 922704ba839SBarry Smith ilink->previous = next; 923704ba839SBarry Smith } 924704ba839SBarry Smith jac->nsplits++; 925704ba839SBarry Smith 926704ba839SBarry Smith PetscFunctionReturn(0); 927704ba839SBarry Smith } 928704ba839SBarry Smith EXTERN_C_END 929704ba839SBarry Smith 9300971522cSBarry Smith #undef __FUNCT__ 9310971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 9320971522cSBarry Smith /*@ 9330971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 9340971522cSBarry Smith 935ad4df100SBarry Smith Logically Collective on PC 9360971522cSBarry Smith 9370971522cSBarry Smith Input Parameters: 9380971522cSBarry Smith + pc - the preconditioner context 939db4c96c1SJed Brown . splitname - name of this split 9400971522cSBarry Smith . n - the number of fields in this split 941db4c96c1SJed Brown - fields - the fields in this split 9420971522cSBarry Smith 9430971522cSBarry Smith Level: intermediate 9440971522cSBarry Smith 945d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 946d32f9abdSBarry Smith 947d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 948d32f9abdSBarry 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 949d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 950d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 951d32f9abdSBarry Smith 952db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 953db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 954db4c96c1SJed Brown 955d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 9560971522cSBarry Smith 9570971522cSBarry Smith @*/ 9587087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 9590971522cSBarry Smith { 9604ac538c5SBarry Smith PetscErrorCode ierr; 9610971522cSBarry Smith 9620971522cSBarry Smith PetscFunctionBegin; 9630700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 964db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 965db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 966db4c96c1SJed Brown PetscValidIntPointer(fields,3); 9674ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr); 9680971522cSBarry Smith PetscFunctionReturn(0); 9690971522cSBarry Smith } 9700971522cSBarry Smith 9710971522cSBarry Smith #undef __FUNCT__ 972704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 973704ba839SBarry Smith /*@ 974704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 975704ba839SBarry Smith 976ad4df100SBarry Smith Logically Collective on PC 977704ba839SBarry Smith 978704ba839SBarry Smith Input Parameters: 979704ba839SBarry Smith + pc - the preconditioner context 980db4c96c1SJed Brown . splitname - name of this split 981db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 982704ba839SBarry Smith 983d32f9abdSBarry Smith 984a6ffb8dbSJed Brown Notes: 985a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 986a6ffb8dbSJed Brown 987db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 988db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 989d32f9abdSBarry Smith 990704ba839SBarry Smith Level: intermediate 991704ba839SBarry Smith 992704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 993704ba839SBarry Smith 994704ba839SBarry Smith @*/ 9957087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 996704ba839SBarry Smith { 9974ac538c5SBarry Smith PetscErrorCode ierr; 998704ba839SBarry Smith 999704ba839SBarry Smith PetscFunctionBegin; 10000700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1001db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1002db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 10034ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1004704ba839SBarry Smith PetscFunctionReturn(0); 1005704ba839SBarry Smith } 1006704ba839SBarry Smith 1007704ba839SBarry Smith #undef __FUNCT__ 100851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 100951f519a2SBarry Smith /*@ 101051f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 101151f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 101251f519a2SBarry Smith 1013ad4df100SBarry Smith Logically Collective on PC 101451f519a2SBarry Smith 101551f519a2SBarry Smith Input Parameters: 101651f519a2SBarry Smith + pc - the preconditioner context 101751f519a2SBarry Smith - bs - the block size 101851f519a2SBarry Smith 101951f519a2SBarry Smith Level: intermediate 102051f519a2SBarry Smith 102151f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 102251f519a2SBarry Smith 102351f519a2SBarry Smith @*/ 10247087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 102551f519a2SBarry Smith { 10264ac538c5SBarry Smith PetscErrorCode ierr; 102751f519a2SBarry Smith 102851f519a2SBarry Smith PetscFunctionBegin; 10290700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1030c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 10314ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 103251f519a2SBarry Smith PetscFunctionReturn(0); 103351f519a2SBarry Smith } 103451f519a2SBarry Smith 103551f519a2SBarry Smith #undef __FUNCT__ 103669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 10370971522cSBarry Smith /*@C 103869a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 10390971522cSBarry Smith 104069a612a9SBarry Smith Collective on KSP 10410971522cSBarry Smith 10420971522cSBarry Smith Input Parameter: 10430971522cSBarry Smith . pc - the preconditioner context 10440971522cSBarry Smith 10450971522cSBarry Smith Output Parameters: 10460971522cSBarry Smith + n - the number of split 104769a612a9SBarry Smith - pc - the array of KSP contexts 10480971522cSBarry Smith 10490971522cSBarry Smith Note: 1050d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1051d32f9abdSBarry Smith (not the KSP just the array that contains them). 10520971522cSBarry Smith 105369a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 10540971522cSBarry Smith 10550971522cSBarry Smith Level: advanced 10560971522cSBarry Smith 10570971522cSBarry Smith .seealso: PCFIELDSPLIT 10580971522cSBarry Smith @*/ 10597087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 10600971522cSBarry Smith { 10614ac538c5SBarry Smith PetscErrorCode ierr; 10620971522cSBarry Smith 10630971522cSBarry Smith PetscFunctionBegin; 10640700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10650971522cSBarry Smith PetscValidIntPointer(n,2); 10664ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 10670971522cSBarry Smith PetscFunctionReturn(0); 10680971522cSBarry Smith } 10690971522cSBarry Smith 1070e69d4d44SBarry Smith #undef __FUNCT__ 1071e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1072e69d4d44SBarry Smith /*@ 1073e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1074e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 1075e69d4d44SBarry Smith 1076e69d4d44SBarry Smith Collective on PC 1077e69d4d44SBarry Smith 1078e69d4d44SBarry Smith Input Parameters: 1079e69d4d44SBarry Smith + pc - the preconditioner context 1080084e4875SJed Brown . ptype - which matrix to use for preconditioning the Schur complement 1081084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1082084e4875SJed Brown 1083084e4875SJed Brown Notes: 1084084e4875SJed Brown The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 1085084e4875SJed Brown There are currently no preconditioners that work directly with the Schur complement so setting 1086084e4875SJed Brown PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 1087e69d4d44SBarry Smith 1088e69d4d44SBarry Smith Options Database: 1089084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1090e69d4d44SBarry Smith 1091e69d4d44SBarry Smith Level: intermediate 1092e69d4d44SBarry Smith 1093084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1094e69d4d44SBarry Smith 1095e69d4d44SBarry Smith @*/ 10967087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1097e69d4d44SBarry Smith { 10984ac538c5SBarry Smith PetscErrorCode ierr; 1099e69d4d44SBarry Smith 1100e69d4d44SBarry Smith PetscFunctionBegin; 11010700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11024ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1103e69d4d44SBarry Smith PetscFunctionReturn(0); 1104e69d4d44SBarry Smith } 1105e69d4d44SBarry Smith 1106e69d4d44SBarry Smith EXTERN_C_BEGIN 1107e69d4d44SBarry Smith #undef __FUNCT__ 1108e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 11097087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1110e69d4d44SBarry Smith { 1111e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1112084e4875SJed Brown PetscErrorCode ierr; 1113e69d4d44SBarry Smith 1114e69d4d44SBarry Smith PetscFunctionBegin; 1115084e4875SJed Brown jac->schurpre = ptype; 1116084e4875SJed Brown if (pre) { 1117084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1118084e4875SJed Brown jac->schur_user = pre; 1119084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1120084e4875SJed Brown } 1121e69d4d44SBarry Smith PetscFunctionReturn(0); 1122e69d4d44SBarry Smith } 1123e69d4d44SBarry Smith EXTERN_C_END 1124e69d4d44SBarry Smith 112530ad9308SMatthew Knepley #undef __FUNCT__ 112630ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 112730ad9308SMatthew Knepley /*@C 11288c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 112930ad9308SMatthew Knepley 113030ad9308SMatthew Knepley Collective on KSP 113130ad9308SMatthew Knepley 113230ad9308SMatthew Knepley Input Parameter: 113330ad9308SMatthew Knepley . pc - the preconditioner context 113430ad9308SMatthew Knepley 113530ad9308SMatthew Knepley Output Parameters: 113630ad9308SMatthew Knepley + A - the (0,0) block 113730ad9308SMatthew Knepley . B - the (0,1) block 113830ad9308SMatthew Knepley . C - the (1,0) block 113930ad9308SMatthew Knepley - D - the (1,1) block 114030ad9308SMatthew Knepley 114130ad9308SMatthew Knepley Level: advanced 114230ad9308SMatthew Knepley 114330ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 114430ad9308SMatthew Knepley @*/ 11457087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 114630ad9308SMatthew Knepley { 114730ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 114830ad9308SMatthew Knepley 114930ad9308SMatthew Knepley PetscFunctionBegin; 11500700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1151c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 115230ad9308SMatthew Knepley if (A) *A = jac->pmat[0]; 115330ad9308SMatthew Knepley if (B) *B = jac->B; 115430ad9308SMatthew Knepley if (C) *C = jac->C; 115530ad9308SMatthew Knepley if (D) *D = jac->pmat[1]; 115630ad9308SMatthew Knepley PetscFunctionReturn(0); 115730ad9308SMatthew Knepley } 115830ad9308SMatthew Knepley 115979416396SBarry Smith EXTERN_C_BEGIN 116079416396SBarry Smith #undef __FUNCT__ 116179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 11627087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 116379416396SBarry Smith { 116479416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1165e69d4d44SBarry Smith PetscErrorCode ierr; 116679416396SBarry Smith 116779416396SBarry Smith PetscFunctionBegin; 116879416396SBarry Smith jac->type = type; 11693b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 11703b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 11713b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1172e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1173e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1174e69d4d44SBarry Smith 11753b224e63SBarry Smith } else { 11763b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 11773b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1178e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 11799e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 11803b224e63SBarry Smith } 118179416396SBarry Smith PetscFunctionReturn(0); 118279416396SBarry Smith } 118379416396SBarry Smith EXTERN_C_END 118479416396SBarry Smith 118551f519a2SBarry Smith EXTERN_C_BEGIN 118651f519a2SBarry Smith #undef __FUNCT__ 118751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 11887087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 118951f519a2SBarry Smith { 119051f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 119151f519a2SBarry Smith 119251f519a2SBarry Smith PetscFunctionBegin; 119365e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 119465e19b50SBarry 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); 119551f519a2SBarry Smith jac->bs = bs; 119651f519a2SBarry Smith PetscFunctionReturn(0); 119751f519a2SBarry Smith } 119851f519a2SBarry Smith EXTERN_C_END 119951f519a2SBarry Smith 120079416396SBarry Smith #undef __FUNCT__ 120179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1202bc08b0f1SBarry Smith /*@ 120379416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 120479416396SBarry Smith 120579416396SBarry Smith Collective on PC 120679416396SBarry Smith 120779416396SBarry Smith Input Parameter: 120879416396SBarry Smith . pc - the preconditioner context 120981540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 121079416396SBarry Smith 121179416396SBarry Smith Options Database Key: 1212a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 121379416396SBarry Smith 121479416396SBarry Smith Level: Developer 121579416396SBarry Smith 121679416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 121779416396SBarry Smith 121879416396SBarry Smith .seealso: PCCompositeSetType() 121979416396SBarry Smith 122079416396SBarry Smith @*/ 12217087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 122279416396SBarry Smith { 12234ac538c5SBarry Smith PetscErrorCode ierr; 122479416396SBarry Smith 122579416396SBarry Smith PetscFunctionBegin; 12260700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12274ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 122879416396SBarry Smith PetscFunctionReturn(0); 122979416396SBarry Smith } 123079416396SBarry Smith 12310971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 12320971522cSBarry Smith /*MC 1233a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 12340971522cSBarry Smith fields or groups of fields 12350971522cSBarry Smith 12360971522cSBarry Smith 1237edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1238edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 12390971522cSBarry Smith 1240a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 124169a612a9SBarry Smith and set the options directly on the resulting KSP object 12420971522cSBarry Smith 12430971522cSBarry Smith Level: intermediate 12440971522cSBarry Smith 124579416396SBarry Smith Options Database Keys: 124681540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 124781540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 124881540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 124981540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 125081540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1251e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 125279416396SBarry Smith 1253edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 12543b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 12553b224e63SBarry Smith 12563b224e63SBarry Smith 1257d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1258d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1259d32f9abdSBarry Smith 1260d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1261d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1262d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1263d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1264d32f9abdSBarry Smith 1265d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1266d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1267d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1268d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1269d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1270d32f9abdSBarry Smith 1271e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1272e69d4d44SBarry Smith ( C D ) 1273e69d4d44SBarry Smith the preconditioner is 1274e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1275e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1276edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1277e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1278edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1279e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1280e69d4d44SBarry Smith option. 1281e69d4d44SBarry Smith 1282edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1283edf189efSBarry Smith is used automatically for a second block. 1284edf189efSBarry Smith 1285a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 12860971522cSBarry Smith 1287a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1288e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 12890971522cSBarry Smith M*/ 12900971522cSBarry Smith 12910971522cSBarry Smith EXTERN_C_BEGIN 12920971522cSBarry Smith #undef __FUNCT__ 12930971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 12947087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 12950971522cSBarry Smith { 12960971522cSBarry Smith PetscErrorCode ierr; 12970971522cSBarry Smith PC_FieldSplit *jac; 12980971522cSBarry Smith 12990971522cSBarry Smith PetscFunctionBegin; 130038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 13010971522cSBarry Smith jac->bs = -1; 13020971522cSBarry Smith jac->nsplits = 0; 13033e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1304e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1305c5d2311dSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL; 130651f519a2SBarry Smith 13070971522cSBarry Smith pc->data = (void*)jac; 13080971522cSBarry Smith 13090971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1310421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 13110971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 13120971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 13130971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 13140971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 13150971522cSBarry Smith pc->ops->applyrichardson = 0; 13160971522cSBarry Smith 131769a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 131869a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 13190971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 13200971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1321704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1322704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 132379416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 132479416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 132551f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 132651f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 13270971522cSBarry Smith PetscFunctionReturn(0); 13280971522cSBarry Smith } 13290971522cSBarry Smith EXTERN_C_END 13300971522cSBarry Smith 13310971522cSBarry Smith 1332a541d17aSBarry Smith /*MC 1333a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1334a541d17aSBarry Smith overview of these methods. 1335a541d17aSBarry Smith 1336a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1337a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1338a541d17aSBarry Smith 1339a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1340a541d17aSBarry Smith B' 0) (x_2) (b_2) 1341a541d17aSBarry Smith 1342a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1343a541d17aSBarry Smith for this block system. 1344a541d17aSBarry Smith 1345a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1346a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1347a541d17aSBarry Smith in the original matrix (for example Ap == A). 1348a541d17aSBarry Smith 1349a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1350a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1351a541d17aSBarry Smith 1352a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1353a541d17aSBarry Smith x_2 = D^ b_2 1354a541d17aSBarry Smith 1355a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1356a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1357a541d17aSBarry Smith 1358a541d17aSBarry Smith Symmetric Gauss-Seidel: x_1 = x_1 + A^(b_1 - A x_1 - B x_2) variant x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2) 1359a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1360a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1361a541d17aSBarry Smith 1362*6ce1633cSBarry Smith Schur complement preconditioner 1363*6ce1633cSBarry Smith the preconditioner is 1364*6ce1633cSBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1365*6ce1633cSBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1366*6ce1633cSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1367*6ce1633cSBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1368*6ce1633cSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1369*6ce1633cSBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1370*6ce1633cSBarry Smith option. 1371*6ce1633cSBarry Smith 13720bc0a719SBarry Smith Level: intermediate 13730bc0a719SBarry Smith 1374a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1375a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1376a541d17aSBarry Smith M*/ 1377