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 81997fe2eSSatish Balay [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 13*e74569cdSMatthew G. Knepley const char *const PCFieldSplitSchurPreTypes[] = {"SELF","A11","USER","FULL","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]; 68*e74569cdSMatthew G. Knepley case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */ 69084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 70084e4875SJed Brown default: 71084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 72084e4875SJed Brown } 73084e4875SJed Brown } 74084e4875SJed Brown 75084e4875SJed Brown 769804daf3SBarry Smith #include <petscdraw.h> 770971522cSBarry Smith #undef __FUNCT__ 780971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 790971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 800971522cSBarry Smith { 810971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 820971522cSBarry Smith PetscErrorCode ierr; 83d9884438SBarry Smith PetscBool iascii,isdraw; 840971522cSBarry Smith PetscInt i,j; 855a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 860971522cSBarry Smith 870971522cSBarry Smith PetscFunctionBegin; 88251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 89d9884438SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 900971522cSBarry Smith if (iascii) { 912c9966d7SBarry Smith if (jac->bs > 0) { 9251f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 932c9966d7SBarry Smith } else { 942c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 952c9966d7SBarry Smith } 96f5236f50SJed Brown if (pc->useAmat) { 97f5236f50SJed Brown ierr = PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr); 98a3df900dSMatthew G Knepley } 9969a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 1000971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1010971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 1021ab39975SBarry Smith if (ilink->fields) { 1030971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 10479416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1055a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 10679416396SBarry Smith if (j > 0) { 10779416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 10879416396SBarry Smith } 1095a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1100971522cSBarry Smith } 1110971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 11279416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1131ab39975SBarry Smith } else { 1141ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1151ab39975SBarry Smith } 1165a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1175a9f2f41SSatish Balay ilink = ilink->next; 1180971522cSBarry Smith } 1190971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1202fa5cd67SKarl Rupp } 1212fa5cd67SKarl Rupp 1222fa5cd67SKarl Rupp if (isdraw) { 123d9884438SBarry Smith PetscDraw draw; 124d9884438SBarry Smith PetscReal x,y,w,wd; 125d9884438SBarry Smith 126d9884438SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 127d9884438SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 128d9884438SBarry Smith w = 2*PetscMin(1.0 - x,x); 129d9884438SBarry Smith wd = w/(jac->nsplits + 1); 130d9884438SBarry Smith x = x - wd*(jac->nsplits-1)/2.0; 131d9884438SBarry Smith for (i=0; i<jac->nsplits; i++) { 132d9884438SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 133d9884438SBarry Smith ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 134d9884438SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 135d9884438SBarry Smith x += wd; 136d9884438SBarry Smith ilink = ilink->next; 137d9884438SBarry Smith } 1380971522cSBarry Smith } 1390971522cSBarry Smith PetscFunctionReturn(0); 1400971522cSBarry Smith } 1410971522cSBarry Smith 1420971522cSBarry Smith #undef __FUNCT__ 1433b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1443b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1453b224e63SBarry Smith { 1463b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1473b224e63SBarry Smith PetscErrorCode ierr; 1484996c5bdSBarry Smith PetscBool iascii,isdraw; 1493b224e63SBarry Smith PetscInt i,j; 1503b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1513b224e63SBarry Smith 1523b224e63SBarry Smith PetscFunctionBegin; 153251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1544996c5bdSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 1553b224e63SBarry Smith if (iascii) { 1562c9966d7SBarry Smith if (jac->bs > 0) { 157c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1582c9966d7SBarry Smith } else { 1592c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1602c9966d7SBarry Smith } 161f5236f50SJed Brown if (pc->useAmat) { 162f5236f50SJed Brown ierr = PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr); 163a3df900dSMatthew G Knepley } 1643e8b8b31SMatthew G Knepley switch (jac->schurpre) { 1653e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1663e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 167e87fab1bSBarry Smith case PC_FIELDSPLIT_SCHUR_PRE_A11: 168e87fab1bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break; 169*e74569cdSMatthew G. Knepley case PC_FIELDSPLIT_SCHUR_PRE_FULL: 170*e74569cdSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the exact Schur complement\n");CHKERRQ(ierr);break; 1713e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1723e8b8b31SMatthew G Knepley if (jac->schur_user) { 1733e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1743e8b8b31SMatthew G Knepley } else { 175e87fab1bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr); 1763e8b8b31SMatthew G Knepley } 1773e8b8b31SMatthew G Knepley break; 1783e8b8b31SMatthew G Knepley default: 17982f516ccSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1803e8b8b31SMatthew G Knepley } 1813b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1823b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1833b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1843b224e63SBarry Smith if (ilink->fields) { 1853b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1863b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1873b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1883b224e63SBarry Smith if (j > 0) { 1893b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1903b224e63SBarry Smith } 1913b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1923b224e63SBarry Smith } 1933b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1943b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1953b224e63SBarry Smith } else { 1963b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1973b224e63SBarry Smith } 1983b224e63SBarry Smith ilink = ilink->next; 1993b224e63SBarry Smith } 200435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");CHKERRQ(ierr); 2013b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 20206de4afeSJed Brown if (jac->head) { 203443836d0SMatthew G Knepley ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 20406de4afeSJed Brown } else {ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr);} 2053b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 20606de4afeSJed Brown if (jac->head && jac->kspupper != jac->head->ksp) { 207443836d0SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr); 208443836d0SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 209b2750c55SJed Brown if (jac->kspupper) {ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);} 210b2750c55SJed Brown else {ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr);} 211443836d0SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 212443836d0SMatthew G Knepley } 213435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 2143b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 21512cae6f2SJed Brown if (jac->kspschur) { 2163b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 21712cae6f2SJed Brown } else { 21812cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 21912cae6f2SJed Brown } 2203b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2213b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 22206de4afeSJed Brown } else if (isdraw && jac->head) { 2234996c5bdSBarry Smith PetscDraw draw; 2244996c5bdSBarry Smith PetscReal x,y,w,wd,h; 2254996c5bdSBarry Smith PetscInt cnt = 2; 2264996c5bdSBarry Smith char str[32]; 2274996c5bdSBarry Smith 2284996c5bdSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 2294996c5bdSBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 230c74581afSBarry Smith if (jac->kspupper != jac->head->ksp) cnt++; 231c74581afSBarry Smith w = 2*PetscMin(1.0 - x,x); 232c74581afSBarry Smith wd = w/(cnt + 1); 233c74581afSBarry Smith 2344996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 2350298fd71SBarry Smith ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 2364996c5bdSBarry Smith y -= h; 2374996c5bdSBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) { 238e87fab1bSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr); 2393b224e63SBarry Smith } else { 2404996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr); 2414996c5bdSBarry Smith } 2420298fd71SBarry Smith ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 2434996c5bdSBarry Smith y -= h; 2444996c5bdSBarry Smith x = x - wd*(cnt-1)/2.0; 2454996c5bdSBarry Smith 2464996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2474996c5bdSBarry Smith ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 2484996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2494996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) { 2504996c5bdSBarry Smith x += wd; 2514996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2524996c5bdSBarry Smith ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 2534996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2544996c5bdSBarry Smith } 2554996c5bdSBarry Smith x += wd; 2564996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2574996c5bdSBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 2584996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2593b224e63SBarry Smith } 2603b224e63SBarry Smith PetscFunctionReturn(0); 2613b224e63SBarry Smith } 2623b224e63SBarry Smith 2633b224e63SBarry Smith #undef __FUNCT__ 2646c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2656c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2666c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2676c924f48SJed Brown { 2686c924f48SJed Brown PetscErrorCode ierr; 2696c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2705d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2715d4c12cdSJungho Lee PetscBool flg,flg_col; 2725d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2736c924f48SJed Brown 2746c924f48SJed Brown PetscFunctionBegin; 275785e854fSJed Brown ierr = PetscMalloc1(jac->bs,&ifields);CHKERRQ(ierr); 276785e854fSJed Brown ierr = PetscMalloc1(jac->bs,&ifields_col);CHKERRQ(ierr); 2776c924f48SJed Brown for (i=0,flg=PETSC_TRUE;; i++) { 2788caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2798caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2808caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2816c924f48SJed Brown nfields = jac->bs; 28229499fbbSJungho Lee nfields_col = jac->bs; 2836c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2845d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2856c924f48SJed Brown if (!flg) break; 2865d4c12cdSJungho Lee else if (flg && !flg_col) { 2875d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2885d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2892fa5cd67SKarl Rupp } else { 2905d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2915d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2925d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2935d4c12cdSJungho Lee } 2946c924f48SJed Brown } 2956c924f48SJed Brown if (i > 0) { 2966c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2976c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2986c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2996c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 3006c924f48SJed Brown } 3016c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 3025d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 3036c924f48SJed Brown PetscFunctionReturn(0); 3046c924f48SJed Brown } 3056c924f48SJed Brown 3066c924f48SJed Brown #undef __FUNCT__ 30769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 30869a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 3090971522cSBarry Smith { 3100971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3110971522cSBarry Smith PetscErrorCode ierr; 3125a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3137287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 3146c924f48SJed Brown PetscInt i; 3150971522cSBarry Smith 3160971522cSBarry Smith PetscFunctionBegin; 3177287d2a3SDmitry Karpeev /* 3187287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 3197287d2a3SDmitry Karpeev Should probably be rewritten. 3207287d2a3SDmitry Karpeev */ 3217287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 3220298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 3234ab8060aSDmitry Karpeev if (pc->dm && jac->dm_splits && !stokes) { 324bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 3250784a22cSJed Brown char **fieldNames; 3267b62db95SJungho Lee IS *fields; 327e7c4fc90SDmitry Karpeev DM *dms; 328bafc1b83SMatthew G Knepley DM subdm[128]; 329bafc1b83SMatthew G Knepley PetscBool flg; 330bafc1b83SMatthew G Knepley 33116621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 332bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 333bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE;; i++) { 334bafc1b83SMatthew G Knepley PetscInt ifields[128]; 335bafc1b83SMatthew G Knepley IS compField; 336bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 337bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 338bafc1b83SMatthew G Knepley 3398caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 340bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 341bafc1b83SMatthew G Knepley if (!flg) break; 34282f516ccSBarry Smith if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 343bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 344bafc1b83SMatthew G Knepley if (nfields == 1) { 345bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 34682f516ccSBarry Smith /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 3470298fd71SBarry Smith ierr = ISView(compField, NULL);CHKERRQ(ierr); */ 348bafc1b83SMatthew G Knepley } else { 3498caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 350bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 35182f516ccSBarry Smith /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", splitname);CHKERRQ(ierr); 3520298fd71SBarry Smith ierr = ISView(compField, NULL);CHKERRQ(ierr); */ 3537287d2a3SDmitry Karpeev } 354bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 355bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 356bafc1b83SMatthew G Knepley f = ifields[j]; 3577b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 3587b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 3597b62db95SJungho Lee } 360bafc1b83SMatthew G Knepley } 361bafc1b83SMatthew G Knepley if (i == 0) { 362bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 363bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 364bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 365bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 366bafc1b83SMatthew G Knepley } 367bafc1b83SMatthew G Knepley } else { 368d724dfffSBarry Smith for (j=0; j<numFields; j++) { 369d724dfffSBarry Smith ierr = DMDestroy(dms+j);CHKERRQ(ierr); 370d724dfffSBarry Smith } 371d724dfffSBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 372785e854fSJed Brown ierr = PetscMalloc1(i, &dms);CHKERRQ(ierr); 3732fa5cd67SKarl Rupp for (j = 0; j < i; ++j) dms[j] = subdm[j]; 374bafc1b83SMatthew G Knepley } 3757b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3767b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 377e7c4fc90SDmitry Karpeev if (dms) { 3788b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 379bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3807287d2a3SDmitry Karpeev const char *prefix; 3817287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr); 3827287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr); 3837b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3847b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3857287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr); 386e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 3872fa5ba8aSJed Brown } 3887b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3898b8307b2SJed Brown } 39066ffff09SJed Brown } else { 391521d7252SBarry Smith if (jac->bs <= 0) { 392704ba839SBarry Smith if (pc->pmat) { 393521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 3942fa5cd67SKarl Rupp } else jac->bs = 1; 395521d7252SBarry Smith } 396d32f9abdSBarry Smith 3976ce1633cSBarry Smith if (stokes) { 3986ce1633cSBarry Smith IS zerodiags,rest; 3996ce1633cSBarry Smith PetscInt nmin,nmax; 4006ce1633cSBarry Smith 4016ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 4026ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 4036ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 4047287d2a3SDmitry Karpeev if (jac->reset) { 4057287d2a3SDmitry Karpeev jac->head->is = rest; 4067287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 4072fa5cd67SKarl Rupp } else { 4086ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 4096ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 4107287d2a3SDmitry Karpeev } 411fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 412fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 4136ce1633cSBarry Smith } else { 414ce94432eSBarry Smith if (jac->reset) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 4159eeaaa73SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr); 4167287d2a3SDmitry Karpeev if (!fieldsplit_default) { 417d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 418d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 4196c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 4206c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 421d32f9abdSBarry Smith } 4227287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 423d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 424db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 4256c924f48SJed Brown char splitname[8]; 4268caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 4275d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 42879416396SBarry Smith } 4295d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 430521d7252SBarry Smith } 43166ffff09SJed Brown } 4326ce1633cSBarry Smith } 433edf189efSBarry Smith } else if (jac->nsplits == 1) { 434edf189efSBarry Smith if (ilink->is) { 435edf189efSBarry Smith IS is2; 436edf189efSBarry Smith PetscInt nmin,nmax; 437edf189efSBarry Smith 438edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 439edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 440db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 441fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 44282f516ccSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 443edf189efSBarry Smith } 444d0af7cd3SBarry Smith 445d0af7cd3SBarry Smith 44682f516ccSBarry Smith if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 44769a612a9SBarry Smith PetscFunctionReturn(0); 44869a612a9SBarry Smith } 44969a612a9SBarry Smith 4505a576424SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 451514bf10dSMatthew G Knepley 45269a612a9SBarry Smith #undef __FUNCT__ 45369a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 45469a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 45569a612a9SBarry Smith { 45669a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 45769a612a9SBarry Smith PetscErrorCode ierr; 4585a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4592c9966d7SBarry Smith PetscInt i,nsplit; 46069a612a9SBarry Smith MatStructure flag = pc->flag; 4615f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 46269a612a9SBarry Smith 46369a612a9SBarry Smith PetscFunctionBegin; 46469a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 46597bbdb24SBarry Smith nsplit = jac->nsplits; 4665a9f2f41SSatish Balay ilink = jac->head; 46797bbdb24SBarry Smith 46897bbdb24SBarry Smith /* get the matrices for each split */ 469704ba839SBarry Smith if (!jac->issetup) { 4701b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 47197bbdb24SBarry Smith 472704ba839SBarry Smith jac->issetup = PETSC_TRUE; 473704ba839SBarry Smith 4745d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4752c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4762c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4772c9966d7SBarry Smith } 47851f519a2SBarry Smith bs = jac->bs; 47997bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4801b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4811b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4821b9fc7fcSBarry Smith if (jac->defaultsplit) { 483ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4845f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 485704ba839SBarry Smith } else if (!ilink->is) { 486ccb205f8SBarry Smith if (ilink->nfields > 1) { 4875f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 488785e854fSJed Brown ierr = PetscMalloc1(ilink->nfields*nslots,&ii);CHKERRQ(ierr); 489785e854fSJed Brown ierr = PetscMalloc1(ilink->nfields*nslots,&jj);CHKERRQ(ierr); 4901b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4911b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4921b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4935f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 49497bbdb24SBarry Smith } 49597bbdb24SBarry Smith } 496ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 497ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 498ccb205f8SBarry Smith } else { 499ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 500ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 501ccb205f8SBarry Smith } 5023e197d65SBarry Smith } 503704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 5045f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 5055f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 506704ba839SBarry Smith ilink = ilink->next; 5071b9fc7fcSBarry Smith } 5081b9fc7fcSBarry Smith } 5091b9fc7fcSBarry Smith 510704ba839SBarry Smith ilink = jac->head; 51197bbdb24SBarry Smith if (!jac->pmat) { 512bdddcaaaSMatthew G Knepley Vec xtmp; 513bdddcaaaSMatthew G Knepley 5140298fd71SBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr); 515785e854fSJed Brown ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr); 516dcca6d9dSJed Brown ierr = PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);CHKERRQ(ierr); 517cf502942SBarry Smith for (i=0; i<nsplit; i++) { 518bdddcaaaSMatthew G Knepley MatNullSpace sp; 519bdddcaaaSMatthew G Knepley 520a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 521a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr); 522a3df900dSMatthew G Knepley if (jac->pmat[i]) { 523a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 524a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 525a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 5262fa5cd67SKarl Rupp 527a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 528a3df900dSMatthew G Knepley } 529a3df900dSMatthew G Knepley } else { 5305f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 531a3df900dSMatthew G Knepley } 532bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 533bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 5340298fd71SBarry Smith ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL; 535bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 5360298fd71SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr); 537ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 538bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr); 539bafc1b83SMatthew G Knepley if (sp) { 540bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 541bafc1b83SMatthew G Knepley } 542ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr); 543ed1f0337SMatthew G Knepley if (sp) { 544ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 545ed1f0337SMatthew G Knepley } 546704ba839SBarry Smith ilink = ilink->next; 547cf502942SBarry Smith } 548bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 54997bbdb24SBarry Smith } else { 550cf502942SBarry Smith for (i=0; i<nsplit; i++) { 551a3df900dSMatthew G Knepley Mat pmat; 552a3df900dSMatthew G Knepley 553a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 554a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr); 555a3df900dSMatthew G Knepley if (!pmat) { 5565f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 557a3df900dSMatthew G Knepley } 558704ba839SBarry Smith ilink = ilink->next; 559cf502942SBarry Smith } 56097bbdb24SBarry Smith } 561f5236f50SJed Brown if (pc->useAmat) { 562519d70e2SJed Brown ilink = jac->head; 563519d70e2SJed Brown if (!jac->mat) { 564785e854fSJed Brown ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr); 565519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5665f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 567519d70e2SJed Brown ilink = ilink->next; 568519d70e2SJed Brown } 569519d70e2SJed Brown } else { 570519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5715f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 572519d70e2SJed Brown ilink = ilink->next; 573519d70e2SJed Brown } 574519d70e2SJed Brown } 575519d70e2SJed Brown } else { 576519d70e2SJed Brown jac->mat = jac->pmat; 577519d70e2SJed Brown } 57897bbdb24SBarry Smith 5796c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 58068dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 58168dd23aaSBarry Smith ilink = jac->head; 58268dd23aaSBarry Smith if (!jac->Afield) { 583785e854fSJed Brown ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr); 58468dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5850298fd71SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58668dd23aaSBarry Smith ilink = ilink->next; 58768dd23aaSBarry Smith } 58868dd23aaSBarry Smith } else { 58968dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5900298fd71SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 59168dd23aaSBarry Smith ilink = ilink->next; 59268dd23aaSBarry Smith } 59368dd23aaSBarry Smith } 59468dd23aaSBarry Smith } 59568dd23aaSBarry Smith 5963b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5973b224e63SBarry Smith IS ccis; 5984aa3045dSJed Brown PetscInt rstart,rend; 599093c86ffSJed Brown char lscname[256]; 600093c86ffSJed Brown PetscObject LSC_L; 601ce94432eSBarry Smith 602ce94432eSBarry Smith if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 60368dd23aaSBarry Smith 604e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 605e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 606e6cab6aaSJed Brown 6073b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 6083b224e63SBarry Smith if (jac->schur) { 6090298fd71SBarry Smith KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper; 610443836d0SMatthew G Knepley 611fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 6123b224e63SBarry Smith ilink = jac->head; 61349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6144aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 615fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6163b224e63SBarry Smith ilink = ilink->next; 61749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6184aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 619fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 620a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 621443836d0SMatthew G Knepley if (kspA != kspInner) { 622443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 623443836d0SMatthew G Knepley } 624443836d0SMatthew G Knepley if (kspUpper != kspA) { 625443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 626443836d0SMatthew G Knepley } 627084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 6283b224e63SBarry Smith } else { 629bafc1b83SMatthew G Knepley const char *Dprefix; 630bafc1b83SMatthew G Knepley char schurprefix[256]; 631514bf10dSMatthew G Knepley char schurtestoption[256]; 632bdddcaaaSMatthew G Knepley MatNullSpace sp; 633514bf10dSMatthew G Knepley PetscBool flg; 6343b224e63SBarry Smith 635a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 6363b224e63SBarry Smith ilink = jac->head; 63749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6384aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 639fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6403b224e63SBarry Smith ilink = ilink->next; 64149bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6424aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 643fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 64420252d06SBarry Smith 645f5236f50SJed Brown /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */ 64620252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 64720252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 64820252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 64920252d06SBarry Smith 650bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 65120252d06SBarry Smith if (sp) { 65220252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 65320252d06SBarry Smith } 65420252d06SBarry Smith 65520252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 6560298fd71SBarry Smith ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr); 657514bf10dSMatthew G Knepley if (flg) { 658514bf10dSMatthew G Knepley DM dmInner; 65921635b76SJed Brown KSP kspInner; 66021635b76SJed Brown 66121635b76SJed Brown ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 66221635b76SJed Brown ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 66321635b76SJed Brown /* Indent this deeper to emphasize the "inner" nature of this solver. */ 66421635b76SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr); 66521635b76SJed Brown ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr); 666514bf10dSMatthew G Knepley 667514bf10dSMatthew G Knepley /* Set DM for new solver */ 668bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 66921635b76SJed Brown ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr); 67021635b76SJed Brown ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr); 671514bf10dSMatthew G Knepley } else { 67221635b76SJed Brown /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or 67321635b76SJed Brown * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact, 67421635b76SJed Brown * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for 67521635b76SJed 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 67621635b76SJed Brown * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used 67721635b76SJed Brown * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */ 67821635b76SJed Brown ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr); 679514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr); 680bafc1b83SMatthew G Knepley } 6815a9f2f41SSatish Balay ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 6825a9f2f41SSatish Balay ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 68397bbdb24SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 68497bbdb24SBarry Smith 685443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 6860298fd71SBarry Smith ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr); 687443836d0SMatthew G Knepley if (flg) { 688443836d0SMatthew G Knepley DM dmInner; 689443836d0SMatthew G Knepley 690443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 69182f516ccSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr); 692443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 693443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 694443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 695443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 696443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 697443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 698443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 699443836d0SMatthew G Knepley } else { 700443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 701443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 702443836d0SMatthew G Knepley } 703443836d0SMatthew G Knepley 704ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr); 70597bbdb24SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 70620252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 70797bbdb24SBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 7087233a360SDmitry Karpeev PC pcschur; 7097233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 7107233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 71197bbdb24SBarry Smith /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 712*e74569cdSMatthew G. Knepley } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) { 713*e74569cdSMatthew G. Knepley ierr = MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);CHKERRQ(ierr); 71497bbdb24SBarry Smith } 715*e74569cdSMatthew G. Knepley ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 71697bbdb24SBarry Smith ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr); 71797bbdb24SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix);CHKERRQ(ierr); 718b20b4189SMatthew G. Knepley /* propogate DM */ 719b20b4189SMatthew G. Knepley { 720b20b4189SMatthew G. Knepley DM sdm; 721b20b4189SMatthew G. Knepley ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr); 722b20b4189SMatthew G. Knepley if (sdm) { 723b20b4189SMatthew G. Knepley ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr); 724b20b4189SMatthew G. Knepley ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr); 725b20b4189SMatthew G. Knepley } 726b20b4189SMatthew G. Knepley } 72797bbdb24SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 72897bbdb24SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 72997bbdb24SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 73097bbdb24SBarry Smith } 73197bbdb24SBarry Smith 7325a9f2f41SSatish Balay /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 7338caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 734519d70e2SJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 7353b224e63SBarry Smith if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 736c1570756SJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 7378caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 73897bbdb24SBarry Smith ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 7395a9f2f41SSatish Balay if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 7400971522cSBarry Smith if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 74197bbdb24SBarry Smith } else { 74268bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 7430971522cSBarry Smith i = 0; 7440971522cSBarry Smith ilink = jac->head; 7450971522cSBarry Smith while (ilink) { 7460971522cSBarry Smith ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 7470971522cSBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 7480971522cSBarry Smith if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 7490971522cSBarry Smith i++; 7500971522cSBarry Smith ilink = ilink->next; 7510971522cSBarry Smith } 7523b224e63SBarry Smith } 7533b224e63SBarry Smith 754c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 7550971522cSBarry Smith PetscFunctionReturn(0); 7560971522cSBarry Smith } 7570971522cSBarry Smith 7585a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 759ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 760ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 7615a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 762ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 763ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 76479416396SBarry Smith 7650971522cSBarry Smith #undef __FUNCT__ 7663b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7673b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7683b224e63SBarry Smith { 7693b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7703b224e63SBarry Smith PetscErrorCode ierr; 7713b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 772443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7733b224e63SBarry Smith 7743b224e63SBarry Smith PetscFunctionBegin; 775c5d2311dSJed Brown switch (jac->schurfactorization) { 776c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 777a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 778c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 779c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 780c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 781443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 782c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 783c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 784c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 785c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 786c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 787c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 788c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 789c5d2311dSJed Brown break; 790c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 791a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 792c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 793c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 794443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 795c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 796c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 797c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 798c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 799c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 800c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 801c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 802c5d2311dSJed Brown ierr = VecScatterEnd(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 break; 805c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 806a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 807c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 808c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 809c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 810c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 811c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 812c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 813c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 814c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 815443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 816c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 817c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 818c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 819c5d2311dSJed Brown break; 820c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 821a04f6461SBarry 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 */ 8223b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8233b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 824443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 8253b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 8263b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 8273b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8283b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8293b224e63SBarry Smith 8303b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 8313b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8323b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8333b224e63SBarry Smith 834443836d0SMatthew G Knepley if (kspUpper == kspA) { 8353b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 8363b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 837443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 838443836d0SMatthew G Knepley } else { 839443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 840443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 841443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 842443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 843443836d0SMatthew G Knepley } 8443b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8453b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 846c5d2311dSJed Brown } 8473b224e63SBarry Smith PetscFunctionReturn(0); 8483b224e63SBarry Smith } 8493b224e63SBarry Smith 8503b224e63SBarry Smith #undef __FUNCT__ 8510971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 8520971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 8530971522cSBarry Smith { 8540971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8550971522cSBarry Smith PetscErrorCode ierr; 8565a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 857939b8a20SBarry Smith PetscInt cnt,bs; 8580971522cSBarry Smith 8590971522cSBarry Smith PetscFunctionBegin; 86079416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8611b9fc7fcSBarry Smith if (jac->defaultsplit) { 862939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 863ce94432eSBarry 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); 864939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 865ce94432eSBarry 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); 8660971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8675a9f2f41SSatish Balay while (ilink) { 8685a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8695a9f2f41SSatish Balay ilink = ilink->next; 8700971522cSBarry Smith } 8710971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8721b9fc7fcSBarry Smith } else { 873efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8745a9f2f41SSatish Balay while (ilink) { 8755a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8765a9f2f41SSatish Balay ilink = ilink->next; 8771b9fc7fcSBarry Smith } 8781b9fc7fcSBarry Smith } 87916913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 88079416396SBarry Smith if (!jac->w1) { 88179416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 88279416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 88379416396SBarry Smith } 884efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8855a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8863e197d65SBarry Smith cnt = 1; 8875a9f2f41SSatish Balay while (ilink->next) { 8885a9f2f41SSatish Balay ilink = ilink->next; 8893e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8903e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8913e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8923e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8933e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8943e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8953e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8963e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8973e197d65SBarry Smith } 89851f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 89911755939SBarry Smith cnt -= 2; 90051f519a2SBarry Smith while (ilink->previous) { 90151f519a2SBarry Smith ilink = ilink->previous; 90211755939SBarry Smith /* compute the residual only over the part of the vector needed */ 90311755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 90411755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 90511755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 90611755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 90711755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 90811755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 90911755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 91051f519a2SBarry Smith } 91111755939SBarry Smith } 912ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 9130971522cSBarry Smith PetscFunctionReturn(0); 9140971522cSBarry Smith } 9150971522cSBarry Smith 916421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 917ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 918ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 919421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 920ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 921ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 922421e10b8SBarry Smith 923421e10b8SBarry Smith #undef __FUNCT__ 9248c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 925421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 926421e10b8SBarry Smith { 927421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 928421e10b8SBarry Smith PetscErrorCode ierr; 929421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 930939b8a20SBarry Smith PetscInt bs; 931421e10b8SBarry Smith 932421e10b8SBarry Smith PetscFunctionBegin; 933421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 934421e10b8SBarry Smith if (jac->defaultsplit) { 935939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 936ce94432eSBarry 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); 937939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 938ce94432eSBarry 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); 939421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 940421e10b8SBarry Smith while (ilink) { 941421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 942421e10b8SBarry Smith ilink = ilink->next; 943421e10b8SBarry Smith } 944421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 945421e10b8SBarry Smith } else { 946421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 947421e10b8SBarry Smith while (ilink) { 948421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 949421e10b8SBarry Smith ilink = ilink->next; 950421e10b8SBarry Smith } 951421e10b8SBarry Smith } 952421e10b8SBarry Smith } else { 953421e10b8SBarry Smith if (!jac->w1) { 954421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 955421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 956421e10b8SBarry Smith } 957421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 958421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 959421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 960421e10b8SBarry Smith while (ilink->next) { 961421e10b8SBarry Smith ilink = ilink->next; 9629989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 963421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 964421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 965421e10b8SBarry Smith } 966421e10b8SBarry Smith while (ilink->previous) { 967421e10b8SBarry Smith ilink = ilink->previous; 9689989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 969421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 970421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 971421e10b8SBarry Smith } 972421e10b8SBarry Smith } else { 973421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 974421e10b8SBarry Smith ilink = ilink->next; 975421e10b8SBarry Smith } 976421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 977421e10b8SBarry Smith while (ilink->previous) { 978421e10b8SBarry Smith ilink = ilink->previous; 9799989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 980421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 981421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 982421e10b8SBarry Smith } 983421e10b8SBarry Smith } 984421e10b8SBarry Smith } 985421e10b8SBarry Smith PetscFunctionReturn(0); 986421e10b8SBarry Smith } 987421e10b8SBarry Smith 9880971522cSBarry Smith #undef __FUNCT__ 989574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 990574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9910971522cSBarry Smith { 9920971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9930971522cSBarry Smith PetscErrorCode ierr; 9945a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9950971522cSBarry Smith 9960971522cSBarry Smith PetscFunctionBegin; 9975a9f2f41SSatish Balay while (ilink) { 998574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 999fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 1000fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 1001443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 1002fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 1003fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1004b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 10055a9f2f41SSatish Balay next = ilink->next; 10065a9f2f41SSatish Balay ilink = next; 10070971522cSBarry Smith } 100805b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1009574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 1010574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 1011574deadeSBarry Smith } else if (jac->mat) { 10120298fd71SBarry Smith jac->mat = NULL; 1013574deadeSBarry Smith } 101497bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 101568dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 10166bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 10176bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 10186bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 10196bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 10206bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 1021d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 10226bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 10236bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 102463ec74ffSBarry Smith jac->reset = PETSC_TRUE; 1025574deadeSBarry Smith PetscFunctionReturn(0); 1026574deadeSBarry Smith } 1027574deadeSBarry Smith 1028574deadeSBarry Smith #undef __FUNCT__ 1029574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1030574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1031574deadeSBarry Smith { 1032574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1033574deadeSBarry Smith PetscErrorCode ierr; 1034574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 1035574deadeSBarry Smith 1036574deadeSBarry Smith PetscFunctionBegin; 1037574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 1038574deadeSBarry Smith while (ilink) { 10396bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 1040574deadeSBarry Smith next = ilink->next; 1041574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 1042574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 10435d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 1044574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 1045574deadeSBarry Smith ilink = next; 1046574deadeSBarry Smith } 1047574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1048c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1049bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr); 1050bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr); 1051bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr); 1052bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr); 1053bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr); 1054bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",NULL);CHKERRQ(ierr); 1055bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr); 10560971522cSBarry Smith PetscFunctionReturn(0); 10570971522cSBarry Smith } 10580971522cSBarry Smith 10590971522cSBarry Smith #undef __FUNCT__ 10600971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10610971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10620971522cSBarry Smith { 10631b9fc7fcSBarry Smith PetscErrorCode ierr; 10646c924f48SJed Brown PetscInt bs; 1065bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10669dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10673b224e63SBarry Smith PCCompositeType ctype; 10681b9fc7fcSBarry Smith 10690971522cSBarry Smith PetscFunctionBegin; 10701b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 10714ab8060aSDmitry Karpeev ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr); 107251f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 107351f519a2SBarry Smith if (flg) { 107451f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 107551f519a2SBarry Smith } 1076704ba839SBarry Smith 10770298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr); 1078c0adfefeSBarry Smith if (stokes) { 1079c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1080c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1081c0adfefeSBarry Smith } 1082c0adfefeSBarry Smith 10833b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10843b224e63SBarry Smith if (flg) { 10853b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10863b224e63SBarry Smith } 1087c30613efSMatthew Knepley /* Only setup fields once */ 1088c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1089d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1090d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10916c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10926c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1093d32f9abdSBarry Smith } 1094c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1095c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1096c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 10970298fd71SBarry 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); 10980298fd71SBarry 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); 1099c5d2311dSJed Brown } 11001b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 11010971522cSBarry Smith PetscFunctionReturn(0); 11020971522cSBarry Smith } 11030971522cSBarry Smith 11040971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 11050971522cSBarry Smith 11060971522cSBarry Smith #undef __FUNCT__ 11070971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 11081e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11090971522cSBarry Smith { 111097bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11110971522cSBarry Smith PetscErrorCode ierr; 11125a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 111369a612a9SBarry Smith char prefix[128]; 11145d4c12cdSJungho Lee PetscInt i; 11150971522cSBarry Smith 11160971522cSBarry Smith PetscFunctionBegin; 11176c924f48SJed Brown if (jac->splitdefined) { 11186c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11196c924f48SJed Brown PetscFunctionReturn(0); 11206c924f48SJed Brown } 112151f519a2SBarry Smith for (i=0; i<n; i++) { 1122e32f2f54SBarry 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); 1123e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 112451f519a2SBarry Smith } 1125b00a9115SJed Brown ierr = PetscNew(&ilink);CHKERRQ(ierr); 1126a04f6461SBarry Smith if (splitname) { 1127db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1128a04f6461SBarry Smith } else { 1129785e854fSJed Brown ierr = PetscMalloc1(3,&ilink->splitname);CHKERRQ(ierr); 1130a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1131a04f6461SBarry Smith } 1132785e854fSJed Brown ierr = PetscMalloc1(n,&ilink->fields);CHKERRQ(ierr); 11335a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 1134785e854fSJed Brown ierr = PetscMalloc1(n,&ilink->fields_col);CHKERRQ(ierr); 11355d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 11362fa5cd67SKarl Rupp 11375a9f2f41SSatish Balay ilink->nfields = n; 11380298fd71SBarry Smith ilink->next = NULL; 1139ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr); 114020252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11415a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11429005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 114369a612a9SBarry Smith 11448caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr); 11455a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11460971522cSBarry Smith 11470971522cSBarry Smith if (!next) { 11485a9f2f41SSatish Balay jac->head = ilink; 11490298fd71SBarry Smith ilink->previous = NULL; 11500971522cSBarry Smith } else { 11510971522cSBarry Smith while (next->next) { 11520971522cSBarry Smith next = next->next; 11530971522cSBarry Smith } 11545a9f2f41SSatish Balay next->next = ilink; 115551f519a2SBarry Smith ilink->previous = next; 11560971522cSBarry Smith } 11570971522cSBarry Smith jac->nsplits++; 11580971522cSBarry Smith PetscFunctionReturn(0); 11590971522cSBarry Smith } 11600971522cSBarry Smith 1161e69d4d44SBarry Smith #undef __FUNCT__ 1162e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11631e6b0712SBarry Smith static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1164e69d4d44SBarry Smith { 1165e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1166e69d4d44SBarry Smith PetscErrorCode ierr; 1167e69d4d44SBarry Smith 1168e69d4d44SBarry Smith PetscFunctionBegin; 1169785e854fSJed Brown ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr); 1170e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 11712fa5cd67SKarl Rupp 1172e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 117313e0d083SBarry Smith if (n) *n = jac->nsplits; 1174e69d4d44SBarry Smith PetscFunctionReturn(0); 1175e69d4d44SBarry Smith } 11760971522cSBarry Smith 11770971522cSBarry Smith #undef __FUNCT__ 117869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11791e6b0712SBarry Smith static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11800971522cSBarry Smith { 11810971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11820971522cSBarry Smith PetscErrorCode ierr; 11830971522cSBarry Smith PetscInt cnt = 0; 11845a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11850971522cSBarry Smith 11860971522cSBarry Smith PetscFunctionBegin; 1187785e854fSJed Brown ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr); 11885a9f2f41SSatish Balay while (ilink) { 11895a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11905a9f2f41SSatish Balay ilink = ilink->next; 11910971522cSBarry Smith } 11925d480477SMatthew 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); 119313e0d083SBarry Smith if (n) *n = jac->nsplits; 11940971522cSBarry Smith PetscFunctionReturn(0); 11950971522cSBarry Smith } 11960971522cSBarry Smith 1197704ba839SBarry Smith #undef __FUNCT__ 1198704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11991e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1200704ba839SBarry Smith { 1201704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1202704ba839SBarry Smith PetscErrorCode ierr; 1203704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1204704ba839SBarry Smith char prefix[128]; 1205704ba839SBarry Smith 1206704ba839SBarry Smith PetscFunctionBegin; 12076c924f48SJed Brown if (jac->splitdefined) { 12086c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 12096c924f48SJed Brown PetscFunctionReturn(0); 12106c924f48SJed Brown } 1211b00a9115SJed Brown ierr = PetscNew(&ilink);CHKERRQ(ierr); 1212a04f6461SBarry Smith if (splitname) { 1213db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1214a04f6461SBarry Smith } else { 1215785e854fSJed Brown ierr = PetscMalloc1(8,&ilink->splitname);CHKERRQ(ierr); 1216b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1217a04f6461SBarry Smith } 12181ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1219b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1220b5787286SJed Brown ilink->is = is; 1221b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1222b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1223b5787286SJed Brown ilink->is_col = is; 12240298fd71SBarry Smith ilink->next = NULL; 1225ce94432eSBarry Smith ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr); 122620252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1227704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 12289005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1229704ba839SBarry Smith 12308caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr); 1231704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1232704ba839SBarry Smith 1233704ba839SBarry Smith if (!next) { 1234704ba839SBarry Smith jac->head = ilink; 12350298fd71SBarry Smith ilink->previous = NULL; 1236704ba839SBarry Smith } else { 1237704ba839SBarry Smith while (next->next) { 1238704ba839SBarry Smith next = next->next; 1239704ba839SBarry Smith } 1240704ba839SBarry Smith next->next = ilink; 1241704ba839SBarry Smith ilink->previous = next; 1242704ba839SBarry Smith } 1243704ba839SBarry Smith jac->nsplits++; 1244704ba839SBarry Smith PetscFunctionReturn(0); 1245704ba839SBarry Smith } 1246704ba839SBarry Smith 12470971522cSBarry Smith #undef __FUNCT__ 12480971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12490971522cSBarry Smith /*@ 12500971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12510971522cSBarry Smith 1252ad4df100SBarry Smith Logically Collective on PC 12530971522cSBarry Smith 12540971522cSBarry Smith Input Parameters: 12550971522cSBarry Smith + pc - the preconditioner context 12560298fd71SBarry Smith . splitname - name of this split, if NULL the number of the split is used 12570971522cSBarry Smith . n - the number of fields in this split 1258db4c96c1SJed Brown - fields - the fields in this split 12590971522cSBarry Smith 12600971522cSBarry Smith Level: intermediate 12610971522cSBarry Smith 1262d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1263d32f9abdSBarry Smith 12647287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1265d32f9abdSBarry 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 1266d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1267d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1268d32f9abdSBarry Smith 1269db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1270db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1271db4c96c1SJed Brown 12725d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12735d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12745d4c12cdSJungho Lee available when this routine is called. 12755d4c12cdSJungho Lee 1276d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12770971522cSBarry Smith 12780971522cSBarry Smith @*/ 12795d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12800971522cSBarry Smith { 12814ac538c5SBarry Smith PetscErrorCode ierr; 12820971522cSBarry Smith 12830971522cSBarry Smith PetscFunctionBegin; 12840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1285db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1286ce94432eSBarry Smith if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1287db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12885d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12890971522cSBarry Smith PetscFunctionReturn(0); 12900971522cSBarry Smith } 12910971522cSBarry Smith 12920971522cSBarry Smith #undef __FUNCT__ 1293704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1294704ba839SBarry Smith /*@ 1295704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1296704ba839SBarry Smith 1297ad4df100SBarry Smith Logically Collective on PC 1298704ba839SBarry Smith 1299704ba839SBarry Smith Input Parameters: 1300704ba839SBarry Smith + pc - the preconditioner context 13010298fd71SBarry Smith . splitname - name of this split, if NULL the number of the split is used 1302db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1303704ba839SBarry Smith 1304d32f9abdSBarry Smith 1305a6ffb8dbSJed Brown Notes: 1306a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1307a6ffb8dbSJed Brown 1308db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1309db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1310d32f9abdSBarry Smith 1311704ba839SBarry Smith Level: intermediate 1312704ba839SBarry Smith 1313704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1314704ba839SBarry Smith 1315704ba839SBarry Smith @*/ 13167087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1317704ba839SBarry Smith { 13184ac538c5SBarry Smith PetscErrorCode ierr; 1319704ba839SBarry Smith 1320704ba839SBarry Smith PetscFunctionBegin; 13210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13227b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1323db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 13244ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1325704ba839SBarry Smith PetscFunctionReturn(0); 1326704ba839SBarry Smith } 1327704ba839SBarry Smith 1328704ba839SBarry Smith #undef __FUNCT__ 132957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 133057a9adfeSMatthew G Knepley /*@ 133157a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 133257a9adfeSMatthew G Knepley 133357a9adfeSMatthew G Knepley Logically Collective on PC 133457a9adfeSMatthew G Knepley 133557a9adfeSMatthew G Knepley Input Parameters: 133657a9adfeSMatthew G Knepley + pc - the preconditioner context 133757a9adfeSMatthew G Knepley - splitname - name of this split 133857a9adfeSMatthew G Knepley 133957a9adfeSMatthew G Knepley Output Parameter: 13400298fd71SBarry Smith - is - the index set that defines the vector elements in this field, or NULL if the field is not found 134157a9adfeSMatthew G Knepley 134257a9adfeSMatthew G Knepley Level: intermediate 134357a9adfeSMatthew G Knepley 134457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 134557a9adfeSMatthew G Knepley 134657a9adfeSMatthew G Knepley @*/ 134757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 134857a9adfeSMatthew G Knepley { 134957a9adfeSMatthew G Knepley PetscErrorCode ierr; 135057a9adfeSMatthew G Knepley 135157a9adfeSMatthew G Knepley PetscFunctionBegin; 135257a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 135357a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 135457a9adfeSMatthew G Knepley PetscValidPointer(is,3); 135557a9adfeSMatthew G Knepley { 135657a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 135757a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 135857a9adfeSMatthew G Knepley PetscBool found; 135957a9adfeSMatthew G Knepley 13600298fd71SBarry Smith *is = NULL; 136157a9adfeSMatthew G Knepley while (ilink) { 136257a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 136357a9adfeSMatthew G Knepley if (found) { 136457a9adfeSMatthew G Knepley *is = ilink->is; 136557a9adfeSMatthew G Knepley break; 136657a9adfeSMatthew G Knepley } 136757a9adfeSMatthew G Knepley ilink = ilink->next; 136857a9adfeSMatthew G Knepley } 136957a9adfeSMatthew G Knepley } 137057a9adfeSMatthew G Knepley PetscFunctionReturn(0); 137157a9adfeSMatthew G Knepley } 137257a9adfeSMatthew G Knepley 137357a9adfeSMatthew G Knepley #undef __FUNCT__ 137451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 137551f519a2SBarry Smith /*@ 137651f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 137751f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 137851f519a2SBarry Smith 1379ad4df100SBarry Smith Logically Collective on PC 138051f519a2SBarry Smith 138151f519a2SBarry Smith Input Parameters: 138251f519a2SBarry Smith + pc - the preconditioner context 138351f519a2SBarry Smith - bs - the block size 138451f519a2SBarry Smith 138551f519a2SBarry Smith Level: intermediate 138651f519a2SBarry Smith 138751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 138851f519a2SBarry Smith 138951f519a2SBarry Smith @*/ 13907087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 139151f519a2SBarry Smith { 13924ac538c5SBarry Smith PetscErrorCode ierr; 139351f519a2SBarry Smith 139451f519a2SBarry Smith PetscFunctionBegin; 13950700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1396c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13974ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 139851f519a2SBarry Smith PetscFunctionReturn(0); 139951f519a2SBarry Smith } 140051f519a2SBarry Smith 140151f519a2SBarry Smith #undef __FUNCT__ 140269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 14030971522cSBarry Smith /*@C 140469a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 14050971522cSBarry Smith 140669a612a9SBarry Smith Collective on KSP 14070971522cSBarry Smith 14080971522cSBarry Smith Input Parameter: 14090971522cSBarry Smith . pc - the preconditioner context 14100971522cSBarry Smith 14110971522cSBarry Smith Output Parameters: 141213e0d083SBarry Smith + n - the number of splits 141369a612a9SBarry Smith - pc - the array of KSP contexts 14140971522cSBarry Smith 14150971522cSBarry Smith Note: 1416d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1417d32f9abdSBarry Smith (not the KSP just the array that contains them). 14180971522cSBarry Smith 141969a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 14200971522cSBarry Smith 1421196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 14220298fd71SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the 1423196cc216SBarry Smith KSP array must be. 1424196cc216SBarry Smith 1425196cc216SBarry Smith 14260971522cSBarry Smith Level: advanced 14270971522cSBarry Smith 14280971522cSBarry Smith .seealso: PCFIELDSPLIT 14290971522cSBarry Smith @*/ 14307087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 14310971522cSBarry Smith { 14324ac538c5SBarry Smith PetscErrorCode ierr; 14330971522cSBarry Smith 14340971522cSBarry Smith PetscFunctionBegin; 14350700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 143613e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14374ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14380971522cSBarry Smith PetscFunctionReturn(0); 14390971522cSBarry Smith } 14400971522cSBarry Smith 1441e69d4d44SBarry Smith #undef __FUNCT__ 1442e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1443e69d4d44SBarry Smith /*@ 1444e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1445a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1446e69d4d44SBarry Smith 1447e69d4d44SBarry Smith Collective on PC 1448e69d4d44SBarry Smith 1449e69d4d44SBarry Smith Input Parameters: 1450e69d4d44SBarry Smith + pc - the preconditioner context 1451e87fab1bSBarry Smith . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default 14520298fd71SBarry Smith - userpre - matrix to use for preconditioning, or NULL 1453084e4875SJed Brown 1454e69d4d44SBarry Smith Options Database: 1455*e74569cdSMatthew G. Knepley . -pc_fieldsplit_schur_precondition <self,user,a11,full> default is a11 1456e69d4d44SBarry Smith 1457fd1303e9SJungho Lee Notes: 1458fd1303e9SJungho Lee $ If ptype is 1459fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1460fd1303e9SJungho Lee $ to this function). 1461e87fab1bSBarry Smith $ a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1462fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1463*e74569cdSMatthew G. Knepley $ full then the preconditioner uses the exact Schur complement (this is expensive) 1464fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1465fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1466fd1303e9SJungho Lee $ preconditioner 1467fd1303e9SJungho Lee 1468e87fab1bSBarry Smith When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense 1469fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1470fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1471fd1303e9SJungho Lee 1472fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1473fd1303e9SJungho Lee the name since it sets a proceedure to use. 1474fd1303e9SJungho Lee 1475e69d4d44SBarry Smith Level: intermediate 1476e69d4d44SBarry Smith 1477fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1478e69d4d44SBarry Smith 1479e69d4d44SBarry Smith @*/ 14807087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1481e69d4d44SBarry Smith { 14824ac538c5SBarry Smith PetscErrorCode ierr; 1483e69d4d44SBarry Smith 1484e69d4d44SBarry Smith PetscFunctionBegin; 14850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14864ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1487e69d4d44SBarry Smith PetscFunctionReturn(0); 1488e69d4d44SBarry Smith } 1489e69d4d44SBarry Smith 1490e69d4d44SBarry Smith #undef __FUNCT__ 1491e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14921e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1493e69d4d44SBarry Smith { 1494e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1495084e4875SJed Brown PetscErrorCode ierr; 1496e69d4d44SBarry Smith 1497e69d4d44SBarry Smith PetscFunctionBegin; 1498084e4875SJed Brown jac->schurpre = ptype; 1499084e4875SJed Brown if (pre) { 15006bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1501084e4875SJed Brown jac->schur_user = pre; 1502084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1503084e4875SJed Brown } 1504e69d4d44SBarry Smith PetscFunctionReturn(0); 1505e69d4d44SBarry Smith } 1506e69d4d44SBarry Smith 150730ad9308SMatthew Knepley #undef __FUNCT__ 1508c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1509ab1df9f5SJed Brown /*@ 1510c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1511ab1df9f5SJed Brown 1512ab1df9f5SJed Brown Collective on PC 1513ab1df9f5SJed Brown 1514ab1df9f5SJed Brown Input Parameters: 1515ab1df9f5SJed Brown + pc - the preconditioner context 1516c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1517ab1df9f5SJed Brown 1518ab1df9f5SJed Brown Options Database: 1519c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1520ab1df9f5SJed Brown 1521ab1df9f5SJed Brown 1522ab1df9f5SJed Brown Level: intermediate 1523ab1df9f5SJed Brown 1524ab1df9f5SJed Brown Notes: 1525ab1df9f5SJed Brown The FULL factorization is 1526ab1df9f5SJed Brown 1527ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1528ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1529ab1df9f5SJed Brown 15306be4592eSBarry 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, 15316be4592eSBarry 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). 1532ab1df9f5SJed Brown 15336be4592eSBarry 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 15346be4592eSBarry 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 15356be4592eSBarry 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, 15366be4592eSBarry 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. 1537ab1df9f5SJed Brown 15386be4592eSBarry 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 15396be4592eSBarry 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). 1540ab1df9f5SJed Brown 1541ab1df9f5SJed Brown References: 1542ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1543ab1df9f5SJed Brown 1544ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1545ab1df9f5SJed Brown 1546ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1547ab1df9f5SJed Brown @*/ 1548c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1549ab1df9f5SJed Brown { 1550ab1df9f5SJed Brown PetscErrorCode ierr; 1551ab1df9f5SJed Brown 1552ab1df9f5SJed Brown PetscFunctionBegin; 1553ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1554c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1555ab1df9f5SJed Brown PetscFunctionReturn(0); 1556ab1df9f5SJed Brown } 1557ab1df9f5SJed Brown 1558ab1df9f5SJed Brown #undef __FUNCT__ 1559c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 15601e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1561ab1df9f5SJed Brown { 1562ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1563ab1df9f5SJed Brown 1564ab1df9f5SJed Brown PetscFunctionBegin; 1565ab1df9f5SJed Brown jac->schurfactorization = ftype; 1566ab1df9f5SJed Brown PetscFunctionReturn(0); 1567ab1df9f5SJed Brown } 1568ab1df9f5SJed Brown 1569ab1df9f5SJed Brown #undef __FUNCT__ 157030ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 157130ad9308SMatthew Knepley /*@C 15728c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 157330ad9308SMatthew Knepley 157430ad9308SMatthew Knepley Collective on KSP 157530ad9308SMatthew Knepley 157630ad9308SMatthew Knepley Input Parameter: 157730ad9308SMatthew Knepley . pc - the preconditioner context 157830ad9308SMatthew Knepley 157930ad9308SMatthew Knepley Output Parameters: 1580a04f6461SBarry Smith + A00 - the (0,0) block 1581a04f6461SBarry Smith . A01 - the (0,1) block 1582a04f6461SBarry Smith . A10 - the (1,0) block 1583a04f6461SBarry Smith - A11 - the (1,1) block 158430ad9308SMatthew Knepley 158530ad9308SMatthew Knepley Level: advanced 158630ad9308SMatthew Knepley 158730ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 158830ad9308SMatthew Knepley @*/ 1589a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 159030ad9308SMatthew Knepley { 159130ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 159230ad9308SMatthew Knepley 159330ad9308SMatthew Knepley PetscFunctionBegin; 15940700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1595ce94432eSBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1596a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1597a04f6461SBarry Smith if (A01) *A01 = jac->B; 1598a04f6461SBarry Smith if (A10) *A10 = jac->C; 1599a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 160030ad9308SMatthew Knepley PetscFunctionReturn(0); 160130ad9308SMatthew Knepley } 160230ad9308SMatthew Knepley 160379416396SBarry Smith #undef __FUNCT__ 160479416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 16051e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 160679416396SBarry Smith { 160779416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1608e69d4d44SBarry Smith PetscErrorCode ierr; 160979416396SBarry Smith 161079416396SBarry Smith PetscFunctionBegin; 161179416396SBarry Smith jac->type = type; 16123b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 16133b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 16143b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 16152fa5cd67SKarl Rupp 1616bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1617bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1618bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1619e69d4d44SBarry Smith 16203b224e63SBarry Smith } else { 16213b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 16223b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 16232fa5cd67SKarl Rupp 1624bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1625bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",0);CHKERRQ(ierr); 1626bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr); 16273b224e63SBarry Smith } 162879416396SBarry Smith PetscFunctionReturn(0); 162979416396SBarry Smith } 163079416396SBarry Smith 163151f519a2SBarry Smith #undef __FUNCT__ 163251f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16331e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 163451f519a2SBarry Smith { 163551f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 163651f519a2SBarry Smith 163751f519a2SBarry Smith PetscFunctionBegin; 1638ce94432eSBarry Smith if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1639ce94432eSBarry 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); 164051f519a2SBarry Smith jac->bs = bs; 164151f519a2SBarry Smith PetscFunctionReturn(0); 164251f519a2SBarry Smith } 164351f519a2SBarry Smith 164479416396SBarry Smith #undef __FUNCT__ 164579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1646bc08b0f1SBarry Smith /*@ 164779416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 164879416396SBarry Smith 164979416396SBarry Smith Collective on PC 165079416396SBarry Smith 165179416396SBarry Smith Input Parameter: 165279416396SBarry Smith . pc - the preconditioner context 165381540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 165479416396SBarry Smith 165579416396SBarry Smith Options Database Key: 1656a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 165779416396SBarry Smith 1658b02e2d75SMatthew G Knepley Level: Intermediate 165979416396SBarry Smith 166079416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 166179416396SBarry Smith 166279416396SBarry Smith .seealso: PCCompositeSetType() 166379416396SBarry Smith 166479416396SBarry Smith @*/ 16657087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 166679416396SBarry Smith { 16674ac538c5SBarry Smith PetscErrorCode ierr; 166879416396SBarry Smith 166979416396SBarry Smith PetscFunctionBegin; 16700700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16714ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 167279416396SBarry Smith PetscFunctionReturn(0); 167379416396SBarry Smith } 167479416396SBarry Smith 1675b02e2d75SMatthew G Knepley #undef __FUNCT__ 1676b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1677b02e2d75SMatthew G Knepley /*@ 1678b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1679b02e2d75SMatthew G Knepley 1680b02e2d75SMatthew G Knepley Not collective 1681b02e2d75SMatthew G Knepley 1682b02e2d75SMatthew G Knepley Input Parameter: 1683b02e2d75SMatthew G Knepley . pc - the preconditioner context 1684b02e2d75SMatthew G Knepley 1685b02e2d75SMatthew G Knepley Output Parameter: 1686b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1687b02e2d75SMatthew G Knepley 1688b02e2d75SMatthew G Knepley Level: Intermediate 1689b02e2d75SMatthew G Knepley 1690b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1691b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1692b02e2d75SMatthew G Knepley @*/ 1693b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1694b02e2d75SMatthew G Knepley { 1695b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit*) pc->data; 1696b02e2d75SMatthew G Knepley 1697b02e2d75SMatthew G Knepley PetscFunctionBegin; 1698b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1699b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1700b02e2d75SMatthew G Knepley *type = jac->type; 1701b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1702b02e2d75SMatthew G Knepley } 1703b02e2d75SMatthew G Knepley 17044ab8060aSDmitry Karpeev #undef __FUNCT__ 17054ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDMSplits" 17064ab8060aSDmitry Karpeev /*@ 17074ab8060aSDmitry Karpeev PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible. 17084ab8060aSDmitry Karpeev 17094ab8060aSDmitry Karpeev Logically Collective 17104ab8060aSDmitry Karpeev 17114ab8060aSDmitry Karpeev Input Parameters: 17124ab8060aSDmitry Karpeev + pc - the preconditioner context 17134ab8060aSDmitry Karpeev - flg - boolean indicating whether to use field splits defined by the DM 17144ab8060aSDmitry Karpeev 17154ab8060aSDmitry Karpeev Options Database Key: 17164ab8060aSDmitry Karpeev . -pc_fieldsplit_dm_splits 17174ab8060aSDmitry Karpeev 17184ab8060aSDmitry Karpeev Level: Intermediate 17194ab8060aSDmitry Karpeev 17204ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative 17214ab8060aSDmitry Karpeev 17224ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits() 17234ab8060aSDmitry Karpeev 17244ab8060aSDmitry Karpeev @*/ 17254ab8060aSDmitry Karpeev PetscErrorCode PCFieldSplitSetDMSplits(PC pc,PetscBool flg) 17264ab8060aSDmitry Karpeev { 17274ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 17284ab8060aSDmitry Karpeev PetscBool isfs; 17294ab8060aSDmitry Karpeev PetscErrorCode ierr; 17304ab8060aSDmitry Karpeev 17314ab8060aSDmitry Karpeev PetscFunctionBegin; 17324ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17334ab8060aSDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 17344ab8060aSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr); 17354ab8060aSDmitry Karpeev if (isfs) { 17364ab8060aSDmitry Karpeev jac->dm_splits = flg; 17374ab8060aSDmitry Karpeev } 17384ab8060aSDmitry Karpeev PetscFunctionReturn(0); 17394ab8060aSDmitry Karpeev } 17404ab8060aSDmitry Karpeev 17414ab8060aSDmitry Karpeev 17424ab8060aSDmitry Karpeev #undef __FUNCT__ 17434ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDMSplits" 17444ab8060aSDmitry Karpeev /*@ 17454ab8060aSDmitry Karpeev PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible. 17464ab8060aSDmitry Karpeev 17474ab8060aSDmitry Karpeev Logically Collective 17484ab8060aSDmitry Karpeev 17494ab8060aSDmitry Karpeev Input Parameter: 17504ab8060aSDmitry Karpeev . pc - the preconditioner context 17514ab8060aSDmitry Karpeev 17524ab8060aSDmitry Karpeev Output Parameter: 17534ab8060aSDmitry Karpeev . flg - boolean indicating whether to use field splits defined by the DM 17544ab8060aSDmitry Karpeev 17554ab8060aSDmitry Karpeev Level: Intermediate 17564ab8060aSDmitry Karpeev 17574ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative 17584ab8060aSDmitry Karpeev 17594ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits() 17604ab8060aSDmitry Karpeev 17614ab8060aSDmitry Karpeev @*/ 17624ab8060aSDmitry Karpeev PetscErrorCode PCFieldSplitGetDMSplits(PC pc,PetscBool* flg) 17634ab8060aSDmitry Karpeev { 17644ab8060aSDmitry Karpeev PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 17654ab8060aSDmitry Karpeev PetscBool isfs; 17664ab8060aSDmitry Karpeev PetscErrorCode ierr; 17674ab8060aSDmitry Karpeev 17684ab8060aSDmitry Karpeev PetscFunctionBegin; 17694ab8060aSDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17704ab8060aSDmitry Karpeev PetscValidPointer(flg,2); 17714ab8060aSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr); 17724ab8060aSDmitry Karpeev if (isfs) { 17734ab8060aSDmitry Karpeev if(flg) *flg = jac->dm_splits; 17744ab8060aSDmitry Karpeev } 17754ab8060aSDmitry Karpeev PetscFunctionReturn(0); 17764ab8060aSDmitry Karpeev } 17774ab8060aSDmitry Karpeev 17780971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 17790971522cSBarry Smith /*MC 1780a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1781a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 17820971522cSBarry Smith 1783edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1784edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 17850971522cSBarry Smith 1786a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 178769a612a9SBarry Smith and set the options directly on the resulting KSP object 17880971522cSBarry Smith 17890971522cSBarry Smith Level: intermediate 17900971522cSBarry Smith 179179416396SBarry Smith Options Database Keys: 179281540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 179381540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 179481540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 179581540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 17960f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 1797*e74569cdSMatthew G. Knepley . -pc_fieldsplit_schur_precondition <self,user,a11,full> - default is a11 1798435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 179979416396SBarry Smith 18005d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 18015d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 18025d4c12cdSJungho Lee 1803c8a0d604SMatthew G Knepley Notes: 1804c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1805d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1806d32f9abdSBarry Smith 1807d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1808d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1809d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1810d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1811d32f9abdSBarry Smith 1812c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1813c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1814c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 181513b05affSJed Brown $ ( I -ksp(A00) A01 ) ( inv(A00) 0 ) ( I 0 ) 1816c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 181713b05affSJed Brown where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. S is the Schur complement 181813b05affSJed Brown $ S = A11 - A10 ksp(A00) A01 181913b05affSJed Brown which is usually dense and not stored explicitly. The action of ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1820c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1821c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1822a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1823c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1824c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 18255668aaf4SBarry Smith diag gives 1826c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1827c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 18285668aaf4SBarry 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 1829c8a0d604SMatthew G Knepley $ ( A00 0 ) 1830c8a0d604SMatthew G Knepley $ ( A10 S ) 1831c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1832c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1833c8a0d604SMatthew G Knepley $ ( 0 S ) 1834c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1835e69d4d44SBarry Smith 1836edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1837edf189efSBarry Smith is used automatically for a second block. 1838edf189efSBarry Smith 1839ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1840ff218e97SBarry Smith Generally it should be used with the AIJ format. 1841ff218e97SBarry Smith 1842ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1843ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1844ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 18450716a85fSBarry Smith 1846a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 18470971522cSBarry Smith 18487e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1849e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 18500971522cSBarry Smith M*/ 18510971522cSBarry Smith 18520971522cSBarry Smith #undef __FUNCT__ 18530971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 18548cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc) 18550971522cSBarry Smith { 18560971522cSBarry Smith PetscErrorCode ierr; 18570971522cSBarry Smith PC_FieldSplit *jac; 18580971522cSBarry Smith 18590971522cSBarry Smith PetscFunctionBegin; 1860b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 18612fa5cd67SKarl Rupp 18620971522cSBarry Smith jac->bs = -1; 18630971522cSBarry Smith jac->nsplits = 0; 18643e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1865e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1866c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 1867fbe7908bSJed Brown jac->dm_splits = PETSC_TRUE; 186851f519a2SBarry Smith 18690971522cSBarry Smith pc->data = (void*)jac; 18700971522cSBarry Smith 18710971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1872421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 18730971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1874574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 18750971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 18760971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 18770971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 18780971522cSBarry Smith pc->ops->applyrichardson = 0; 18790971522cSBarry Smith 1880bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 1881bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1882bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 1883bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 1884bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 18850971522cSBarry Smith PetscFunctionReturn(0); 18860971522cSBarry Smith } 18870971522cSBarry Smith 18880971522cSBarry Smith 1889a541d17aSBarry Smith 1890