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); 199*06de4afeSJed Brown if (jac->head) { 200443836d0SMatthew G Knepley ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 201*06de4afeSJed Brown } else {ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr);} 2023b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 203*06de4afeSJed Brown if (jac->head && jac->kspupper != jac->head->ksp) { 204443836d0SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr); 205443836d0SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 206b2750c55SJed Brown if (jac->kspupper) {ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);} 207b2750c55SJed Brown else {ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr);} 208443836d0SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 209443836d0SMatthew G Knepley } 210435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 2113b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 21212cae6f2SJed Brown if (jac->kspschur) { 2133b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 21412cae6f2SJed Brown } else { 21512cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 21612cae6f2SJed Brown } 2173b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2183b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 219*06de4afeSJed Brown } else if (isdraw && jac->head) { 2204996c5bdSBarry Smith PetscDraw draw; 2214996c5bdSBarry Smith PetscReal x,y,w,wd,h; 2224996c5bdSBarry Smith PetscInt cnt = 2; 2234996c5bdSBarry Smith char str[32]; 2244996c5bdSBarry Smith 2254996c5bdSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 2264996c5bdSBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 227c74581afSBarry Smith if (jac->kspupper != jac->head->ksp) cnt++; 228c74581afSBarry Smith w = 2*PetscMin(1.0 - x,x); 229c74581afSBarry Smith wd = w/(cnt + 1); 230c74581afSBarry Smith 2314996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 2320298fd71SBarry Smith ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 2334996c5bdSBarry Smith y -= h; 2344996c5bdSBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) { 235e87fab1bSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr); 2363b224e63SBarry Smith } else { 2374996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr); 2384996c5bdSBarry Smith } 2390298fd71SBarry Smith ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 2404996c5bdSBarry Smith y -= h; 2414996c5bdSBarry Smith x = x - wd*(cnt-1)/2.0; 2424996c5bdSBarry Smith 2434996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2444996c5bdSBarry Smith ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 2454996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2464996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) { 2474996c5bdSBarry Smith x += wd; 2484996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2494996c5bdSBarry Smith ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 2504996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2514996c5bdSBarry Smith } 2524996c5bdSBarry Smith x += wd; 2534996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2544996c5bdSBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 2554996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2563b224e63SBarry Smith } 2573b224e63SBarry Smith PetscFunctionReturn(0); 2583b224e63SBarry Smith } 2593b224e63SBarry Smith 2603b224e63SBarry Smith #undef __FUNCT__ 2616c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2626c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2636c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2646c924f48SJed Brown { 2656c924f48SJed Brown PetscErrorCode ierr; 2666c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2675d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2685d4c12cdSJungho Lee PetscBool flg,flg_col; 2695d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2706c924f48SJed Brown 2716c924f48SJed Brown PetscFunctionBegin; 2726c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2735d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2746c924f48SJed Brown for (i=0,flg=PETSC_TRUE;; i++) { 2758caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2768caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2778caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2786c924f48SJed Brown nfields = jac->bs; 27929499fbbSJungho Lee nfields_col = jac->bs; 2806c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2815d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2826c924f48SJed Brown if (!flg) break; 2835d4c12cdSJungho Lee else if (flg && !flg_col) { 2845d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2855d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2862fa5cd67SKarl Rupp } else { 2875d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2885d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2895d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2905d4c12cdSJungho Lee } 2916c924f48SJed Brown } 2926c924f48SJed Brown if (i > 0) { 2936c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2946c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2956c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2966c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2976c924f48SJed Brown } 2986c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2995d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 3006c924f48SJed Brown PetscFunctionReturn(0); 3016c924f48SJed Brown } 3026c924f48SJed Brown 3036c924f48SJed Brown #undef __FUNCT__ 30469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 30569a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 3060971522cSBarry Smith { 3070971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3080971522cSBarry Smith PetscErrorCode ierr; 3095a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3107287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 3116c924f48SJed Brown PetscInt i; 3120971522cSBarry Smith 3130971522cSBarry Smith PetscFunctionBegin; 3147287d2a3SDmitry Karpeev /* 3157287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 3167287d2a3SDmitry Karpeev Should probably be rewritten. 3177287d2a3SDmitry Karpeev */ 3187287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 3190298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 3204ab8060aSDmitry Karpeev if (pc->dm && jac->dm_splits && !stokes) { 321bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 3220784a22cSJed Brown char **fieldNames; 3237b62db95SJungho Lee IS *fields; 324e7c4fc90SDmitry Karpeev DM *dms; 325bafc1b83SMatthew G Knepley DM subdm[128]; 326bafc1b83SMatthew G Knepley PetscBool flg; 327bafc1b83SMatthew G Knepley 32816621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 329bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 330bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE;; i++) { 331bafc1b83SMatthew G Knepley PetscInt ifields[128]; 332bafc1b83SMatthew G Knepley IS compField; 333bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 334bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 335bafc1b83SMatthew G Knepley 3368caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 337bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 338bafc1b83SMatthew G Knepley if (!flg) break; 33982f516ccSBarry Smith if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 340bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 341bafc1b83SMatthew G Knepley if (nfields == 1) { 342bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 34382f516ccSBarry Smith /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 3440298fd71SBarry Smith ierr = ISView(compField, NULL);CHKERRQ(ierr); */ 345bafc1b83SMatthew G Knepley } else { 3468caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 347bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 34882f516ccSBarry Smith /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", splitname);CHKERRQ(ierr); 3490298fd71SBarry Smith ierr = ISView(compField, NULL);CHKERRQ(ierr); */ 3507287d2a3SDmitry Karpeev } 351bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 352bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 353bafc1b83SMatthew G Knepley f = ifields[j]; 3547b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 3557b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 3567b62db95SJungho Lee } 357bafc1b83SMatthew G Knepley } 358bafc1b83SMatthew G Knepley if (i == 0) { 359bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 360bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 361bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 362bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 363bafc1b83SMatthew G Knepley } 364bafc1b83SMatthew G Knepley } else { 365d724dfffSBarry Smith for (j=0; j<numFields; j++) { 366d724dfffSBarry Smith ierr = DMDestroy(dms+j);CHKERRQ(ierr); 367d724dfffSBarry Smith } 368d724dfffSBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 369bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 3702fa5cd67SKarl Rupp for (j = 0; j < i; ++j) dms[j] = subdm[j]; 371bafc1b83SMatthew G Knepley } 3727b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3737b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 374e7c4fc90SDmitry Karpeev if (dms) { 3758b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 376bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3777287d2a3SDmitry Karpeev const char *prefix; 3787287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr); 3797287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr); 3807b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3817b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3827287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr); 383e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 3842fa5ba8aSJed Brown } 3857b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3868b8307b2SJed Brown } 38766ffff09SJed Brown } else { 388521d7252SBarry Smith if (jac->bs <= 0) { 389704ba839SBarry Smith if (pc->pmat) { 390521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 3912fa5cd67SKarl Rupp } else jac->bs = 1; 392521d7252SBarry Smith } 393d32f9abdSBarry Smith 3946ce1633cSBarry Smith if (stokes) { 3956ce1633cSBarry Smith IS zerodiags,rest; 3966ce1633cSBarry Smith PetscInt nmin,nmax; 3976ce1633cSBarry Smith 3986ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3996ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 4006ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 4017287d2a3SDmitry Karpeev if (jac->reset) { 4027287d2a3SDmitry Karpeev jac->head->is = rest; 4037287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 4042fa5cd67SKarl Rupp } else { 4056ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 4066ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 4077287d2a3SDmitry Karpeev } 408fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 409fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 4106ce1633cSBarry Smith } else { 411ce94432eSBarry Smith if (jac->reset) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 4129eeaaa73SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr); 4137287d2a3SDmitry Karpeev if (!fieldsplit_default) { 414d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 415d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 4166c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 4176c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 418d32f9abdSBarry Smith } 4197287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 420d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 421db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 4226c924f48SJed Brown char splitname[8]; 4238caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 4245d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 42579416396SBarry Smith } 4265d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 427521d7252SBarry Smith } 42866ffff09SJed Brown } 4296ce1633cSBarry Smith } 430edf189efSBarry Smith } else if (jac->nsplits == 1) { 431edf189efSBarry Smith if (ilink->is) { 432edf189efSBarry Smith IS is2; 433edf189efSBarry Smith PetscInt nmin,nmax; 434edf189efSBarry Smith 435edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 436edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 437db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 438fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 43982f516ccSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 440edf189efSBarry Smith } 441d0af7cd3SBarry Smith 442d0af7cd3SBarry Smith 44382f516ccSBarry Smith if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 44469a612a9SBarry Smith PetscFunctionReturn(0); 44569a612a9SBarry Smith } 44669a612a9SBarry Smith 4475a576424SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 448514bf10dSMatthew G Knepley 44969a612a9SBarry Smith #undef __FUNCT__ 45069a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 45169a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 45269a612a9SBarry Smith { 45369a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 45469a612a9SBarry Smith PetscErrorCode ierr; 4555a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4562c9966d7SBarry Smith PetscInt i,nsplit; 45769a612a9SBarry Smith MatStructure flag = pc->flag; 4585f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 45969a612a9SBarry Smith 46069a612a9SBarry Smith PetscFunctionBegin; 46169a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 46297bbdb24SBarry Smith nsplit = jac->nsplits; 4635a9f2f41SSatish Balay ilink = jac->head; 46497bbdb24SBarry Smith 46597bbdb24SBarry Smith /* get the matrices for each split */ 466704ba839SBarry Smith if (!jac->issetup) { 4671b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 46897bbdb24SBarry Smith 469704ba839SBarry Smith jac->issetup = PETSC_TRUE; 470704ba839SBarry Smith 4715d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4722c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4732c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4742c9966d7SBarry Smith } 47551f519a2SBarry Smith bs = jac->bs; 47697bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4771b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4781b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4791b9fc7fcSBarry Smith if (jac->defaultsplit) { 480ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4815f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 482704ba839SBarry Smith } else if (!ilink->is) { 483ccb205f8SBarry Smith if (ilink->nfields > 1) { 4845f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4855a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4865f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4871b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4881b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4891b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4905f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 49197bbdb24SBarry Smith } 49297bbdb24SBarry Smith } 493ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 494ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 495ccb205f8SBarry Smith } else { 496ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 497ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 498ccb205f8SBarry Smith } 4993e197d65SBarry Smith } 500704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 5015f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 5025f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 503704ba839SBarry Smith ilink = ilink->next; 5041b9fc7fcSBarry Smith } 5051b9fc7fcSBarry Smith } 5061b9fc7fcSBarry Smith 507704ba839SBarry Smith ilink = jac->head; 50897bbdb24SBarry Smith if (!jac->pmat) { 509bdddcaaaSMatthew G Knepley Vec xtmp; 510bdddcaaaSMatthew G Knepley 5110298fd71SBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr); 512cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 513bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 514cf502942SBarry Smith for (i=0; i<nsplit; i++) { 515bdddcaaaSMatthew G Knepley MatNullSpace sp; 516bdddcaaaSMatthew G Knepley 517a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 518a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr); 519a3df900dSMatthew G Knepley if (jac->pmat[i]) { 520a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 521a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 522a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 5232fa5cd67SKarl Rupp 524a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 525a3df900dSMatthew G Knepley } 526a3df900dSMatthew G Knepley } else { 5275f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 528a3df900dSMatthew G Knepley } 529bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 530bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 5310298fd71SBarry Smith ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL; 532bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5330298fd71SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr); 534ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 535bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr); 536bafc1b83SMatthew G Knepley if (sp) { 537bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 538bafc1b83SMatthew G Knepley } 539ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr); 540ed1f0337SMatthew G Knepley if (sp) { 541ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 542ed1f0337SMatthew G Knepley } 543704ba839SBarry Smith ilink = ilink->next; 544cf502942SBarry Smith } 545bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 54697bbdb24SBarry Smith } else { 547cf502942SBarry Smith for (i=0; i<nsplit; i++) { 548a3df900dSMatthew G Knepley Mat pmat; 549a3df900dSMatthew G Knepley 550a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 551a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr); 552a3df900dSMatthew G Knepley if (!pmat) { 5535f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 554a3df900dSMatthew G Knepley } 555704ba839SBarry Smith ilink = ilink->next; 556cf502942SBarry Smith } 55797bbdb24SBarry Smith } 558f5236f50SJed Brown if (pc->useAmat) { 559519d70e2SJed Brown ilink = jac->head; 560519d70e2SJed Brown if (!jac->mat) { 561519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 562519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5635f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 564519d70e2SJed Brown ilink = ilink->next; 565519d70e2SJed Brown } 566519d70e2SJed Brown } else { 567519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5685f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 569519d70e2SJed Brown ilink = ilink->next; 570519d70e2SJed Brown } 571519d70e2SJed Brown } 572519d70e2SJed Brown } else { 573519d70e2SJed Brown jac->mat = jac->pmat; 574519d70e2SJed Brown } 57597bbdb24SBarry Smith 5766c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 57768dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 57868dd23aaSBarry Smith ilink = jac->head; 57968dd23aaSBarry Smith if (!jac->Afield) { 58068dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 58168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5820298fd71SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58368dd23aaSBarry Smith ilink = ilink->next; 58468dd23aaSBarry Smith } 58568dd23aaSBarry Smith } else { 58668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5870298fd71SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58868dd23aaSBarry Smith ilink = ilink->next; 58968dd23aaSBarry Smith } 59068dd23aaSBarry Smith } 59168dd23aaSBarry Smith } 59268dd23aaSBarry Smith 5933b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5943b224e63SBarry Smith IS ccis; 5954aa3045dSJed Brown PetscInt rstart,rend; 596093c86ffSJed Brown char lscname[256]; 597093c86ffSJed Brown PetscObject LSC_L; 598ce94432eSBarry Smith 599ce94432eSBarry Smith if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 60068dd23aaSBarry Smith 601e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 602e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 603e6cab6aaSJed Brown 6043b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 6053b224e63SBarry Smith if (jac->schur) { 6060298fd71SBarry Smith KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper; 607443836d0SMatthew G Knepley 608fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 6093b224e63SBarry Smith ilink = jac->head; 61049bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6114aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 612fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6133b224e63SBarry Smith ilink = ilink->next; 61449bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6154aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 616fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 617a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 618443836d0SMatthew G Knepley if (kspA != kspInner) { 619443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 620443836d0SMatthew G Knepley } 621443836d0SMatthew G Knepley if (kspUpper != kspA) { 622443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 623443836d0SMatthew G Knepley } 624084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 6253b224e63SBarry Smith } else { 626bafc1b83SMatthew G Knepley const char *Dprefix; 627bafc1b83SMatthew G Knepley char schurprefix[256]; 628514bf10dSMatthew G Knepley char schurtestoption[256]; 629bdddcaaaSMatthew G Knepley MatNullSpace sp; 630514bf10dSMatthew G Knepley PetscBool flg; 6313b224e63SBarry Smith 632a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 6333b224e63SBarry Smith ilink = jac->head; 63449bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6354aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 636fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6373b224e63SBarry Smith ilink = ilink->next; 63849bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6394aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 640fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 64120252d06SBarry Smith 642f5236f50SJed Brown /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */ 64320252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 64420252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 64520252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 64620252d06SBarry Smith 647bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 64820252d06SBarry Smith if (sp) { 64920252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 65020252d06SBarry Smith } 65120252d06SBarry Smith 65220252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 6530298fd71SBarry Smith ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr); 654514bf10dSMatthew G Knepley if (flg) { 655514bf10dSMatthew G Knepley DM dmInner; 65621635b76SJed Brown KSP kspInner; 65721635b76SJed Brown 65821635b76SJed Brown ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 65921635b76SJed Brown ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 66021635b76SJed Brown /* Indent this deeper to emphasize the "inner" nature of this solver. */ 66121635b76SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr); 66221635b76SJed Brown ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr); 663514bf10dSMatthew G Knepley 664514bf10dSMatthew G Knepley /* Set DM for new solver */ 665bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 66621635b76SJed Brown ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr); 66721635b76SJed Brown ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr); 668514bf10dSMatthew G Knepley } else { 66921635b76SJed Brown /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or 67021635b76SJed Brown * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact, 67121635b76SJed Brown * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for 67221635b76SJed 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 67321635b76SJed Brown * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used 67421635b76SJed Brown * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */ 67521635b76SJed Brown ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr); 676514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr); 677bafc1b83SMatthew G Knepley } 6785a9f2f41SSatish Balay ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 6795a9f2f41SSatish Balay ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 68097bbdb24SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 68197bbdb24SBarry Smith 682443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 6830298fd71SBarry Smith ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr); 684443836d0SMatthew G Knepley if (flg) { 685443836d0SMatthew G Knepley DM dmInner; 686443836d0SMatthew G Knepley 687443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 68882f516ccSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr); 689443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 690443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 691443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 692443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 693443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 694443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 695443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 696443836d0SMatthew G Knepley } else { 697443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 698443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 699443836d0SMatthew G Knepley } 700443836d0SMatthew G Knepley 701ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr); 70297bbdb24SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 70320252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 70497bbdb24SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 70597bbdb24SBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 7067233a360SDmitry Karpeev PC pcschur; 7077233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 7087233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 70997bbdb24SBarry Smith /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 71097bbdb24SBarry Smith } 71197bbdb24SBarry Smith ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr); 71297bbdb24SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix);CHKERRQ(ierr); 71397bbdb24SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 71497bbdb24SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 71597bbdb24SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 71697bbdb24SBarry Smith } 71797bbdb24SBarry Smith 7185a9f2f41SSatish Balay /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 7198caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 720519d70e2SJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 7213b224e63SBarry Smith if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 722c1570756SJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 7238caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 72497bbdb24SBarry Smith ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 7255a9f2f41SSatish Balay if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 7260971522cSBarry Smith if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 72797bbdb24SBarry Smith } else { 72868bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 7290971522cSBarry Smith i = 0; 7300971522cSBarry Smith ilink = jac->head; 7310971522cSBarry Smith while (ilink) { 7320971522cSBarry Smith ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 7330971522cSBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 7340971522cSBarry Smith if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 7350971522cSBarry Smith i++; 7360971522cSBarry Smith ilink = ilink->next; 7370971522cSBarry Smith } 7383b224e63SBarry Smith } 7393b224e63SBarry Smith 740c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 7410971522cSBarry Smith PetscFunctionReturn(0); 7420971522cSBarry Smith } 7430971522cSBarry Smith 7445a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 745ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 746ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 7475a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 748ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 749ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 75079416396SBarry Smith 7510971522cSBarry Smith #undef __FUNCT__ 7523b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7533b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7543b224e63SBarry Smith { 7553b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7563b224e63SBarry Smith PetscErrorCode ierr; 7573b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 758443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7593b224e63SBarry Smith 7603b224e63SBarry Smith PetscFunctionBegin; 761c5d2311dSJed Brown switch (jac->schurfactorization) { 762c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 763a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 764c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 765c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 766c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 767443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 768c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 769c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 770c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 771c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 772c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 773c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 774c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 775c5d2311dSJed Brown break; 776c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 777a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 778c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 779c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 780443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 781c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 782c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 783c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 784c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 785c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 786c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 787c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 788c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 789c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 790c5d2311dSJed Brown break; 791c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 792a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 793c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 794c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 795c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 796c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 797c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 798c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 799c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 800c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 801443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 802c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 803c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 804c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 805c5d2311dSJed Brown break; 806c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 807a04f6461SBarry 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 */ 8083b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8093b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 810443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 8113b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 8123b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 8133b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8143b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8153b224e63SBarry Smith 8163b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 8173b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8183b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8193b224e63SBarry Smith 820443836d0SMatthew G Knepley if (kspUpper == kspA) { 8213b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 8223b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 823443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 824443836d0SMatthew G Knepley } else { 825443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 826443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 827443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 828443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 829443836d0SMatthew G Knepley } 8303b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8313b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 832c5d2311dSJed Brown } 8333b224e63SBarry Smith PetscFunctionReturn(0); 8343b224e63SBarry Smith } 8353b224e63SBarry Smith 8363b224e63SBarry Smith #undef __FUNCT__ 8370971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 8380971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 8390971522cSBarry Smith { 8400971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8410971522cSBarry Smith PetscErrorCode ierr; 8425a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 843939b8a20SBarry Smith PetscInt cnt,bs; 8440971522cSBarry Smith 8450971522cSBarry Smith PetscFunctionBegin; 84679416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8471b9fc7fcSBarry Smith if (jac->defaultsplit) { 848939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 849ce94432eSBarry 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); 850939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 851ce94432eSBarry 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); 8520971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8535a9f2f41SSatish Balay while (ilink) { 8545a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8555a9f2f41SSatish Balay ilink = ilink->next; 8560971522cSBarry Smith } 8570971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8581b9fc7fcSBarry Smith } else { 859efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8605a9f2f41SSatish Balay while (ilink) { 8615a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8625a9f2f41SSatish Balay ilink = ilink->next; 8631b9fc7fcSBarry Smith } 8641b9fc7fcSBarry Smith } 86516913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 86679416396SBarry Smith if (!jac->w1) { 86779416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 86879416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 86979416396SBarry Smith } 870efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8715a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8723e197d65SBarry Smith cnt = 1; 8735a9f2f41SSatish Balay while (ilink->next) { 8745a9f2f41SSatish Balay ilink = ilink->next; 8753e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8763e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8773e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8783e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8793e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8803e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8813e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8823e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8833e197d65SBarry Smith } 88451f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 88511755939SBarry Smith cnt -= 2; 88651f519a2SBarry Smith while (ilink->previous) { 88751f519a2SBarry Smith ilink = ilink->previous; 88811755939SBarry Smith /* compute the residual only over the part of the vector needed */ 88911755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 89011755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 89111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 89211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 89311755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 89411755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 89511755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 89651f519a2SBarry Smith } 89711755939SBarry Smith } 898ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 8990971522cSBarry Smith PetscFunctionReturn(0); 9000971522cSBarry Smith } 9010971522cSBarry Smith 902421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 903ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 904ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 905421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 906ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 907ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 908421e10b8SBarry Smith 909421e10b8SBarry Smith #undef __FUNCT__ 9108c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 911421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 912421e10b8SBarry Smith { 913421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 914421e10b8SBarry Smith PetscErrorCode ierr; 915421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 916939b8a20SBarry Smith PetscInt bs; 917421e10b8SBarry Smith 918421e10b8SBarry Smith PetscFunctionBegin; 919421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 920421e10b8SBarry Smith if (jac->defaultsplit) { 921939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 922ce94432eSBarry 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); 923939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 924ce94432eSBarry 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); 925421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 926421e10b8SBarry Smith while (ilink) { 927421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 928421e10b8SBarry Smith ilink = ilink->next; 929421e10b8SBarry Smith } 930421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 931421e10b8SBarry Smith } else { 932421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 933421e10b8SBarry Smith while (ilink) { 934421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 935421e10b8SBarry Smith ilink = ilink->next; 936421e10b8SBarry Smith } 937421e10b8SBarry Smith } 938421e10b8SBarry Smith } else { 939421e10b8SBarry Smith if (!jac->w1) { 940421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 941421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 942421e10b8SBarry Smith } 943421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 944421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 945421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 946421e10b8SBarry Smith while (ilink->next) { 947421e10b8SBarry Smith ilink = ilink->next; 9489989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 949421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 950421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 951421e10b8SBarry Smith } 952421e10b8SBarry Smith while (ilink->previous) { 953421e10b8SBarry Smith ilink = ilink->previous; 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 } else { 959421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 960421e10b8SBarry Smith ilink = ilink->next; 961421e10b8SBarry Smith } 962421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 963421e10b8SBarry Smith while (ilink->previous) { 964421e10b8SBarry Smith ilink = ilink->previous; 9659989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 966421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 967421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 968421e10b8SBarry Smith } 969421e10b8SBarry Smith } 970421e10b8SBarry Smith } 971421e10b8SBarry Smith PetscFunctionReturn(0); 972421e10b8SBarry Smith } 973421e10b8SBarry Smith 9740971522cSBarry Smith #undef __FUNCT__ 975574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 976574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9770971522cSBarry Smith { 9780971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9790971522cSBarry Smith PetscErrorCode ierr; 9805a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9810971522cSBarry Smith 9820971522cSBarry Smith PetscFunctionBegin; 9835a9f2f41SSatish Balay while (ilink) { 984574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 985fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 986fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 987443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 988fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 989fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 990b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9915a9f2f41SSatish Balay next = ilink->next; 9925a9f2f41SSatish Balay ilink = next; 9930971522cSBarry Smith } 99405b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 995574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 996574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 997574deadeSBarry Smith } else if (jac->mat) { 9980298fd71SBarry Smith jac->mat = NULL; 999574deadeSBarry Smith } 100097bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 100168dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 10026bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 10036bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 10046bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 10056bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 10066bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 1007d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 10086bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 10096bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 101063ec74ffSBarry Smith jac->reset = PETSC_TRUE; 1011574deadeSBarry Smith PetscFunctionReturn(0); 1012574deadeSBarry Smith } 1013574deadeSBarry Smith 1014574deadeSBarry Smith #undef __FUNCT__ 1015574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1016574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1017574deadeSBarry Smith { 1018574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1019574deadeSBarry Smith PetscErrorCode ierr; 1020574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 1021574deadeSBarry Smith 1022574deadeSBarry Smith PetscFunctionBegin; 1023574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 1024574deadeSBarry Smith while (ilink) { 10256bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 1026574deadeSBarry Smith next = ilink->next; 1027574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 1028574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 10295d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 1030574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 1031574deadeSBarry Smith ilink = next; 1032574deadeSBarry Smith } 1033574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1034c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1035bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr); 1036bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr); 1037bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr); 1038bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr); 1039bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr); 1040bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",NULL);CHKERRQ(ierr); 1041bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr); 10420971522cSBarry Smith PetscFunctionReturn(0); 10430971522cSBarry Smith } 10440971522cSBarry Smith 10450971522cSBarry Smith #undef __FUNCT__ 10460971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10470971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10480971522cSBarry Smith { 10491b9fc7fcSBarry Smith PetscErrorCode ierr; 10506c924f48SJed Brown PetscInt bs; 1051bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10529dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10533b224e63SBarry Smith PCCompositeType ctype; 10541b9fc7fcSBarry Smith 10550971522cSBarry Smith PetscFunctionBegin; 10561b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 10574ab8060aSDmitry Karpeev ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr); 105851f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 105951f519a2SBarry Smith if (flg) { 106051f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 106151f519a2SBarry Smith } 1062704ba839SBarry Smith 10630298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 1064c0adfefeSBarry Smith if (stokes) { 1065c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1066c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1067c0adfefeSBarry Smith } 1068c0adfefeSBarry Smith 10693b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10703b224e63SBarry Smith if (flg) { 10713b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10723b224e63SBarry Smith } 1073c30613efSMatthew Knepley /* Only setup fields once */ 1074c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1075d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1076d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10776c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10786c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1079d32f9abdSBarry Smith } 1080c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1081c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1082c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 10830298fd71SBarry 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); 10840298fd71SBarry 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); 1085c5d2311dSJed Brown } 10861b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10870971522cSBarry Smith PetscFunctionReturn(0); 10880971522cSBarry Smith } 10890971522cSBarry Smith 10900971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10910971522cSBarry Smith 10920971522cSBarry Smith #undef __FUNCT__ 10930971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 10941e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10950971522cSBarry Smith { 109697bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10970971522cSBarry Smith PetscErrorCode ierr; 10985a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 109969a612a9SBarry Smith char prefix[128]; 11005d4c12cdSJungho Lee PetscInt i; 11010971522cSBarry Smith 11020971522cSBarry Smith PetscFunctionBegin; 11036c924f48SJed Brown if (jac->splitdefined) { 11046c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11056c924f48SJed Brown PetscFunctionReturn(0); 11066c924f48SJed Brown } 110751f519a2SBarry Smith for (i=0; i<n; i++) { 1108e32f2f54SBarry 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); 1109e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 111051f519a2SBarry Smith } 1111704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1112a04f6461SBarry Smith if (splitname) { 1113db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1114a04f6461SBarry Smith } else { 1115a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1116a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1117a04f6461SBarry Smith } 1118704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 11195a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 11205d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 11215d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 11222fa5cd67SKarl Rupp 11235a9f2f41SSatish Balay ilink->nfields = n; 11240298fd71SBarry Smith ilink->next = NULL; 1125ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr); 112620252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11275a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11289005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 112969a612a9SBarry Smith 11308caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr); 11315a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11320971522cSBarry Smith 11330971522cSBarry Smith if (!next) { 11345a9f2f41SSatish Balay jac->head = ilink; 11350298fd71SBarry Smith ilink->previous = NULL; 11360971522cSBarry Smith } else { 11370971522cSBarry Smith while (next->next) { 11380971522cSBarry Smith next = next->next; 11390971522cSBarry Smith } 11405a9f2f41SSatish Balay next->next = ilink; 114151f519a2SBarry Smith ilink->previous = next; 11420971522cSBarry Smith } 11430971522cSBarry Smith jac->nsplits++; 11440971522cSBarry Smith PetscFunctionReturn(0); 11450971522cSBarry Smith } 11460971522cSBarry Smith 1147e69d4d44SBarry Smith #undef __FUNCT__ 1148e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11491e6b0712SBarry Smith static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1150e69d4d44SBarry Smith { 1151e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1152e69d4d44SBarry Smith PetscErrorCode ierr; 1153e69d4d44SBarry Smith 1154e69d4d44SBarry Smith PetscFunctionBegin; 1155519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1156e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 11572fa5cd67SKarl Rupp 1158e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 115913e0d083SBarry Smith if (n) *n = jac->nsplits; 1160e69d4d44SBarry Smith PetscFunctionReturn(0); 1161e69d4d44SBarry Smith } 11620971522cSBarry Smith 11630971522cSBarry Smith #undef __FUNCT__ 116469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11651e6b0712SBarry Smith static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11660971522cSBarry Smith { 11670971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11680971522cSBarry Smith PetscErrorCode ierr; 11690971522cSBarry Smith PetscInt cnt = 0; 11705a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11710971522cSBarry Smith 11720971522cSBarry Smith PetscFunctionBegin; 11735d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11745a9f2f41SSatish Balay while (ilink) { 11755a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11765a9f2f41SSatish Balay ilink = ilink->next; 11770971522cSBarry Smith } 11785d480477SMatthew 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); 117913e0d083SBarry Smith if (n) *n = jac->nsplits; 11800971522cSBarry Smith PetscFunctionReturn(0); 11810971522cSBarry Smith } 11820971522cSBarry Smith 1183704ba839SBarry Smith #undef __FUNCT__ 1184704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11851e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1186704ba839SBarry Smith { 1187704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1188704ba839SBarry Smith PetscErrorCode ierr; 1189704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1190704ba839SBarry Smith char prefix[128]; 1191704ba839SBarry Smith 1192704ba839SBarry Smith PetscFunctionBegin; 11936c924f48SJed Brown if (jac->splitdefined) { 11946c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11956c924f48SJed Brown PetscFunctionReturn(0); 11966c924f48SJed Brown } 119716913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1198a04f6461SBarry Smith if (splitname) { 1199db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1200a04f6461SBarry Smith } else { 1201b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1202b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1203a04f6461SBarry Smith } 12041ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1205b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1206b5787286SJed Brown ilink->is = is; 1207b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1208b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1209b5787286SJed Brown ilink->is_col = is; 12100298fd71SBarry Smith ilink->next = NULL; 1211ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr); 121220252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1213704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 12149005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1215704ba839SBarry Smith 12168caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr); 1217704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1218704ba839SBarry Smith 1219704ba839SBarry Smith if (!next) { 1220704ba839SBarry Smith jac->head = ilink; 12210298fd71SBarry Smith ilink->previous = NULL; 1222704ba839SBarry Smith } else { 1223704ba839SBarry Smith while (next->next) { 1224704ba839SBarry Smith next = next->next; 1225704ba839SBarry Smith } 1226704ba839SBarry Smith next->next = ilink; 1227704ba839SBarry Smith ilink->previous = next; 1228704ba839SBarry Smith } 1229704ba839SBarry Smith jac->nsplits++; 1230704ba839SBarry Smith PetscFunctionReturn(0); 1231704ba839SBarry Smith } 1232704ba839SBarry Smith 12330971522cSBarry Smith #undef __FUNCT__ 12340971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12350971522cSBarry Smith /*@ 12360971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12370971522cSBarry Smith 1238ad4df100SBarry Smith Logically Collective on PC 12390971522cSBarry Smith 12400971522cSBarry Smith Input Parameters: 12410971522cSBarry Smith + pc - the preconditioner context 12420298fd71SBarry Smith . splitname - name of this split, if NULL the number of the split is used 12430971522cSBarry Smith . n - the number of fields in this split 1244db4c96c1SJed Brown - fields - the fields in this split 12450971522cSBarry Smith 12460971522cSBarry Smith Level: intermediate 12470971522cSBarry Smith 1248d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1249d32f9abdSBarry Smith 12507287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1251d32f9abdSBarry 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 1252d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1253d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1254d32f9abdSBarry Smith 1255db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1256db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1257db4c96c1SJed Brown 12585d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12595d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12605d4c12cdSJungho Lee available when this routine is called. 12615d4c12cdSJungho Lee 1262d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12630971522cSBarry Smith 12640971522cSBarry Smith @*/ 12655d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12660971522cSBarry Smith { 12674ac538c5SBarry Smith PetscErrorCode ierr; 12680971522cSBarry Smith 12690971522cSBarry Smith PetscFunctionBegin; 12700700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1271db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1272ce94432eSBarry Smith if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1273db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12745d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12750971522cSBarry Smith PetscFunctionReturn(0); 12760971522cSBarry Smith } 12770971522cSBarry Smith 12780971522cSBarry Smith #undef __FUNCT__ 1279704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1280704ba839SBarry Smith /*@ 1281704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1282704ba839SBarry Smith 1283ad4df100SBarry Smith Logically Collective on PC 1284704ba839SBarry Smith 1285704ba839SBarry Smith Input Parameters: 1286704ba839SBarry Smith + pc - the preconditioner context 12870298fd71SBarry Smith . splitname - name of this split, if NULL the number of the split is used 1288db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1289704ba839SBarry Smith 1290d32f9abdSBarry Smith 1291a6ffb8dbSJed Brown Notes: 1292a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1293a6ffb8dbSJed Brown 1294db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1295db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1296d32f9abdSBarry Smith 1297704ba839SBarry Smith Level: intermediate 1298704ba839SBarry Smith 1299704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1300704ba839SBarry Smith 1301704ba839SBarry Smith @*/ 13027087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1303704ba839SBarry Smith { 13044ac538c5SBarry Smith PetscErrorCode ierr; 1305704ba839SBarry Smith 1306704ba839SBarry Smith PetscFunctionBegin; 13070700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13087b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1309db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 13104ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1311704ba839SBarry Smith PetscFunctionReturn(0); 1312704ba839SBarry Smith } 1313704ba839SBarry Smith 1314704ba839SBarry Smith #undef __FUNCT__ 131557a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 131657a9adfeSMatthew G Knepley /*@ 131757a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 131857a9adfeSMatthew G Knepley 131957a9adfeSMatthew G Knepley Logically Collective on PC 132057a9adfeSMatthew G Knepley 132157a9adfeSMatthew G Knepley Input Parameters: 132257a9adfeSMatthew G Knepley + pc - the preconditioner context 132357a9adfeSMatthew G Knepley - splitname - name of this split 132457a9adfeSMatthew G Knepley 132557a9adfeSMatthew G Knepley Output Parameter: 13260298fd71SBarry Smith - is - the index set that defines the vector elements in this field, or NULL if the field is not found 132757a9adfeSMatthew G Knepley 132857a9adfeSMatthew G Knepley Level: intermediate 132957a9adfeSMatthew G Knepley 133057a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 133157a9adfeSMatthew G Knepley 133257a9adfeSMatthew G Knepley @*/ 133357a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 133457a9adfeSMatthew G Knepley { 133557a9adfeSMatthew G Knepley PetscErrorCode ierr; 133657a9adfeSMatthew G Knepley 133757a9adfeSMatthew G Knepley PetscFunctionBegin; 133857a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 133957a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 134057a9adfeSMatthew G Knepley PetscValidPointer(is,3); 134157a9adfeSMatthew G Knepley { 134257a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 134357a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 134457a9adfeSMatthew G Knepley PetscBool found; 134557a9adfeSMatthew G Knepley 13460298fd71SBarry Smith *is = NULL; 134757a9adfeSMatthew G Knepley while (ilink) { 134857a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 134957a9adfeSMatthew G Knepley if (found) { 135057a9adfeSMatthew G Knepley *is = ilink->is; 135157a9adfeSMatthew G Knepley break; 135257a9adfeSMatthew G Knepley } 135357a9adfeSMatthew G Knepley ilink = ilink->next; 135457a9adfeSMatthew G Knepley } 135557a9adfeSMatthew G Knepley } 135657a9adfeSMatthew G Knepley PetscFunctionReturn(0); 135757a9adfeSMatthew G Knepley } 135857a9adfeSMatthew G Knepley 135957a9adfeSMatthew G Knepley #undef __FUNCT__ 136051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 136151f519a2SBarry Smith /*@ 136251f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 136351f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 136451f519a2SBarry Smith 1365ad4df100SBarry Smith Logically Collective on PC 136651f519a2SBarry Smith 136751f519a2SBarry Smith Input Parameters: 136851f519a2SBarry Smith + pc - the preconditioner context 136951f519a2SBarry Smith - bs - the block size 137051f519a2SBarry Smith 137151f519a2SBarry Smith Level: intermediate 137251f519a2SBarry Smith 137351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 137451f519a2SBarry Smith 137551f519a2SBarry Smith @*/ 13767087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 137751f519a2SBarry Smith { 13784ac538c5SBarry Smith PetscErrorCode ierr; 137951f519a2SBarry Smith 138051f519a2SBarry Smith PetscFunctionBegin; 13810700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1382c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13834ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 138451f519a2SBarry Smith PetscFunctionReturn(0); 138551f519a2SBarry Smith } 138651f519a2SBarry Smith 138751f519a2SBarry Smith #undef __FUNCT__ 138869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 13890971522cSBarry Smith /*@C 139069a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 13910971522cSBarry Smith 139269a612a9SBarry Smith Collective on KSP 13930971522cSBarry Smith 13940971522cSBarry Smith Input Parameter: 13950971522cSBarry Smith . pc - the preconditioner context 13960971522cSBarry Smith 13970971522cSBarry Smith Output Parameters: 139813e0d083SBarry Smith + n - the number of splits 139969a612a9SBarry Smith - pc - the array of KSP contexts 14000971522cSBarry Smith 14010971522cSBarry Smith Note: 1402d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1403d32f9abdSBarry Smith (not the KSP just the array that contains them). 14040971522cSBarry Smith 140569a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 14060971522cSBarry Smith 1407196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 14080298fd71SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the 1409196cc216SBarry Smith KSP array must be. 1410196cc216SBarry Smith 1411196cc216SBarry Smith 14120971522cSBarry Smith Level: advanced 14130971522cSBarry Smith 14140971522cSBarry Smith .seealso: PCFIELDSPLIT 14150971522cSBarry Smith @*/ 14167087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 14170971522cSBarry Smith { 14184ac538c5SBarry Smith PetscErrorCode ierr; 14190971522cSBarry Smith 14200971522cSBarry Smith PetscFunctionBegin; 14210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 142213e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14234ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14240971522cSBarry Smith PetscFunctionReturn(0); 14250971522cSBarry Smith } 14260971522cSBarry Smith 1427e69d4d44SBarry Smith #undef __FUNCT__ 1428e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1429e69d4d44SBarry Smith /*@ 1430e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1431a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1432e69d4d44SBarry Smith 1433e69d4d44SBarry Smith Collective on PC 1434e69d4d44SBarry Smith 1435e69d4d44SBarry Smith Input Parameters: 1436e69d4d44SBarry Smith + pc - the preconditioner context 1437e87fab1bSBarry Smith . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default 14380298fd71SBarry Smith - userpre - matrix to use for preconditioning, or NULL 1439084e4875SJed Brown 1440e69d4d44SBarry Smith Options Database: 1441e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> default is a11 1442e69d4d44SBarry Smith 1443fd1303e9SJungho Lee Notes: 1444fd1303e9SJungho Lee $ If ptype is 1445fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1446fd1303e9SJungho Lee $ to this function). 1447e87fab1bSBarry Smith $ a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1448fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1449fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1450fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1451fd1303e9SJungho Lee $ preconditioner 1452fd1303e9SJungho Lee 1453e87fab1bSBarry Smith When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense 1454fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1455fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1456fd1303e9SJungho Lee 1457fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1458fd1303e9SJungho Lee the name since it sets a proceedure to use. 1459fd1303e9SJungho Lee 1460e69d4d44SBarry Smith Level: intermediate 1461e69d4d44SBarry Smith 1462fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1463e69d4d44SBarry Smith 1464e69d4d44SBarry Smith @*/ 14657087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1466e69d4d44SBarry Smith { 14674ac538c5SBarry Smith PetscErrorCode ierr; 1468e69d4d44SBarry Smith 1469e69d4d44SBarry Smith PetscFunctionBegin; 14700700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14714ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1472e69d4d44SBarry Smith PetscFunctionReturn(0); 1473e69d4d44SBarry Smith } 1474e69d4d44SBarry Smith 1475e69d4d44SBarry Smith #undef __FUNCT__ 1476e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14771e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1478e69d4d44SBarry Smith { 1479e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1480084e4875SJed Brown PetscErrorCode ierr; 1481e69d4d44SBarry Smith 1482e69d4d44SBarry Smith PetscFunctionBegin; 1483084e4875SJed Brown jac->schurpre = ptype; 1484084e4875SJed Brown if (pre) { 14856bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1486084e4875SJed Brown jac->schur_user = pre; 1487084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1488084e4875SJed Brown } 1489e69d4d44SBarry Smith PetscFunctionReturn(0); 1490e69d4d44SBarry Smith } 1491e69d4d44SBarry Smith 149230ad9308SMatthew Knepley #undef __FUNCT__ 1493c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1494ab1df9f5SJed Brown /*@ 1495c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1496ab1df9f5SJed Brown 1497ab1df9f5SJed Brown Collective on PC 1498ab1df9f5SJed Brown 1499ab1df9f5SJed Brown Input Parameters: 1500ab1df9f5SJed Brown + pc - the preconditioner context 1501c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1502ab1df9f5SJed Brown 1503ab1df9f5SJed Brown Options Database: 1504c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1505ab1df9f5SJed Brown 1506ab1df9f5SJed Brown 1507ab1df9f5SJed Brown Level: intermediate 1508ab1df9f5SJed Brown 1509ab1df9f5SJed Brown Notes: 1510ab1df9f5SJed Brown The FULL factorization is 1511ab1df9f5SJed Brown 1512ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1513ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1514ab1df9f5SJed Brown 15156be4592eSBarry 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, 15166be4592eSBarry 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). 1517ab1df9f5SJed Brown 15186be4592eSBarry 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 15196be4592eSBarry 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 15206be4592eSBarry 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, 15216be4592eSBarry 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. 1522ab1df9f5SJed Brown 15236be4592eSBarry 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 15246be4592eSBarry 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). 1525ab1df9f5SJed Brown 1526ab1df9f5SJed Brown References: 1527ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1528ab1df9f5SJed Brown 1529ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1530ab1df9f5SJed Brown 1531ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1532ab1df9f5SJed Brown @*/ 1533c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1534ab1df9f5SJed Brown { 1535ab1df9f5SJed Brown PetscErrorCode ierr; 1536ab1df9f5SJed Brown 1537ab1df9f5SJed Brown PetscFunctionBegin; 1538ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1539c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1540ab1df9f5SJed Brown PetscFunctionReturn(0); 1541ab1df9f5SJed Brown } 1542ab1df9f5SJed Brown 1543ab1df9f5SJed Brown #undef __FUNCT__ 1544c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 15451e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1546ab1df9f5SJed Brown { 1547ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1548ab1df9f5SJed Brown 1549ab1df9f5SJed Brown PetscFunctionBegin; 1550ab1df9f5SJed Brown jac->schurfactorization = ftype; 1551ab1df9f5SJed Brown PetscFunctionReturn(0); 1552ab1df9f5SJed Brown } 1553ab1df9f5SJed Brown 1554ab1df9f5SJed Brown #undef __FUNCT__ 155530ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 155630ad9308SMatthew Knepley /*@C 15578c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 155830ad9308SMatthew Knepley 155930ad9308SMatthew Knepley Collective on KSP 156030ad9308SMatthew Knepley 156130ad9308SMatthew Knepley Input Parameter: 156230ad9308SMatthew Knepley . pc - the preconditioner context 156330ad9308SMatthew Knepley 156430ad9308SMatthew Knepley Output Parameters: 1565a04f6461SBarry Smith + A00 - the (0,0) block 1566a04f6461SBarry Smith . A01 - the (0,1) block 1567a04f6461SBarry Smith . A10 - the (1,0) block 1568a04f6461SBarry Smith - A11 - the (1,1) block 156930ad9308SMatthew Knepley 157030ad9308SMatthew Knepley Level: advanced 157130ad9308SMatthew Knepley 157230ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 157330ad9308SMatthew Knepley @*/ 1574a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 157530ad9308SMatthew Knepley { 157630ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 157730ad9308SMatthew Knepley 157830ad9308SMatthew Knepley PetscFunctionBegin; 15790700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1580ce94432eSBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1581a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1582a04f6461SBarry Smith if (A01) *A01 = jac->B; 1583a04f6461SBarry Smith if (A10) *A10 = jac->C; 1584a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 158530ad9308SMatthew Knepley PetscFunctionReturn(0); 158630ad9308SMatthew Knepley } 158730ad9308SMatthew Knepley 158879416396SBarry Smith #undef __FUNCT__ 158979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 15901e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 159179416396SBarry Smith { 159279416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1593e69d4d44SBarry Smith PetscErrorCode ierr; 159479416396SBarry Smith 159579416396SBarry Smith PetscFunctionBegin; 159679416396SBarry Smith jac->type = type; 15973b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 15983b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 15993b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 16002fa5cd67SKarl Rupp 1601bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1602bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1603bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1604e69d4d44SBarry Smith 16053b224e63SBarry Smith } else { 16063b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 16073b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 16082fa5cd67SKarl Rupp 1609bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1610bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",0);CHKERRQ(ierr); 1611bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr); 16123b224e63SBarry Smith } 161379416396SBarry Smith PetscFunctionReturn(0); 161479416396SBarry Smith } 161579416396SBarry Smith 161651f519a2SBarry Smith #undef __FUNCT__ 161751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16181e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 161951f519a2SBarry Smith { 162051f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 162151f519a2SBarry Smith 162251f519a2SBarry Smith PetscFunctionBegin; 1623ce94432eSBarry Smith if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1624ce94432eSBarry 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); 162551f519a2SBarry Smith jac->bs = bs; 162651f519a2SBarry Smith PetscFunctionReturn(0); 162751f519a2SBarry Smith } 162851f519a2SBarry Smith 162979416396SBarry Smith #undef __FUNCT__ 163079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1631bc08b0f1SBarry Smith /*@ 163279416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 163379416396SBarry Smith 163479416396SBarry Smith Collective on PC 163579416396SBarry Smith 163679416396SBarry Smith Input Parameter: 163779416396SBarry Smith . pc - the preconditioner context 163881540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 163979416396SBarry Smith 164079416396SBarry Smith Options Database Key: 1641a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 164279416396SBarry Smith 1643b02e2d75SMatthew G Knepley Level: Intermediate 164479416396SBarry Smith 164579416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 164679416396SBarry Smith 164779416396SBarry Smith .seealso: PCCompositeSetType() 164879416396SBarry Smith 164979416396SBarry Smith @*/ 16507087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 165179416396SBarry Smith { 16524ac538c5SBarry Smith PetscErrorCode ierr; 165379416396SBarry Smith 165479416396SBarry Smith PetscFunctionBegin; 16550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16564ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 165779416396SBarry Smith PetscFunctionReturn(0); 165879416396SBarry Smith } 165979416396SBarry Smith 1660b02e2d75SMatthew G Knepley #undef __FUNCT__ 1661b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1662b02e2d75SMatthew G Knepley /*@ 1663b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1664b02e2d75SMatthew G Knepley 1665b02e2d75SMatthew G Knepley Not collective 1666b02e2d75SMatthew G Knepley 1667b02e2d75SMatthew G Knepley Input Parameter: 1668b02e2d75SMatthew G Knepley . pc - the preconditioner context 1669b02e2d75SMatthew G Knepley 1670b02e2d75SMatthew G Knepley Output Parameter: 1671b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1672b02e2d75SMatthew G Knepley 1673b02e2d75SMatthew G Knepley Level: Intermediate 1674b02e2d75SMatthew G Knepley 1675b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1676b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1677b02e2d75SMatthew G Knepley @*/ 1678b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1679b02e2d75SMatthew G Knepley { 1680b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 1681b02e2d75SMatthew G Knepley 1682b02e2d75SMatthew G Knepley PetscFunctionBegin; 1683b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1684b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1685b02e2d75SMatthew G Knepley *type = jac->type; 1686b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1687b02e2d75SMatthew G Knepley } 1688b02e2d75SMatthew G Knepley 16894ab8060aSDmitry Karpeev #undef __FUNCT__ 16904ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDMSplits" 16914ab8060aSDmitry Karpeev /*@ 16924ab8060aSDmitry Karpeev PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible. 16934ab8060aSDmitry Karpeev 16944ab8060aSDmitry Karpeev Logically Collective 16954ab8060aSDmitry Karpeev 16964ab8060aSDmitry Karpeev Input Parameters: 16974ab8060aSDmitry Karpeev + pc - the preconditioner context 16984ab8060aSDmitry Karpeev - flg - boolean indicating whether to use field splits defined by the DM 16994ab8060aSDmitry Karpeev 17004ab8060aSDmitry Karpeev Options Database Key: 17014ab8060aSDmitry Karpeev . -pc_fieldsplit_dm_splits 17024ab8060aSDmitry Karpeev 17034ab8060aSDmitry Karpeev Level: Intermediate 17044ab8060aSDmitry Karpeev 17054ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative 17064ab8060aSDmitry Karpeev 17074ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits() 17084ab8060aSDmitry Karpeev 17094ab8060aSDmitry Karpeev @*/ 17104ab8060aSDmitry Karpeev PetscErrorCode PCFieldSplitSetDMSplits(PC pc,PetscBool flg) 17114ab8060aSDmitry Karpeev { 17124ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 17134ab8060aSDmitry Karpeev PetscBool isfs; 17144ab8060aSDmitry Karpeev PetscErrorCode ierr; 17154ab8060aSDmitry Karpeev 17164ab8060aSDmitry Karpeev PetscFunctionBegin; 17174ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17184ab8060aSDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 17194ab8060aSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr); 17204ab8060aSDmitry Karpeev if (isfs) { 17214ab8060aSDmitry Karpeev jac->dm_splits = flg; 17224ab8060aSDmitry Karpeev } 17234ab8060aSDmitry Karpeev PetscFunctionReturn(0); 17244ab8060aSDmitry Karpeev } 17254ab8060aSDmitry Karpeev 17264ab8060aSDmitry Karpeev 17274ab8060aSDmitry Karpeev #undef __FUNCT__ 17284ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDMSplits" 17294ab8060aSDmitry Karpeev /*@ 17304ab8060aSDmitry Karpeev PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible. 17314ab8060aSDmitry Karpeev 17324ab8060aSDmitry Karpeev Logically Collective 17334ab8060aSDmitry Karpeev 17344ab8060aSDmitry Karpeev Input Parameter: 17354ab8060aSDmitry Karpeev . pc - the preconditioner context 17364ab8060aSDmitry Karpeev 17374ab8060aSDmitry Karpeev Output Parameter: 17384ab8060aSDmitry Karpeev . flg - boolean indicating whether to use field splits defined by the DM 17394ab8060aSDmitry Karpeev 17404ab8060aSDmitry Karpeev Level: Intermediate 17414ab8060aSDmitry Karpeev 17424ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative 17434ab8060aSDmitry Karpeev 17444ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits() 17454ab8060aSDmitry Karpeev 17464ab8060aSDmitry Karpeev @*/ 17474ab8060aSDmitry Karpeev PetscErrorCode PCFieldSplitGetDMSplits(PC pc,PetscBool* flg) 17484ab8060aSDmitry Karpeev { 17494ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 17504ab8060aSDmitry Karpeev PetscBool isfs; 17514ab8060aSDmitry Karpeev PetscErrorCode ierr; 17524ab8060aSDmitry Karpeev 17534ab8060aSDmitry Karpeev PetscFunctionBegin; 17544ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17554ab8060aSDmitry Karpeev PetscValidPointer(flg,2); 17564ab8060aSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr); 17574ab8060aSDmitry Karpeev if (isfs) { 17584ab8060aSDmitry Karpeev if(flg) *flg = jac->dm_splits; 17594ab8060aSDmitry Karpeev } 17604ab8060aSDmitry Karpeev PetscFunctionReturn(0); 17614ab8060aSDmitry Karpeev } 17624ab8060aSDmitry Karpeev 17630971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 17640971522cSBarry Smith /*MC 1765a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1766a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 17670971522cSBarry Smith 1768edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1769edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 17700971522cSBarry Smith 1771a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 177269a612a9SBarry Smith and set the options directly on the resulting KSP object 17730971522cSBarry Smith 17740971522cSBarry Smith Level: intermediate 17750971522cSBarry Smith 177679416396SBarry Smith Options Database Keys: 177781540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 177881540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 177981540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 178081540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 17810f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 1782e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11 1783435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 178479416396SBarry Smith 17855d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 17865d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 17875d4c12cdSJungho Lee 1788c8a0d604SMatthew G Knepley Notes: 1789c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1790d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1791d32f9abdSBarry Smith 1792d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1793d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1794d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1795d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1796d32f9abdSBarry Smith 1797c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1798c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1799c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1800c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1801c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1802a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1803c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1804c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1805c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1806a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1807c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1808c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 18095668aaf4SBarry Smith diag gives 1810c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1811c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 18125668aaf4SBarry 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 1813c8a0d604SMatthew G Knepley $ ( A00 0 ) 1814c8a0d604SMatthew G Knepley $ ( A10 S ) 1815c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1816c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1817c8a0d604SMatthew G Knepley $ ( 0 S ) 1818c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1819e69d4d44SBarry Smith 1820edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1821edf189efSBarry Smith is used automatically for a second block. 1822edf189efSBarry Smith 1823ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1824ff218e97SBarry Smith Generally it should be used with the AIJ format. 1825ff218e97SBarry Smith 1826ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1827ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1828ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 18290716a85fSBarry Smith 1830a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 18310971522cSBarry Smith 18327e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1833e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 18340971522cSBarry Smith M*/ 18350971522cSBarry Smith 18360971522cSBarry Smith #undef __FUNCT__ 18370971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 18388cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc) 18390971522cSBarry Smith { 18400971522cSBarry Smith PetscErrorCode ierr; 18410971522cSBarry Smith PC_FieldSplit *jac; 18420971522cSBarry Smith 18430971522cSBarry Smith PetscFunctionBegin; 184438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 18452fa5cd67SKarl Rupp 18460971522cSBarry Smith jac->bs = -1; 18470971522cSBarry Smith jac->nsplits = 0; 18483e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1849e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1850c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 1851fbe7908bSJed Brown jac->dm_splits = PETSC_TRUE; 185251f519a2SBarry Smith 18530971522cSBarry Smith pc->data = (void*)jac; 18540971522cSBarry Smith 18550971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1856421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 18570971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1858574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 18590971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 18600971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 18610971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 18620971522cSBarry Smith pc->ops->applyrichardson = 0; 18630971522cSBarry Smith 1864bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1865bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1866bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 1867bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 1868bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 18690971522cSBarry Smith PetscFunctionReturn(0); 18700971522cSBarry Smith } 18710971522cSBarry Smith 18720971522cSBarry Smith 1873a541d17aSBarry Smith 1874