1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 31e25c274SJed Brown #include <petscdm.h> 40971522cSBarry Smith 5443836d0SMatthew G Knepley /* 6443836d0SMatthew G Knepley There is a nice discussion of block preconditioners in 7443836d0SMatthew G Knepley 8443836d0SMatthew G Knepley [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier–Stokes equations 9443836d0SMatthew G Knepley Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808 101b8e4c5fSJed Brown http://chess.cs.umd.edu/~elman/papers/tax.pdf 11443836d0SMatthew G Knepley */ 12443836d0SMatthew G Knepley 13e87fab1bSBarry Smith const char *const PCFieldSplitSchurPreTypes[] = {"SELF","A11","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 14c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0}; 15c5d2311dSJed Brown 160971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 170971522cSBarry Smith struct _PC_FieldSplitLink { 1869a612a9SBarry Smith KSP ksp; 19443836d0SMatthew G Knepley Vec x,y,z; 20db4c96c1SJed Brown char *splitname; 210971522cSBarry Smith PetscInt nfields; 225d4c12cdSJungho Lee PetscInt *fields,*fields_col; 231b9fc7fcSBarry Smith VecScatter sctx; 245d4c12cdSJungho Lee IS is,is_col; 2551f519a2SBarry Smith PC_FieldSplitLink next,previous; 260971522cSBarry Smith }; 270971522cSBarry Smith 280971522cSBarry Smith typedef struct { 2968dd23aaSBarry Smith PCCompositeType type; 30ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 31ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 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; 392fa5cd67SKarl Rupp 4030ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 4130ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 4230ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 43443836d0SMatthew G Knepley Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */ 44084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 45084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 46c9c6ffaaSJed Brown PCFieldSplitSchurFactType schurfactorization; 4730ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 48443836d0SMatthew G Knepley KSP kspupper; /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */ 4997bbdb24SBarry Smith PC_FieldSplitLink head; 5063ec74ffSBarry Smith PetscBool reset; /* indicates PCReset() has been last called on this object, hack */ 51c1570756SJed Brown PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */ 524ab8060aSDmitry Karpeev PetscBool dm_splits; /* Whether to use DMCreateFieldDecomposition() whenever possible */ 530971522cSBarry Smith } PC_FieldSplit; 540971522cSBarry Smith 5516913363SBarry Smith /* 5616913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5716913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5816913363SBarry Smith PC you could change this. 5916913363SBarry Smith */ 60084e4875SJed Brown 61e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 62084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 63084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 64084e4875SJed Brown { 65084e4875SJed Brown switch (jac->schurpre) { 66084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 67e87fab1bSBarry Smith case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1]; 68084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 69084e4875SJed Brown default: 70084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 71084e4875SJed Brown } 72084e4875SJed Brown } 73084e4875SJed Brown 74084e4875SJed Brown 759804daf3SBarry Smith #include <petscdraw.h> 760971522cSBarry Smith #undef __FUNCT__ 770971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 780971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 790971522cSBarry Smith { 800971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 810971522cSBarry Smith PetscErrorCode ierr; 82d9884438SBarry Smith PetscBool iascii,isdraw; 830971522cSBarry Smith PetscInt i,j; 845a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 850971522cSBarry Smith 860971522cSBarry Smith PetscFunctionBegin; 87251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 88d9884438SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 890971522cSBarry Smith if (iascii) { 902c9966d7SBarry Smith if (jac->bs > 0) { 9151f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 922c9966d7SBarry Smith } else { 932c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 942c9966d7SBarry Smith } 95f5236f50SJed Brown if (pc->useAmat) { 96f5236f50SJed Brown ierr = PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr); 97a3df900dSMatthew G Knepley } 9869a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 990971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1000971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 1011ab39975SBarry Smith if (ilink->fields) { 1020971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 10379416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1045a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 10579416396SBarry Smith if (j > 0) { 10679416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 10779416396SBarry Smith } 1085a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1090971522cSBarry Smith } 1100971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 11179416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1121ab39975SBarry Smith } else { 1131ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1141ab39975SBarry Smith } 1155a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1165a9f2f41SSatish Balay ilink = ilink->next; 1170971522cSBarry Smith } 1180971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1192fa5cd67SKarl Rupp } 1202fa5cd67SKarl Rupp 1212fa5cd67SKarl Rupp if (isdraw) { 122d9884438SBarry Smith PetscDraw draw; 123d9884438SBarry Smith PetscReal x,y,w,wd; 124d9884438SBarry Smith 125d9884438SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 126d9884438SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 127d9884438SBarry Smith w = 2*PetscMin(1.0 - x,x); 128d9884438SBarry Smith wd = w/(jac->nsplits + 1); 129d9884438SBarry Smith x = x - wd*(jac->nsplits-1)/2.0; 130d9884438SBarry Smith for (i=0; i<jac->nsplits; i++) { 131d9884438SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 132d9884438SBarry Smith ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 133d9884438SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 134d9884438SBarry Smith x += wd; 135d9884438SBarry Smith ilink = ilink->next; 136d9884438SBarry Smith } 1370971522cSBarry Smith } 1380971522cSBarry Smith PetscFunctionReturn(0); 1390971522cSBarry Smith } 1400971522cSBarry Smith 1410971522cSBarry Smith #undef __FUNCT__ 1423b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1433b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1443b224e63SBarry Smith { 1453b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1463b224e63SBarry Smith PetscErrorCode ierr; 1474996c5bdSBarry Smith PetscBool iascii,isdraw; 1483b224e63SBarry Smith PetscInt i,j; 1493b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1503b224e63SBarry Smith 1513b224e63SBarry Smith PetscFunctionBegin; 152251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1534996c5bdSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 1543b224e63SBarry Smith if (iascii) { 1552c9966d7SBarry Smith if (jac->bs > 0) { 156c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1572c9966d7SBarry Smith } else { 1582c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1592c9966d7SBarry Smith } 160f5236f50SJed Brown if (pc->useAmat) { 161f5236f50SJed Brown ierr = PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr); 162a3df900dSMatthew G Knepley } 1633e8b8b31SMatthew G Knepley switch (jac->schurpre) { 1643e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1653e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 166e87fab1bSBarry Smith case PC_FIELDSPLIT_SCHUR_PRE_A11: 167e87fab1bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break; 1683e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1693e8b8b31SMatthew G Knepley if (jac->schur_user) { 1703e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1713e8b8b31SMatthew G Knepley } else { 172e87fab1bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr); 1733e8b8b31SMatthew G Knepley } 1743e8b8b31SMatthew G Knepley break; 1753e8b8b31SMatthew G Knepley default: 17682f516ccSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1773e8b8b31SMatthew G Knepley } 1783b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1793b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1803b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1813b224e63SBarry Smith if (ilink->fields) { 1823b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1833b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1843b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1853b224e63SBarry Smith if (j > 0) { 1863b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1873b224e63SBarry Smith } 1883b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1893b224e63SBarry Smith } 1903b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1913b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1923b224e63SBarry Smith } else { 1933b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1943b224e63SBarry Smith } 1953b224e63SBarry Smith ilink = ilink->next; 1963b224e63SBarry Smith } 197435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1983b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 199443836d0SMatthew G Knepley ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 2003b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 201443836d0SMatthew G Knepley if (jac->kspupper != jac->head->ksp) { 202443836d0SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr); 203443836d0SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 204443836d0SMatthew G Knepley ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 205443836d0SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 206443836d0SMatthew G Knepley } 207435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 2083b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 20912cae6f2SJed Brown if (jac->kspschur) { 2103b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 21112cae6f2SJed Brown } else { 21212cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 21312cae6f2SJed Brown } 2143b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2153b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2164996c5bdSBarry Smith } else if (isdraw) { 2174996c5bdSBarry Smith PetscDraw draw; 2184996c5bdSBarry Smith PetscReal x,y,w,wd,h; 2194996c5bdSBarry Smith PetscInt cnt = 2; 2204996c5bdSBarry Smith char str[32]; 2214996c5bdSBarry Smith 2224996c5bdSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 2234996c5bdSBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 224c74581afSBarry Smith if (jac->kspupper != jac->head->ksp) cnt++; 225c74581afSBarry Smith w = 2*PetscMin(1.0 - x,x); 226c74581afSBarry Smith wd = w/(cnt + 1); 227c74581afSBarry Smith 2284996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 2290298fd71SBarry Smith ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 2304996c5bdSBarry Smith y -= h; 2314996c5bdSBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) { 232e87fab1bSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr); 2333b224e63SBarry Smith } else { 2344996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr); 2354996c5bdSBarry Smith } 2360298fd71SBarry Smith ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 2374996c5bdSBarry Smith y -= h; 2384996c5bdSBarry Smith x = x - wd*(cnt-1)/2.0; 2394996c5bdSBarry Smith 2404996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2414996c5bdSBarry Smith ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 2424996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2434996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) { 2444996c5bdSBarry Smith x += wd; 2454996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2464996c5bdSBarry Smith ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 2474996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2484996c5bdSBarry Smith } 2494996c5bdSBarry Smith x += wd; 2504996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2514996c5bdSBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 2524996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2533b224e63SBarry Smith } 2543b224e63SBarry Smith PetscFunctionReturn(0); 2553b224e63SBarry Smith } 2563b224e63SBarry Smith 2573b224e63SBarry Smith #undef __FUNCT__ 2586c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2596c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2606c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2616c924f48SJed Brown { 2626c924f48SJed Brown PetscErrorCode ierr; 2636c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2645d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2655d4c12cdSJungho Lee PetscBool flg,flg_col; 2665d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2676c924f48SJed Brown 2686c924f48SJed Brown PetscFunctionBegin; 2696c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2705d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2716c924f48SJed Brown for (i=0,flg=PETSC_TRUE;; i++) { 2728caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2738caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2748caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2756c924f48SJed Brown nfields = jac->bs; 27629499fbbSJungho Lee nfields_col = jac->bs; 2776c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2785d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2796c924f48SJed Brown if (!flg) break; 2805d4c12cdSJungho Lee else if (flg && !flg_col) { 2815d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2825d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2832fa5cd67SKarl Rupp } else { 2845d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2855d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2865d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2875d4c12cdSJungho Lee } 2886c924f48SJed Brown } 2896c924f48SJed Brown if (i > 0) { 2906c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2916c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2926c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2936c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2946c924f48SJed Brown } 2956c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2965d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2976c924f48SJed Brown PetscFunctionReturn(0); 2986c924f48SJed Brown } 2996c924f48SJed Brown 3006c924f48SJed Brown #undef __FUNCT__ 30169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 30269a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 3030971522cSBarry Smith { 3040971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3050971522cSBarry Smith PetscErrorCode ierr; 3065a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3077287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 3086c924f48SJed Brown PetscInt i; 3090971522cSBarry Smith 3100971522cSBarry Smith PetscFunctionBegin; 3117287d2a3SDmitry Karpeev /* 3127287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 3137287d2a3SDmitry Karpeev Should probably be rewritten. 3147287d2a3SDmitry Karpeev */ 3157287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 3160298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 3174ab8060aSDmitry Karpeev if (pc->dm && jac->dm_splits && !stokes) { 318bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 3190784a22cSJed Brown char **fieldNames; 3207b62db95SJungho Lee IS *fields; 321e7c4fc90SDmitry Karpeev DM *dms; 322bafc1b83SMatthew G Knepley DM subdm[128]; 323bafc1b83SMatthew G Knepley PetscBool flg; 324bafc1b83SMatthew G Knepley 32516621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 326bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 327bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE;; i++) { 328bafc1b83SMatthew G Knepley PetscInt ifields[128]; 329bafc1b83SMatthew G Knepley IS compField; 330bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 331bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 332bafc1b83SMatthew G Knepley 3338caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 334bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 335bafc1b83SMatthew G Knepley if (!flg) break; 33682f516ccSBarry Smith if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 337bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 338bafc1b83SMatthew G Knepley if (nfields == 1) { 339bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 34082f516ccSBarry Smith /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 3410298fd71SBarry Smith ierr = ISView(compField, NULL);CHKERRQ(ierr); */ 342bafc1b83SMatthew G Knepley } else { 3438caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 344bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 34582f516ccSBarry Smith /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", splitname);CHKERRQ(ierr); 3460298fd71SBarry Smith ierr = ISView(compField, NULL);CHKERRQ(ierr); */ 3477287d2a3SDmitry Karpeev } 348bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 349bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 350bafc1b83SMatthew G Knepley f = ifields[j]; 3517b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 3527b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 3537b62db95SJungho Lee } 354bafc1b83SMatthew G Knepley } 355bafc1b83SMatthew G Knepley if (i == 0) { 356bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 357bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 358bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 359bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 360bafc1b83SMatthew G Knepley } 361bafc1b83SMatthew G Knepley } else { 362d724dfffSBarry Smith for (j=0; j<numFields; j++) { 363d724dfffSBarry Smith ierr = DMDestroy(dms+j);CHKERRQ(ierr); 364d724dfffSBarry Smith } 365d724dfffSBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 366bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 3672fa5cd67SKarl Rupp for (j = 0; j < i; ++j) dms[j] = subdm[j]; 368bafc1b83SMatthew G Knepley } 3697b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3707b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 371e7c4fc90SDmitry Karpeev if (dms) { 3728b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 373bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3747287d2a3SDmitry Karpeev const char *prefix; 3757287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr); 3767287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr); 3777b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3787b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3797287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr); 380e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 3812fa5ba8aSJed Brown } 3827b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3838b8307b2SJed Brown } 38466ffff09SJed Brown } else { 385521d7252SBarry Smith if (jac->bs <= 0) { 386704ba839SBarry Smith if (pc->pmat) { 387521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 3882fa5cd67SKarl Rupp } else jac->bs = 1; 389521d7252SBarry Smith } 390d32f9abdSBarry Smith 3916ce1633cSBarry Smith if (stokes) { 3926ce1633cSBarry Smith IS zerodiags,rest; 3936ce1633cSBarry Smith PetscInt nmin,nmax; 3946ce1633cSBarry Smith 3956ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3966ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3976ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 3987287d2a3SDmitry Karpeev if (jac->reset) { 3997287d2a3SDmitry Karpeev jac->head->is = rest; 4007287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 4012fa5cd67SKarl Rupp } else { 4026ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 4036ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 4047287d2a3SDmitry Karpeev } 405fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 406fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 4076ce1633cSBarry Smith } else { 408ce94432eSBarry Smith if (jac->reset) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 4099eeaaa73SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr); 4107287d2a3SDmitry Karpeev if (!fieldsplit_default) { 411d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 412d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 4136c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 4146c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 415d32f9abdSBarry Smith } 4167287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 417d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 418db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 4196c924f48SJed Brown char splitname[8]; 4208caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 4215d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 42279416396SBarry Smith } 4235d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 424521d7252SBarry Smith } 42566ffff09SJed Brown } 4266ce1633cSBarry Smith } 427edf189efSBarry Smith } else if (jac->nsplits == 1) { 428edf189efSBarry Smith if (ilink->is) { 429edf189efSBarry Smith IS is2; 430edf189efSBarry Smith PetscInt nmin,nmax; 431edf189efSBarry Smith 432edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 433edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 434db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 435fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 43682f516ccSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 437edf189efSBarry Smith } 438d0af7cd3SBarry Smith 439d0af7cd3SBarry Smith 44082f516ccSBarry Smith if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 44169a612a9SBarry Smith PetscFunctionReturn(0); 44269a612a9SBarry Smith } 44369a612a9SBarry Smith 4445a576424SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 445514bf10dSMatthew G Knepley 44669a612a9SBarry Smith #undef __FUNCT__ 44769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 44869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 44969a612a9SBarry Smith { 45069a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 45169a612a9SBarry Smith PetscErrorCode ierr; 4525a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4532c9966d7SBarry Smith PetscInt i,nsplit; 45469a612a9SBarry Smith MatStructure flag = pc->flag; 4555f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 45669a612a9SBarry Smith 45769a612a9SBarry Smith PetscFunctionBegin; 45869a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 45997bbdb24SBarry Smith nsplit = jac->nsplits; 4605a9f2f41SSatish Balay ilink = jac->head; 46197bbdb24SBarry Smith 46297bbdb24SBarry Smith /* get the matrices for each split */ 463704ba839SBarry Smith if (!jac->issetup) { 4641b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 46597bbdb24SBarry Smith 466704ba839SBarry Smith jac->issetup = PETSC_TRUE; 467704ba839SBarry Smith 4685d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4692c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4702c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4712c9966d7SBarry Smith } 47251f519a2SBarry Smith bs = jac->bs; 47397bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4741b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4751b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4761b9fc7fcSBarry Smith if (jac->defaultsplit) { 477ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4785f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 479704ba839SBarry Smith } else if (!ilink->is) { 480ccb205f8SBarry Smith if (ilink->nfields > 1) { 4815f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4825a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4835f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4841b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4851b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4861b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4875f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 48897bbdb24SBarry Smith } 48997bbdb24SBarry Smith } 490ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 491ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 492ccb205f8SBarry Smith } else { 493ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 494ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 495ccb205f8SBarry Smith } 4963e197d65SBarry Smith } 497704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4985f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4995f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 500704ba839SBarry Smith ilink = ilink->next; 5011b9fc7fcSBarry Smith } 5021b9fc7fcSBarry Smith } 5031b9fc7fcSBarry Smith 504704ba839SBarry Smith ilink = jac->head; 50597bbdb24SBarry Smith if (!jac->pmat) { 506bdddcaaaSMatthew G Knepley Vec xtmp; 507bdddcaaaSMatthew G Knepley 5080298fd71SBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr); 509cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 510bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 511cf502942SBarry Smith for (i=0; i<nsplit; i++) { 512bdddcaaaSMatthew G Knepley MatNullSpace sp; 513bdddcaaaSMatthew G Knepley 514a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 515a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr); 516a3df900dSMatthew G Knepley if (jac->pmat[i]) { 517a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 518a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 519a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 5202fa5cd67SKarl Rupp 521a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 522a3df900dSMatthew G Knepley } 523a3df900dSMatthew G Knepley } else { 5245f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 525a3df900dSMatthew G Knepley } 526bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 527bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 5280298fd71SBarry Smith ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL; 529bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5300298fd71SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr); 531ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 532bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr); 533bafc1b83SMatthew G Knepley if (sp) { 534bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 535bafc1b83SMatthew G Knepley } 536ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr); 537ed1f0337SMatthew G Knepley if (sp) { 538ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 539ed1f0337SMatthew G Knepley } 540704ba839SBarry Smith ilink = ilink->next; 541cf502942SBarry Smith } 542bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 54397bbdb24SBarry Smith } else { 544cf502942SBarry Smith for (i=0; i<nsplit; i++) { 545a3df900dSMatthew G Knepley Mat pmat; 546a3df900dSMatthew G Knepley 547a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 548a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr); 549a3df900dSMatthew G Knepley if (!pmat) { 5505f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 551a3df900dSMatthew G Knepley } 552704ba839SBarry Smith ilink = ilink->next; 553cf502942SBarry Smith } 55497bbdb24SBarry Smith } 555f5236f50SJed Brown if (pc->useAmat) { 556519d70e2SJed Brown ilink = jac->head; 557519d70e2SJed Brown if (!jac->mat) { 558519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 559519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5605f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 561519d70e2SJed Brown ilink = ilink->next; 562519d70e2SJed Brown } 563519d70e2SJed Brown } else { 564519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5655f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 566519d70e2SJed Brown ilink = ilink->next; 567519d70e2SJed Brown } 568519d70e2SJed Brown } 569519d70e2SJed Brown } else { 570519d70e2SJed Brown jac->mat = jac->pmat; 571519d70e2SJed Brown } 57297bbdb24SBarry Smith 5736c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 57468dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 57568dd23aaSBarry Smith ilink = jac->head; 57668dd23aaSBarry Smith if (!jac->Afield) { 57768dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 57868dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5790298fd71SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58068dd23aaSBarry Smith ilink = ilink->next; 58168dd23aaSBarry Smith } 58268dd23aaSBarry Smith } else { 58368dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5840298fd71SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58568dd23aaSBarry Smith ilink = ilink->next; 58668dd23aaSBarry Smith } 58768dd23aaSBarry Smith } 58868dd23aaSBarry Smith } 58968dd23aaSBarry Smith 5903b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5913b224e63SBarry Smith IS ccis; 5924aa3045dSJed Brown PetscInt rstart,rend; 593093c86ffSJed Brown char lscname[256]; 594093c86ffSJed Brown PetscObject LSC_L; 595ce94432eSBarry Smith 596ce94432eSBarry Smith if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 59768dd23aaSBarry Smith 598e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 599e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 600e6cab6aaSJed Brown 6013b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 6023b224e63SBarry Smith if (jac->schur) { 6030298fd71SBarry Smith KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper; 604443836d0SMatthew G Knepley 605fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 6063b224e63SBarry Smith ilink = jac->head; 60749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6084aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 609fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6103b224e63SBarry Smith ilink = ilink->next; 61149bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6124aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 613fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 614a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 615443836d0SMatthew G Knepley if (kspA != kspInner) { 616443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 617443836d0SMatthew G Knepley } 618443836d0SMatthew G Knepley if (kspUpper != kspA) { 619443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 620443836d0SMatthew G Knepley } 621084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 6223b224e63SBarry Smith } else { 623bafc1b83SMatthew G Knepley const char *Dprefix; 624bafc1b83SMatthew G Knepley char schurprefix[256]; 625514bf10dSMatthew G Knepley char schurtestoption[256]; 626bdddcaaaSMatthew G Knepley MatNullSpace sp; 627514bf10dSMatthew G Knepley PetscBool flg; 6283b224e63SBarry Smith 629a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 6303b224e63SBarry Smith ilink = jac->head; 63149bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6324aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 633fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6343b224e63SBarry Smith ilink = ilink->next; 63549bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6364aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 637fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 63820252d06SBarry Smith 639f5236f50SJed Brown /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */ 64020252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 64120252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 64220252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 64320252d06SBarry Smith 644bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 64520252d06SBarry Smith if (sp) { 64620252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 64720252d06SBarry Smith } 64820252d06SBarry Smith 64920252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 6500298fd71SBarry Smith ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr); 651514bf10dSMatthew G Knepley if (flg) { 652514bf10dSMatthew G Knepley DM dmInner; 65321635b76SJed Brown KSP kspInner; 65421635b76SJed Brown 65521635b76SJed Brown ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 65621635b76SJed Brown ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 65721635b76SJed Brown /* Indent this deeper to emphasize the "inner" nature of this solver. */ 65821635b76SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr); 65921635b76SJed Brown ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr); 660514bf10dSMatthew G Knepley 661514bf10dSMatthew G Knepley /* Set DM for new solver */ 662bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 66321635b76SJed Brown ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr); 66421635b76SJed Brown ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr); 665514bf10dSMatthew G Knepley } else { 66621635b76SJed Brown /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or 66721635b76SJed Brown * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact, 66821635b76SJed Brown * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for 66921635b76SJed Brown * S = D - C A_inner^{-1} B, we expect S to be defined using an accurate definition of A_inner^{-1}, so we make 67021635b76SJed Brown * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used 67121635b76SJed Brown * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */ 67221635b76SJed Brown ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr); 673514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr); 674bafc1b83SMatthew G Knepley } 6755a9f2f41SSatish Balay ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 6765a9f2f41SSatish Balay ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 67797bbdb24SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 67897bbdb24SBarry Smith 679443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 6800298fd71SBarry Smith ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr); 681443836d0SMatthew G Knepley if (flg) { 682443836d0SMatthew G Knepley DM dmInner; 683443836d0SMatthew G Knepley 684443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 68582f516ccSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr); 686443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 687443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 688443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 689443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 690443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 691443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 692443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 693443836d0SMatthew G Knepley } else { 694443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 695443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 696443836d0SMatthew G Knepley } 697443836d0SMatthew G Knepley 698ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr); 69997bbdb24SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 70020252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 70197bbdb24SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 70297bbdb24SBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 7037233a360SDmitry Karpeev PC pcschur; 7047233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 7057233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 70697bbdb24SBarry Smith /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 70797bbdb24SBarry Smith } 70897bbdb24SBarry Smith ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr); 70997bbdb24SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix);CHKERRQ(ierr); 710*b20b4189SMatthew G. Knepley /* propogate DM */ 711*b20b4189SMatthew G. Knepley { 712*b20b4189SMatthew G. Knepley DM sdm; 713*b20b4189SMatthew G. Knepley ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr); 714*b20b4189SMatthew G. Knepley if (sdm) { 715*b20b4189SMatthew G. Knepley ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr); 716*b20b4189SMatthew G. Knepley ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr); 717*b20b4189SMatthew G. Knepley } 718*b20b4189SMatthew G. Knepley } 71997bbdb24SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 72097bbdb24SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 72197bbdb24SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 72297bbdb24SBarry Smith } 72397bbdb24SBarry Smith 7245a9f2f41SSatish Balay /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 7258caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 726519d70e2SJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 7273b224e63SBarry Smith if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 728c1570756SJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 7298caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 73097bbdb24SBarry Smith ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 7315a9f2f41SSatish Balay if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 7320971522cSBarry Smith if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 73397bbdb24SBarry Smith } else { 73468bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 7350971522cSBarry Smith i = 0; 7360971522cSBarry Smith ilink = jac->head; 7370971522cSBarry Smith while (ilink) { 7380971522cSBarry Smith ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 7390971522cSBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 7400971522cSBarry Smith if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 7410971522cSBarry Smith i++; 7420971522cSBarry Smith ilink = ilink->next; 7430971522cSBarry Smith } 7443b224e63SBarry Smith } 7453b224e63SBarry Smith 746c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 7470971522cSBarry Smith PetscFunctionReturn(0); 7480971522cSBarry Smith } 7490971522cSBarry Smith 7505a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 751ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 752ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 7535a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 754ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 755ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 75679416396SBarry Smith 7570971522cSBarry Smith #undef __FUNCT__ 7583b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7593b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7603b224e63SBarry Smith { 7613b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7623b224e63SBarry Smith PetscErrorCode ierr; 7633b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 764443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7653b224e63SBarry Smith 7663b224e63SBarry Smith PetscFunctionBegin; 767c5d2311dSJed Brown switch (jac->schurfactorization) { 768c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 769a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 770c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 771c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 772c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 773443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 774c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 775c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 776c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 777c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 778c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 779c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 780c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 781c5d2311dSJed Brown break; 782c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 783a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 784c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 785c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 786443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 787c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 788c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 789c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 790c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 791c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 792c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 793c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 794c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 795c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 796c5d2311dSJed Brown break; 797c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 798a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 799c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 800c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 801c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 802c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 803c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 804c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 805c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 806c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 807443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 808c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 809c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 810c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 811c5d2311dSJed Brown break; 812c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 813a04f6461SBarry 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 */ 8143b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8153b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 816443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 8173b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 8183b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 8193b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8203b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8213b224e63SBarry Smith 8223b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 8233b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8243b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8253b224e63SBarry Smith 826443836d0SMatthew G Knepley if (kspUpper == kspA) { 8273b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 8283b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 829443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 830443836d0SMatthew G Knepley } else { 831443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 832443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 833443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 834443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 835443836d0SMatthew G Knepley } 8363b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8373b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 838c5d2311dSJed Brown } 8393b224e63SBarry Smith PetscFunctionReturn(0); 8403b224e63SBarry Smith } 8413b224e63SBarry Smith 8423b224e63SBarry Smith #undef __FUNCT__ 8430971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 8440971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 8450971522cSBarry Smith { 8460971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8470971522cSBarry Smith PetscErrorCode ierr; 8485a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 849939b8a20SBarry Smith PetscInt cnt,bs; 8500971522cSBarry Smith 8510971522cSBarry Smith PetscFunctionBegin; 85279416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8531b9fc7fcSBarry Smith if (jac->defaultsplit) { 854939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 855ce94432eSBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 856939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 857ce94432eSBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 8580971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8595a9f2f41SSatish Balay while (ilink) { 8605a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8615a9f2f41SSatish Balay ilink = ilink->next; 8620971522cSBarry Smith } 8630971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8641b9fc7fcSBarry Smith } else { 865efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8665a9f2f41SSatish Balay while (ilink) { 8675a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8685a9f2f41SSatish Balay ilink = ilink->next; 8691b9fc7fcSBarry Smith } 8701b9fc7fcSBarry Smith } 87116913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 87279416396SBarry Smith if (!jac->w1) { 87379416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 87479416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 87579416396SBarry Smith } 876efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8775a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8783e197d65SBarry Smith cnt = 1; 8795a9f2f41SSatish Balay while (ilink->next) { 8805a9f2f41SSatish Balay ilink = ilink->next; 8813e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8823e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8833e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8843e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8853e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8863e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8873e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8883e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8893e197d65SBarry Smith } 89051f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 89111755939SBarry Smith cnt -= 2; 89251f519a2SBarry Smith while (ilink->previous) { 89351f519a2SBarry Smith ilink = ilink->previous; 89411755939SBarry Smith /* compute the residual only over the part of the vector needed */ 89511755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 89611755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 89711755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 89811755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 89911755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 90011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 90111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 90251f519a2SBarry Smith } 90311755939SBarry Smith } 904ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 9050971522cSBarry Smith PetscFunctionReturn(0); 9060971522cSBarry Smith } 9070971522cSBarry Smith 908421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 909ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 910ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 911421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 912ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 913ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 914421e10b8SBarry Smith 915421e10b8SBarry Smith #undef __FUNCT__ 9168c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 917421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 918421e10b8SBarry Smith { 919421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 920421e10b8SBarry Smith PetscErrorCode ierr; 921421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 922939b8a20SBarry Smith PetscInt bs; 923421e10b8SBarry Smith 924421e10b8SBarry Smith PetscFunctionBegin; 925421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 926421e10b8SBarry Smith if (jac->defaultsplit) { 927939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 928ce94432eSBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 929939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 930ce94432eSBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 931421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 932421e10b8SBarry Smith while (ilink) { 933421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 934421e10b8SBarry Smith ilink = ilink->next; 935421e10b8SBarry Smith } 936421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 937421e10b8SBarry Smith } else { 938421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 939421e10b8SBarry Smith while (ilink) { 940421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 941421e10b8SBarry Smith ilink = ilink->next; 942421e10b8SBarry Smith } 943421e10b8SBarry Smith } 944421e10b8SBarry Smith } else { 945421e10b8SBarry Smith if (!jac->w1) { 946421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 947421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 948421e10b8SBarry Smith } 949421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 950421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 951421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 952421e10b8SBarry Smith while (ilink->next) { 953421e10b8SBarry Smith ilink = ilink->next; 9549989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 955421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 956421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 957421e10b8SBarry Smith } 958421e10b8SBarry Smith while (ilink->previous) { 959421e10b8SBarry Smith ilink = ilink->previous; 9609989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 961421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 962421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 963421e10b8SBarry Smith } 964421e10b8SBarry Smith } else { 965421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 966421e10b8SBarry Smith ilink = ilink->next; 967421e10b8SBarry Smith } 968421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 969421e10b8SBarry Smith while (ilink->previous) { 970421e10b8SBarry Smith ilink = ilink->previous; 9719989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 972421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 973421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 974421e10b8SBarry Smith } 975421e10b8SBarry Smith } 976421e10b8SBarry Smith } 977421e10b8SBarry Smith PetscFunctionReturn(0); 978421e10b8SBarry Smith } 979421e10b8SBarry Smith 9800971522cSBarry Smith #undef __FUNCT__ 981574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 982574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9830971522cSBarry Smith { 9840971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9850971522cSBarry Smith PetscErrorCode ierr; 9865a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9870971522cSBarry Smith 9880971522cSBarry Smith PetscFunctionBegin; 9895a9f2f41SSatish Balay while (ilink) { 990574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 991fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 992fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 993443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 994fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 995fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 996b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9975a9f2f41SSatish Balay next = ilink->next; 9985a9f2f41SSatish Balay ilink = next; 9990971522cSBarry Smith } 100005b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1001574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 1002574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 1003574deadeSBarry Smith } else if (jac->mat) { 10040298fd71SBarry Smith jac->mat = NULL; 1005574deadeSBarry Smith } 100697bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 100768dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 10086bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 10096bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 10106bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 10116bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 10126bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 1013d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 10146bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 10156bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 101663ec74ffSBarry Smith jac->reset = PETSC_TRUE; 1017574deadeSBarry Smith PetscFunctionReturn(0); 1018574deadeSBarry Smith } 1019574deadeSBarry Smith 1020574deadeSBarry Smith #undef __FUNCT__ 1021574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1022574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1023574deadeSBarry Smith { 1024574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1025574deadeSBarry Smith PetscErrorCode ierr; 1026574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 1027574deadeSBarry Smith 1028574deadeSBarry Smith PetscFunctionBegin; 1029574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 1030574deadeSBarry Smith while (ilink) { 10316bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 1032574deadeSBarry Smith next = ilink->next; 1033574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 1034574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 10355d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 1036574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 1037574deadeSBarry Smith ilink = next; 1038574deadeSBarry Smith } 1039574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1040c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1041bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr); 1042bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr); 1043bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr); 1044bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr); 1045bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr); 1046bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",NULL);CHKERRQ(ierr); 1047bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr); 10480971522cSBarry Smith PetscFunctionReturn(0); 10490971522cSBarry Smith } 10500971522cSBarry Smith 10510971522cSBarry Smith #undef __FUNCT__ 10520971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10530971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10540971522cSBarry Smith { 10551b9fc7fcSBarry Smith PetscErrorCode ierr; 10566c924f48SJed Brown PetscInt bs; 1057bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10589dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10593b224e63SBarry Smith PCCompositeType ctype; 10601b9fc7fcSBarry Smith 10610971522cSBarry Smith PetscFunctionBegin; 10621b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 10634ab8060aSDmitry Karpeev ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr); 106451f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 106551f519a2SBarry Smith if (flg) { 106651f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 106751f519a2SBarry Smith } 1068704ba839SBarry Smith 10690298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 1070c0adfefeSBarry Smith if (stokes) { 1071c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1072c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1073c0adfefeSBarry Smith } 1074c0adfefeSBarry Smith 10753b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10763b224e63SBarry Smith if (flg) { 10773b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10783b224e63SBarry Smith } 1079c30613efSMatthew Knepley /* Only setup fields once */ 1080c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1081d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1082d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10836c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10846c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1085d32f9abdSBarry Smith } 1086c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1087c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1088c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 10890298fd71SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_schur_fact_type","Which off-diagonal parts of the block factorization to use","PCFieldSplitSetSchurFactType",PCFieldSplitSchurFactTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,NULL);CHKERRQ(ierr); 10900298fd71SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);CHKERRQ(ierr); 1091c5d2311dSJed Brown } 10921b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10930971522cSBarry Smith PetscFunctionReturn(0); 10940971522cSBarry Smith } 10950971522cSBarry Smith 10960971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10970971522cSBarry Smith 10980971522cSBarry Smith #undef __FUNCT__ 10990971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 11001e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11010971522cSBarry Smith { 110297bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11030971522cSBarry Smith PetscErrorCode ierr; 11045a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 110569a612a9SBarry Smith char prefix[128]; 11065d4c12cdSJungho Lee PetscInt i; 11070971522cSBarry Smith 11080971522cSBarry Smith PetscFunctionBegin; 11096c924f48SJed Brown if (jac->splitdefined) { 11106c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11116c924f48SJed Brown PetscFunctionReturn(0); 11126c924f48SJed Brown } 111351f519a2SBarry Smith for (i=0; i<n; i++) { 1114e32f2f54SBarry 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); 1115e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 111651f519a2SBarry Smith } 1117704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1118a04f6461SBarry Smith if (splitname) { 1119db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1120a04f6461SBarry Smith } else { 1121a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1122a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1123a04f6461SBarry Smith } 1124704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 11255a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 11265d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 11275d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 11282fa5cd67SKarl Rupp 11295a9f2f41SSatish Balay ilink->nfields = n; 11300298fd71SBarry Smith ilink->next = NULL; 1131ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr); 113220252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11335a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11349005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 113569a612a9SBarry Smith 11368caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr); 11375a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11380971522cSBarry Smith 11390971522cSBarry Smith if (!next) { 11405a9f2f41SSatish Balay jac->head = ilink; 11410298fd71SBarry Smith ilink->previous = NULL; 11420971522cSBarry Smith } else { 11430971522cSBarry Smith while (next->next) { 11440971522cSBarry Smith next = next->next; 11450971522cSBarry Smith } 11465a9f2f41SSatish Balay next->next = ilink; 114751f519a2SBarry Smith ilink->previous = next; 11480971522cSBarry Smith } 11490971522cSBarry Smith jac->nsplits++; 11500971522cSBarry Smith PetscFunctionReturn(0); 11510971522cSBarry Smith } 11520971522cSBarry Smith 1153e69d4d44SBarry Smith #undef __FUNCT__ 1154e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11551e6b0712SBarry Smith static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1156e69d4d44SBarry Smith { 1157e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1158e69d4d44SBarry Smith PetscErrorCode ierr; 1159e69d4d44SBarry Smith 1160e69d4d44SBarry Smith PetscFunctionBegin; 1161519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1162e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 11632fa5cd67SKarl Rupp 1164e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 116513e0d083SBarry Smith if (n) *n = jac->nsplits; 1166e69d4d44SBarry Smith PetscFunctionReturn(0); 1167e69d4d44SBarry Smith } 11680971522cSBarry Smith 11690971522cSBarry Smith #undef __FUNCT__ 117069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11711e6b0712SBarry Smith static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11720971522cSBarry Smith { 11730971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11740971522cSBarry Smith PetscErrorCode ierr; 11750971522cSBarry Smith PetscInt cnt = 0; 11765a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11770971522cSBarry Smith 11780971522cSBarry Smith PetscFunctionBegin; 11795d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11805a9f2f41SSatish Balay while (ilink) { 11815a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11825a9f2f41SSatish Balay ilink = ilink->next; 11830971522cSBarry Smith } 11845d480477SMatthew 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); 118513e0d083SBarry Smith if (n) *n = jac->nsplits; 11860971522cSBarry Smith PetscFunctionReturn(0); 11870971522cSBarry Smith } 11880971522cSBarry Smith 1189704ba839SBarry Smith #undef __FUNCT__ 1190704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11911e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1192704ba839SBarry Smith { 1193704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1194704ba839SBarry Smith PetscErrorCode ierr; 1195704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1196704ba839SBarry Smith char prefix[128]; 1197704ba839SBarry Smith 1198704ba839SBarry Smith PetscFunctionBegin; 11996c924f48SJed Brown if (jac->splitdefined) { 12006c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 12016c924f48SJed Brown PetscFunctionReturn(0); 12026c924f48SJed Brown } 120316913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1204a04f6461SBarry Smith if (splitname) { 1205db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1206a04f6461SBarry Smith } else { 1207b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1208b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1209a04f6461SBarry Smith } 12101ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1211b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1212b5787286SJed Brown ilink->is = is; 1213b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1214b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1215b5787286SJed Brown ilink->is_col = is; 12160298fd71SBarry Smith ilink->next = NULL; 1217ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr); 121820252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1219704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 12209005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1221704ba839SBarry Smith 12228caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr); 1223704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1224704ba839SBarry Smith 1225704ba839SBarry Smith if (!next) { 1226704ba839SBarry Smith jac->head = ilink; 12270298fd71SBarry Smith ilink->previous = NULL; 1228704ba839SBarry Smith } else { 1229704ba839SBarry Smith while (next->next) { 1230704ba839SBarry Smith next = next->next; 1231704ba839SBarry Smith } 1232704ba839SBarry Smith next->next = ilink; 1233704ba839SBarry Smith ilink->previous = next; 1234704ba839SBarry Smith } 1235704ba839SBarry Smith jac->nsplits++; 1236704ba839SBarry Smith PetscFunctionReturn(0); 1237704ba839SBarry Smith } 1238704ba839SBarry Smith 12390971522cSBarry Smith #undef __FUNCT__ 12400971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12410971522cSBarry Smith /*@ 12420971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12430971522cSBarry Smith 1244ad4df100SBarry Smith Logically Collective on PC 12450971522cSBarry Smith 12460971522cSBarry Smith Input Parameters: 12470971522cSBarry Smith + pc - the preconditioner context 12480298fd71SBarry Smith . splitname - name of this split, if NULL the number of the split is used 12490971522cSBarry Smith . n - the number of fields in this split 1250db4c96c1SJed Brown - fields - the fields in this split 12510971522cSBarry Smith 12520971522cSBarry Smith Level: intermediate 12530971522cSBarry Smith 1254d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1255d32f9abdSBarry Smith 12567287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1257d32f9abdSBarry 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 1258d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1259d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1260d32f9abdSBarry Smith 1261db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1262db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1263db4c96c1SJed Brown 12645d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12655d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12665d4c12cdSJungho Lee available when this routine is called. 12675d4c12cdSJungho Lee 1268d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12690971522cSBarry Smith 12700971522cSBarry Smith @*/ 12715d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12720971522cSBarry Smith { 12734ac538c5SBarry Smith PetscErrorCode ierr; 12740971522cSBarry Smith 12750971522cSBarry Smith PetscFunctionBegin; 12760700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1277db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1278ce94432eSBarry Smith if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1279db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12805d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12810971522cSBarry Smith PetscFunctionReturn(0); 12820971522cSBarry Smith } 12830971522cSBarry Smith 12840971522cSBarry Smith #undef __FUNCT__ 1285704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1286704ba839SBarry Smith /*@ 1287704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1288704ba839SBarry Smith 1289ad4df100SBarry Smith Logically Collective on PC 1290704ba839SBarry Smith 1291704ba839SBarry Smith Input Parameters: 1292704ba839SBarry Smith + pc - the preconditioner context 12930298fd71SBarry Smith . splitname - name of this split, if NULL the number of the split is used 1294db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1295704ba839SBarry Smith 1296d32f9abdSBarry Smith 1297a6ffb8dbSJed Brown Notes: 1298a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1299a6ffb8dbSJed Brown 1300db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1301db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1302d32f9abdSBarry Smith 1303704ba839SBarry Smith Level: intermediate 1304704ba839SBarry Smith 1305704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1306704ba839SBarry Smith 1307704ba839SBarry Smith @*/ 13087087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1309704ba839SBarry Smith { 13104ac538c5SBarry Smith PetscErrorCode ierr; 1311704ba839SBarry Smith 1312704ba839SBarry Smith PetscFunctionBegin; 13130700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13147b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1315db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 13164ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1317704ba839SBarry Smith PetscFunctionReturn(0); 1318704ba839SBarry Smith } 1319704ba839SBarry Smith 1320704ba839SBarry Smith #undef __FUNCT__ 132157a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 132257a9adfeSMatthew G Knepley /*@ 132357a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 132457a9adfeSMatthew G Knepley 132557a9adfeSMatthew G Knepley Logically Collective on PC 132657a9adfeSMatthew G Knepley 132757a9adfeSMatthew G Knepley Input Parameters: 132857a9adfeSMatthew G Knepley + pc - the preconditioner context 132957a9adfeSMatthew G Knepley - splitname - name of this split 133057a9adfeSMatthew G Knepley 133157a9adfeSMatthew G Knepley Output Parameter: 13320298fd71SBarry Smith - is - the index set that defines the vector elements in this field, or NULL if the field is not found 133357a9adfeSMatthew G Knepley 133457a9adfeSMatthew G Knepley Level: intermediate 133557a9adfeSMatthew G Knepley 133657a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 133757a9adfeSMatthew G Knepley 133857a9adfeSMatthew G Knepley @*/ 133957a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 134057a9adfeSMatthew G Knepley { 134157a9adfeSMatthew G Knepley PetscErrorCode ierr; 134257a9adfeSMatthew G Knepley 134357a9adfeSMatthew G Knepley PetscFunctionBegin; 134457a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 134557a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 134657a9adfeSMatthew G Knepley PetscValidPointer(is,3); 134757a9adfeSMatthew G Knepley { 134857a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 134957a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 135057a9adfeSMatthew G Knepley PetscBool found; 135157a9adfeSMatthew G Knepley 13520298fd71SBarry Smith *is = NULL; 135357a9adfeSMatthew G Knepley while (ilink) { 135457a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 135557a9adfeSMatthew G Knepley if (found) { 135657a9adfeSMatthew G Knepley *is = ilink->is; 135757a9adfeSMatthew G Knepley break; 135857a9adfeSMatthew G Knepley } 135957a9adfeSMatthew G Knepley ilink = ilink->next; 136057a9adfeSMatthew G Knepley } 136157a9adfeSMatthew G Knepley } 136257a9adfeSMatthew G Knepley PetscFunctionReturn(0); 136357a9adfeSMatthew G Knepley } 136457a9adfeSMatthew G Knepley 136557a9adfeSMatthew G Knepley #undef __FUNCT__ 136651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 136751f519a2SBarry Smith /*@ 136851f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 136951f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 137051f519a2SBarry Smith 1371ad4df100SBarry Smith Logically Collective on PC 137251f519a2SBarry Smith 137351f519a2SBarry Smith Input Parameters: 137451f519a2SBarry Smith + pc - the preconditioner context 137551f519a2SBarry Smith - bs - the block size 137651f519a2SBarry Smith 137751f519a2SBarry Smith Level: intermediate 137851f519a2SBarry Smith 137951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 138051f519a2SBarry Smith 138151f519a2SBarry Smith @*/ 13827087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 138351f519a2SBarry Smith { 13844ac538c5SBarry Smith PetscErrorCode ierr; 138551f519a2SBarry Smith 138651f519a2SBarry Smith PetscFunctionBegin; 13870700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1388c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13894ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 139051f519a2SBarry Smith PetscFunctionReturn(0); 139151f519a2SBarry Smith } 139251f519a2SBarry Smith 139351f519a2SBarry Smith #undef __FUNCT__ 139469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 13950971522cSBarry Smith /*@C 139669a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 13970971522cSBarry Smith 139869a612a9SBarry Smith Collective on KSP 13990971522cSBarry Smith 14000971522cSBarry Smith Input Parameter: 14010971522cSBarry Smith . pc - the preconditioner context 14020971522cSBarry Smith 14030971522cSBarry Smith Output Parameters: 140413e0d083SBarry Smith + n - the number of splits 140569a612a9SBarry Smith - pc - the array of KSP contexts 14060971522cSBarry Smith 14070971522cSBarry Smith Note: 1408d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1409d32f9abdSBarry Smith (not the KSP just the array that contains them). 14100971522cSBarry Smith 141169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 14120971522cSBarry Smith 1413196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 14140298fd71SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the 1415196cc216SBarry Smith KSP array must be. 1416196cc216SBarry Smith 1417196cc216SBarry Smith 14180971522cSBarry Smith Level: advanced 14190971522cSBarry Smith 14200971522cSBarry Smith .seealso: PCFIELDSPLIT 14210971522cSBarry Smith @*/ 14227087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 14230971522cSBarry Smith { 14244ac538c5SBarry Smith PetscErrorCode ierr; 14250971522cSBarry Smith 14260971522cSBarry Smith PetscFunctionBegin; 14270700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 142813e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14294ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14300971522cSBarry Smith PetscFunctionReturn(0); 14310971522cSBarry Smith } 14320971522cSBarry Smith 1433e69d4d44SBarry Smith #undef __FUNCT__ 1434e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1435e69d4d44SBarry Smith /*@ 1436e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1437a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1438e69d4d44SBarry Smith 1439e69d4d44SBarry Smith Collective on PC 1440e69d4d44SBarry Smith 1441e69d4d44SBarry Smith Input Parameters: 1442e69d4d44SBarry Smith + pc - the preconditioner context 1443e87fab1bSBarry Smith . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default 14440298fd71SBarry Smith - userpre - matrix to use for preconditioning, or NULL 1445084e4875SJed Brown 1446e69d4d44SBarry Smith Options Database: 1447e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> default is a11 1448e69d4d44SBarry Smith 1449fd1303e9SJungho Lee Notes: 1450fd1303e9SJungho Lee $ If ptype is 1451fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1452fd1303e9SJungho Lee $ to this function). 1453e87fab1bSBarry Smith $ a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1454fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1455fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1456fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1457fd1303e9SJungho Lee $ preconditioner 1458fd1303e9SJungho Lee 1459e87fab1bSBarry Smith When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense 1460fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1461fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1462fd1303e9SJungho Lee 1463fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1464fd1303e9SJungho Lee the name since it sets a proceedure to use. 1465fd1303e9SJungho Lee 1466e69d4d44SBarry Smith Level: intermediate 1467e69d4d44SBarry Smith 1468fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1469e69d4d44SBarry Smith 1470e69d4d44SBarry Smith @*/ 14717087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1472e69d4d44SBarry Smith { 14734ac538c5SBarry Smith PetscErrorCode ierr; 1474e69d4d44SBarry Smith 1475e69d4d44SBarry Smith PetscFunctionBegin; 14760700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14774ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1478e69d4d44SBarry Smith PetscFunctionReturn(0); 1479e69d4d44SBarry Smith } 1480e69d4d44SBarry Smith 1481e69d4d44SBarry Smith #undef __FUNCT__ 1482e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14831e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1484e69d4d44SBarry Smith { 1485e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1486084e4875SJed Brown PetscErrorCode ierr; 1487e69d4d44SBarry Smith 1488e69d4d44SBarry Smith PetscFunctionBegin; 1489084e4875SJed Brown jac->schurpre = ptype; 1490084e4875SJed Brown if (pre) { 14916bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1492084e4875SJed Brown jac->schur_user = pre; 1493084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1494084e4875SJed Brown } 1495e69d4d44SBarry Smith PetscFunctionReturn(0); 1496e69d4d44SBarry Smith } 1497e69d4d44SBarry Smith 149830ad9308SMatthew Knepley #undef __FUNCT__ 1499c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1500ab1df9f5SJed Brown /*@ 1501c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1502ab1df9f5SJed Brown 1503ab1df9f5SJed Brown Collective on PC 1504ab1df9f5SJed Brown 1505ab1df9f5SJed Brown Input Parameters: 1506ab1df9f5SJed Brown + pc - the preconditioner context 1507c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1508ab1df9f5SJed Brown 1509ab1df9f5SJed Brown Options Database: 1510c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1511ab1df9f5SJed Brown 1512ab1df9f5SJed Brown 1513ab1df9f5SJed Brown Level: intermediate 1514ab1df9f5SJed Brown 1515ab1df9f5SJed Brown Notes: 1516ab1df9f5SJed Brown The FULL factorization is 1517ab1df9f5SJed Brown 1518ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1519ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1520ab1df9f5SJed Brown 15216be4592eSBarry Smith where S = D - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D, 15226be4592eSBarry Smith and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES). 1523ab1df9f5SJed Brown 15246be4592eSBarry Smith If applied exactly, FULL factorization is a direct solver. The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial 15256be4592eSBarry Smith of degree 2, so KSPGMRES converges in 2 iterations. If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner 15266be4592eSBarry Smith application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice. With DIAG, 15276be4592eSBarry Smith the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations. 1528ab1df9f5SJed Brown 15296be4592eSBarry Smith For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES. Note that a flexible method like KSPFGMRES 15306be4592eSBarry Smith or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used inside a split). 1531ab1df9f5SJed Brown 1532ab1df9f5SJed Brown References: 1533ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1534ab1df9f5SJed Brown 1535ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1536ab1df9f5SJed Brown 1537ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1538ab1df9f5SJed Brown @*/ 1539c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1540ab1df9f5SJed Brown { 1541ab1df9f5SJed Brown PetscErrorCode ierr; 1542ab1df9f5SJed Brown 1543ab1df9f5SJed Brown PetscFunctionBegin; 1544ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1545c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1546ab1df9f5SJed Brown PetscFunctionReturn(0); 1547ab1df9f5SJed Brown } 1548ab1df9f5SJed Brown 1549ab1df9f5SJed Brown #undef __FUNCT__ 1550c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 15511e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1552ab1df9f5SJed Brown { 1553ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1554ab1df9f5SJed Brown 1555ab1df9f5SJed Brown PetscFunctionBegin; 1556ab1df9f5SJed Brown jac->schurfactorization = ftype; 1557ab1df9f5SJed Brown PetscFunctionReturn(0); 1558ab1df9f5SJed Brown } 1559ab1df9f5SJed Brown 1560ab1df9f5SJed Brown #undef __FUNCT__ 156130ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 156230ad9308SMatthew Knepley /*@C 15638c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 156430ad9308SMatthew Knepley 156530ad9308SMatthew Knepley Collective on KSP 156630ad9308SMatthew Knepley 156730ad9308SMatthew Knepley Input Parameter: 156830ad9308SMatthew Knepley . pc - the preconditioner context 156930ad9308SMatthew Knepley 157030ad9308SMatthew Knepley Output Parameters: 1571a04f6461SBarry Smith + A00 - the (0,0) block 1572a04f6461SBarry Smith . A01 - the (0,1) block 1573a04f6461SBarry Smith . A10 - the (1,0) block 1574a04f6461SBarry Smith - A11 - the (1,1) block 157530ad9308SMatthew Knepley 157630ad9308SMatthew Knepley Level: advanced 157730ad9308SMatthew Knepley 157830ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 157930ad9308SMatthew Knepley @*/ 1580a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 158130ad9308SMatthew Knepley { 158230ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 158330ad9308SMatthew Knepley 158430ad9308SMatthew Knepley PetscFunctionBegin; 15850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1586ce94432eSBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1587a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1588a04f6461SBarry Smith if (A01) *A01 = jac->B; 1589a04f6461SBarry Smith if (A10) *A10 = jac->C; 1590a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 159130ad9308SMatthew Knepley PetscFunctionReturn(0); 159230ad9308SMatthew Knepley } 159330ad9308SMatthew Knepley 159479416396SBarry Smith #undef __FUNCT__ 159579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 15961e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 159779416396SBarry Smith { 159879416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1599e69d4d44SBarry Smith PetscErrorCode ierr; 160079416396SBarry Smith 160179416396SBarry Smith PetscFunctionBegin; 160279416396SBarry Smith jac->type = type; 16033b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 16043b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 16053b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 16062fa5cd67SKarl Rupp 1607bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1608bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1609bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1610e69d4d44SBarry Smith 16113b224e63SBarry Smith } else { 16123b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 16133b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 16142fa5cd67SKarl Rupp 1615bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1616bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",0);CHKERRQ(ierr); 1617bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr); 16183b224e63SBarry Smith } 161979416396SBarry Smith PetscFunctionReturn(0); 162079416396SBarry Smith } 162179416396SBarry Smith 162251f519a2SBarry Smith #undef __FUNCT__ 162351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16241e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 162551f519a2SBarry Smith { 162651f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 162751f519a2SBarry Smith 162851f519a2SBarry Smith PetscFunctionBegin; 1629ce94432eSBarry Smith if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1630ce94432eSBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 163151f519a2SBarry Smith jac->bs = bs; 163251f519a2SBarry Smith PetscFunctionReturn(0); 163351f519a2SBarry Smith } 163451f519a2SBarry Smith 163579416396SBarry Smith #undef __FUNCT__ 163679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1637bc08b0f1SBarry Smith /*@ 163879416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 163979416396SBarry Smith 164079416396SBarry Smith Collective on PC 164179416396SBarry Smith 164279416396SBarry Smith Input Parameter: 164379416396SBarry Smith . pc - the preconditioner context 164481540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 164579416396SBarry Smith 164679416396SBarry Smith Options Database Key: 1647a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 164879416396SBarry Smith 1649b02e2d75SMatthew G Knepley Level: Intermediate 165079416396SBarry Smith 165179416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 165279416396SBarry Smith 165379416396SBarry Smith .seealso: PCCompositeSetType() 165479416396SBarry Smith 165579416396SBarry Smith @*/ 16567087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 165779416396SBarry Smith { 16584ac538c5SBarry Smith PetscErrorCode ierr; 165979416396SBarry Smith 166079416396SBarry Smith PetscFunctionBegin; 16610700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16624ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 166379416396SBarry Smith PetscFunctionReturn(0); 166479416396SBarry Smith } 166579416396SBarry Smith 1666b02e2d75SMatthew G Knepley #undef __FUNCT__ 1667b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1668b02e2d75SMatthew G Knepley /*@ 1669b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1670b02e2d75SMatthew G Knepley 1671b02e2d75SMatthew G Knepley Not collective 1672b02e2d75SMatthew G Knepley 1673b02e2d75SMatthew G Knepley Input Parameter: 1674b02e2d75SMatthew G Knepley . pc - the preconditioner context 1675b02e2d75SMatthew G Knepley 1676b02e2d75SMatthew G Knepley Output Parameter: 1677b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1678b02e2d75SMatthew G Knepley 1679b02e2d75SMatthew G Knepley Level: Intermediate 1680b02e2d75SMatthew G Knepley 1681b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1682b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1683b02e2d75SMatthew G Knepley @*/ 1684b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1685b02e2d75SMatthew G Knepley { 1686b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 1687b02e2d75SMatthew G Knepley 1688b02e2d75SMatthew G Knepley PetscFunctionBegin; 1689b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1690b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1691b02e2d75SMatthew G Knepley *type = jac->type; 1692b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1693b02e2d75SMatthew G Knepley } 1694b02e2d75SMatthew G Knepley 16954ab8060aSDmitry Karpeev #undef __FUNCT__ 16964ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDMSplits" 16974ab8060aSDmitry Karpeev /*@ 16984ab8060aSDmitry Karpeev PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible. 16994ab8060aSDmitry Karpeev 17004ab8060aSDmitry Karpeev Logically Collective 17014ab8060aSDmitry Karpeev 17024ab8060aSDmitry Karpeev Input Parameters: 17034ab8060aSDmitry Karpeev + pc - the preconditioner context 17044ab8060aSDmitry Karpeev - flg - boolean indicating whether to use field splits defined by the DM 17054ab8060aSDmitry Karpeev 17064ab8060aSDmitry Karpeev Options Database Key: 17074ab8060aSDmitry Karpeev . -pc_fieldsplit_dm_splits 17084ab8060aSDmitry Karpeev 17094ab8060aSDmitry Karpeev Level: Intermediate 17104ab8060aSDmitry Karpeev 17114ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative 17124ab8060aSDmitry Karpeev 17134ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits() 17144ab8060aSDmitry Karpeev 17154ab8060aSDmitry Karpeev @*/ 17164ab8060aSDmitry Karpeev PetscErrorCode PCFieldSplitSetDMSplits(PC pc,PetscBool flg) 17174ab8060aSDmitry Karpeev { 17184ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 17194ab8060aSDmitry Karpeev PetscBool isfs; 17204ab8060aSDmitry Karpeev PetscErrorCode ierr; 17214ab8060aSDmitry Karpeev 17224ab8060aSDmitry Karpeev PetscFunctionBegin; 17234ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17244ab8060aSDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 17254ab8060aSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr); 17264ab8060aSDmitry Karpeev if (isfs) { 17274ab8060aSDmitry Karpeev jac->dm_splits = flg; 17284ab8060aSDmitry Karpeev } 17294ab8060aSDmitry Karpeev PetscFunctionReturn(0); 17304ab8060aSDmitry Karpeev } 17314ab8060aSDmitry Karpeev 17324ab8060aSDmitry Karpeev 17334ab8060aSDmitry Karpeev #undef __FUNCT__ 17344ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDMSplits" 17354ab8060aSDmitry Karpeev /*@ 17364ab8060aSDmitry Karpeev PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible. 17374ab8060aSDmitry Karpeev 17384ab8060aSDmitry Karpeev Logically Collective 17394ab8060aSDmitry Karpeev 17404ab8060aSDmitry Karpeev Input Parameter: 17414ab8060aSDmitry Karpeev . pc - the preconditioner context 17424ab8060aSDmitry Karpeev 17434ab8060aSDmitry Karpeev Output Parameter: 17444ab8060aSDmitry Karpeev . flg - boolean indicating whether to use field splits defined by the DM 17454ab8060aSDmitry Karpeev 17464ab8060aSDmitry Karpeev Level: Intermediate 17474ab8060aSDmitry Karpeev 17484ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative 17494ab8060aSDmitry Karpeev 17504ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits() 17514ab8060aSDmitry Karpeev 17524ab8060aSDmitry Karpeev @*/ 17534ab8060aSDmitry Karpeev PetscErrorCode PCFieldSplitGetDMSplits(PC pc,PetscBool* flg) 17544ab8060aSDmitry Karpeev { 17554ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 17564ab8060aSDmitry Karpeev PetscBool isfs; 17574ab8060aSDmitry Karpeev PetscErrorCode ierr; 17584ab8060aSDmitry Karpeev 17594ab8060aSDmitry Karpeev PetscFunctionBegin; 17604ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17614ab8060aSDmitry Karpeev PetscValidPointer(flg,2); 17624ab8060aSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr); 17634ab8060aSDmitry Karpeev if (isfs) { 17644ab8060aSDmitry Karpeev if(flg) *flg = jac->dm_splits; 17654ab8060aSDmitry Karpeev } 17664ab8060aSDmitry Karpeev PetscFunctionReturn(0); 17674ab8060aSDmitry Karpeev } 17684ab8060aSDmitry Karpeev 17690971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 17700971522cSBarry Smith /*MC 1771a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1772a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 17730971522cSBarry Smith 1774edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1775edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 17760971522cSBarry Smith 1777a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 177869a612a9SBarry Smith and set the options directly on the resulting KSP object 17790971522cSBarry Smith 17800971522cSBarry Smith Level: intermediate 17810971522cSBarry Smith 178279416396SBarry Smith Options Database Keys: 178381540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 178481540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 178581540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 178681540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 17870f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 1788e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11 1789435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 179079416396SBarry Smith 17915d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 17925d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 17935d4c12cdSJungho Lee 1794c8a0d604SMatthew G Knepley Notes: 1795c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1796d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1797d32f9abdSBarry Smith 1798d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1799d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1800d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1801d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1802d32f9abdSBarry Smith 1803c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1804c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1805c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1806c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1807c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1808a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1809c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1810c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1811c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1812a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1813c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1814c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 18155668aaf4SBarry Smith diag gives 1816c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1817c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 18185668aaf4SBarry Smith note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of 1819c8a0d604SMatthew G Knepley $ ( A00 0 ) 1820c8a0d604SMatthew G Knepley $ ( A10 S ) 1821c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1822c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1823c8a0d604SMatthew G Knepley $ ( 0 S ) 1824c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1825e69d4d44SBarry Smith 1826edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1827edf189efSBarry Smith is used automatically for a second block. 1828edf189efSBarry Smith 1829ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1830ff218e97SBarry Smith Generally it should be used with the AIJ format. 1831ff218e97SBarry Smith 1832ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1833ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1834ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 18350716a85fSBarry Smith 1836a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 18370971522cSBarry Smith 18387e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1839e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 18400971522cSBarry Smith M*/ 18410971522cSBarry Smith 18420971522cSBarry Smith #undef __FUNCT__ 18430971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 18448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc) 18450971522cSBarry Smith { 18460971522cSBarry Smith PetscErrorCode ierr; 18470971522cSBarry Smith PC_FieldSplit *jac; 18480971522cSBarry Smith 18490971522cSBarry Smith PetscFunctionBegin; 185038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 18512fa5cd67SKarl Rupp 18520971522cSBarry Smith jac->bs = -1; 18530971522cSBarry Smith jac->nsplits = 0; 18543e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1855e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1856c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 1857fbe7908bSJed Brown jac->dm_splits = PETSC_TRUE; 185851f519a2SBarry Smith 18590971522cSBarry Smith pc->data = (void*)jac; 18600971522cSBarry Smith 18610971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1862421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 18630971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1864574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 18650971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 18660971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 18670971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 18680971522cSBarry Smith pc->ops->applyrichardson = 0; 18690971522cSBarry Smith 1870bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1871bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1872bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 1873bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 1874bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 18750971522cSBarry Smith PetscFunctionReturn(0); 18760971522cSBarry Smith } 18770971522cSBarry Smith 18780971522cSBarry Smith 1879a541d17aSBarry Smith 1880