1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 307475bc1SBarry Smith #include <petscdmcomposite.h> 40971522cSBarry Smith 5443836d0SMatthew G Knepley /* 6443836d0SMatthew G Knepley There is a nice discussion of block preconditioners in 7443836d0SMatthew G Knepley 8443836d0SMatthew G Knepley [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier–Stokes equations 9443836d0SMatthew G Knepley Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808 101b8e4c5fSJed Brown http://chess.cs.umd.edu/~elman/papers/tax.pdf 11443836d0SMatthew G Knepley */ 12443836d0SMatthew G Knepley 13e87fab1bSBarry Smith const char *const PCFieldSplitSchurPreTypes[] = {"SELF","A11","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 14c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0}; 15c5d2311dSJed Brown 160971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 170971522cSBarry Smith struct _PC_FieldSplitLink { 1869a612a9SBarry Smith KSP ksp; 19443836d0SMatthew G Knepley Vec x,y,z; 20db4c96c1SJed Brown char *splitname; 210971522cSBarry Smith PetscInt nfields; 225d4c12cdSJungho Lee PetscInt *fields,*fields_col; 231b9fc7fcSBarry Smith VecScatter sctx; 245d4c12cdSJungho Lee IS is,is_col; 2551f519a2SBarry Smith PC_FieldSplitLink next,previous; 260971522cSBarry Smith }; 270971522cSBarry Smith 280971522cSBarry Smith typedef struct { 2968dd23aaSBarry Smith PCCompositeType type; 30ace3abfcSBarry Smith PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 31ace3abfcSBarry Smith PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 32ace3abfcSBarry Smith PetscBool realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 3330ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 3430ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 3579416396SBarry Smith Vec *x,*y,w1,w2; 36519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 37519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 3830ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 39ace3abfcSBarry Smith PetscBool issetup; 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 */ 520971522cSBarry Smith } PC_FieldSplit; 530971522cSBarry Smith 5416913363SBarry Smith /* 5516913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 5616913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 5716913363SBarry Smith PC you could change this. 5816913363SBarry Smith */ 59084e4875SJed Brown 60e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 61084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 62084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 63084e4875SJed Brown { 64084e4875SJed Brown switch (jac->schurpre) { 65084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 66e87fab1bSBarry Smith case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1]; 67084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 68084e4875SJed Brown default: 69084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 70084e4875SJed Brown } 71084e4875SJed Brown } 72084e4875SJed Brown 73084e4875SJed Brown 740971522cSBarry Smith #undef __FUNCT__ 750971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 760971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 770971522cSBarry Smith { 780971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 790971522cSBarry Smith PetscErrorCode ierr; 80d9884438SBarry Smith PetscBool iascii,isdraw; 810971522cSBarry Smith PetscInt i,j; 825a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 830971522cSBarry Smith 840971522cSBarry Smith PetscFunctionBegin; 85251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 86d9884438SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 870971522cSBarry Smith if (iascii) { 882c9966d7SBarry Smith if (jac->bs > 0) { 8951f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 902c9966d7SBarry Smith } else { 912c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr); 922c9966d7SBarry Smith } 93a3df900dSMatthew G Knepley if (jac->realdiagonal) { 94a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 95a3df900dSMatthew G Knepley } 9669a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 970971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 980971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 991ab39975SBarry Smith if (ilink->fields) { 1000971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 10179416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1025a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 10379416396SBarry Smith if (j > 0) { 10479416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 10579416396SBarry Smith } 1065a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1070971522cSBarry Smith } 1080971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 10979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1101ab39975SBarry Smith } else { 1111ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1121ab39975SBarry Smith } 1135a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 1145a9f2f41SSatish Balay ilink = ilink->next; 1150971522cSBarry Smith } 1160971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 117d9884438SBarry Smith } if (isdraw) { 118d9884438SBarry Smith PetscDraw draw; 119d9884438SBarry Smith PetscReal x,y,w,wd; 120d9884438SBarry Smith 121d9884438SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 122d9884438SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 123d9884438SBarry Smith w = 2*PetscMin(1.0 - x,x); 124d9884438SBarry Smith wd = w/(jac->nsplits + 1); 125d9884438SBarry Smith x = x - wd*(jac->nsplits-1)/2.0; 126d9884438SBarry Smith for (i=0; i<jac->nsplits; i++) { 127d9884438SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 128d9884438SBarry Smith ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 129d9884438SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 130d9884438SBarry Smith x += wd; 131d9884438SBarry Smith ilink = ilink->next; 132d9884438SBarry Smith } 1330971522cSBarry Smith } 1340971522cSBarry Smith PetscFunctionReturn(0); 1350971522cSBarry Smith } 1360971522cSBarry Smith 1370971522cSBarry Smith #undef __FUNCT__ 1383b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1393b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1403b224e63SBarry Smith { 1413b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1423b224e63SBarry Smith PetscErrorCode ierr; 1434996c5bdSBarry Smith PetscBool iascii,isdraw; 1443b224e63SBarry Smith PetscInt i,j; 1453b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1463b224e63SBarry Smith 1473b224e63SBarry Smith PetscFunctionBegin; 148251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1494996c5bdSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 1503b224e63SBarry Smith if (iascii) { 1512c9966d7SBarry Smith if (jac->bs > 0) { 152c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1532c9966d7SBarry Smith } else { 1542c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1552c9966d7SBarry Smith } 156a3df900dSMatthew G Knepley if (jac->realdiagonal) { 157a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 158a3df900dSMatthew G Knepley } 1593e8b8b31SMatthew G Knepley switch (jac->schurpre) { 1603e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1613e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 162e87fab1bSBarry Smith case PC_FIELDSPLIT_SCHUR_PRE_A11: 163e87fab1bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break; 1643e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1653e8b8b31SMatthew G Knepley if (jac->schur_user) { 1663e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1673e8b8b31SMatthew G Knepley } else { 168e87fab1bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr); 1693e8b8b31SMatthew G Knepley } 1703e8b8b31SMatthew G Knepley break; 1713e8b8b31SMatthew G Knepley default: 1723e8b8b31SMatthew G Knepley SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1733e8b8b31SMatthew G Knepley } 1743b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1753b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1763b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1773b224e63SBarry Smith if (ilink->fields) { 1783b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1793b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1803b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1813b224e63SBarry Smith if (j > 0) { 1823b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1833b224e63SBarry Smith } 1843b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1853b224e63SBarry Smith } 1863b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1873b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1883b224e63SBarry Smith } else { 1893b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1903b224e63SBarry Smith } 1913b224e63SBarry Smith ilink = ilink->next; 1923b224e63SBarry Smith } 193435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1943b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 195443836d0SMatthew G Knepley ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 1963b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 197443836d0SMatthew G Knepley if (jac->kspupper != jac->head->ksp) { 198443836d0SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr); 199443836d0SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 200443836d0SMatthew G Knepley ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 201443836d0SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 202443836d0SMatthew G Knepley } 203435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 2043b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 20512cae6f2SJed Brown if (jac->kspschur) { 2063b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 20712cae6f2SJed Brown } else { 20812cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 20912cae6f2SJed Brown } 2103b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2113b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2124996c5bdSBarry Smith } else if (isdraw) { 2134996c5bdSBarry Smith PetscDraw draw; 2144996c5bdSBarry Smith PetscReal x,y,w,wd,h; 2154996c5bdSBarry Smith PetscInt cnt = 2; 2164996c5bdSBarry Smith char str[32]; 2174996c5bdSBarry Smith 2184996c5bdSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 2194996c5bdSBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 220c74581afSBarry Smith if (jac->kspupper != jac->head->ksp) cnt++; 221c74581afSBarry Smith w = 2*PetscMin(1.0 - x,x); 222c74581afSBarry Smith wd = w/(cnt + 1); 223c74581afSBarry Smith 2244996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 2254996c5bdSBarry Smith ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,PETSC_NULL,&h);CHKERRQ(ierr); 2264996c5bdSBarry Smith y -= h; 2274996c5bdSBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) { 228e87fab1bSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr); 2293b224e63SBarry Smith } else { 2304996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr); 2314996c5bdSBarry Smith } 232c74581afSBarry Smith ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,PETSC_NULL,&h);CHKERRQ(ierr); 2334996c5bdSBarry Smith y -= h; 2344996c5bdSBarry Smith x = x - wd*(cnt-1)/2.0; 2354996c5bdSBarry Smith 2364996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2374996c5bdSBarry Smith ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 2384996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2394996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) { 2404996c5bdSBarry Smith x += wd; 2414996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2424996c5bdSBarry Smith ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 2434996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2444996c5bdSBarry Smith } 2454996c5bdSBarry Smith x += wd; 2464996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 2474996c5bdSBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 2484996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2493b224e63SBarry Smith } 2503b224e63SBarry Smith PetscFunctionReturn(0); 2513b224e63SBarry Smith } 2523b224e63SBarry Smith 2533b224e63SBarry Smith #undef __FUNCT__ 2546c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2556c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2566c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2576c924f48SJed Brown { 2586c924f48SJed Brown PetscErrorCode ierr; 2596c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2605d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2615d4c12cdSJungho Lee PetscBool flg,flg_col; 2625d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2636c924f48SJed Brown 2646c924f48SJed Brown PetscFunctionBegin; 2656c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2665d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2676c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2688caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2698caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2708caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2716c924f48SJed Brown nfields = jac->bs; 27229499fbbSJungho Lee nfields_col = jac->bs; 2736c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2745d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2756c924f48SJed Brown if (!flg) break; 2765d4c12cdSJungho Lee else if (flg && !flg_col) { 2775d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2785d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2795d4c12cdSJungho Lee } 2805d4c12cdSJungho Lee else { 2815d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2825d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2835d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2845d4c12cdSJungho Lee } 2856c924f48SJed Brown } 2866c924f48SJed Brown if (i > 0) { 2876c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2886c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2896c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2906c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2916c924f48SJed Brown } 2926c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2935d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2946c924f48SJed Brown PetscFunctionReturn(0); 2956c924f48SJed Brown } 2966c924f48SJed Brown 2976c924f48SJed Brown #undef __FUNCT__ 29869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 29969a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 3000971522cSBarry Smith { 3010971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3020971522cSBarry Smith PetscErrorCode ierr; 3035a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3047287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 3056c924f48SJed Brown PetscInt i; 3060971522cSBarry Smith 3070971522cSBarry Smith PetscFunctionBegin; 3087287d2a3SDmitry Karpeev /* 3097287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 3107287d2a3SDmitry Karpeev Should probably be rewritten. 3117287d2a3SDmitry Karpeev */ 3127287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 313d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 314d0af7cd3SBarry Smith if (pc->dm && !stokes) { 315bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 3160784a22cSJed Brown char **fieldNames; 3177b62db95SJungho Lee IS *fields; 318e7c4fc90SDmitry Karpeev DM *dms; 319bafc1b83SMatthew G Knepley DM subdm[128]; 320bafc1b83SMatthew G Knepley PetscBool flg; 321bafc1b83SMatthew G Knepley 32216621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 323bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 324bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE; ; i++) { 325bafc1b83SMatthew G Knepley PetscInt ifields[128]; 326bafc1b83SMatthew G Knepley IS compField; 327bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 328bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 329bafc1b83SMatthew G Knepley 3308caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 331bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 332bafc1b83SMatthew G Knepley if (!flg) break; 333bafc1b83SMatthew G Knepley if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 334bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 335bafc1b83SMatthew G Knepley if (nfields == 1) { 336bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 3374d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 3384d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 339bafc1b83SMatthew G Knepley } else { 3408caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 341bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 3424d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr); 3434d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 3447287d2a3SDmitry Karpeev } 345bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 346bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 347bafc1b83SMatthew G Knepley f = ifields[j]; 3487b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 3497b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 3507b62db95SJungho Lee } 351bafc1b83SMatthew G Knepley } 352bafc1b83SMatthew G Knepley if (i == 0) { 353bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 354bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 355bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 356bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 357bafc1b83SMatthew G Knepley } 358bafc1b83SMatthew G Knepley } else { 359d724dfffSBarry Smith for (j=0; j<numFields; j++) { 360d724dfffSBarry Smith ierr = DMDestroy(dms+j);CHKERRQ(ierr); 361d724dfffSBarry Smith } 362d724dfffSBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 363bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 364bafc1b83SMatthew G Knepley for (j = 0; j < i; ++j) { 365bafc1b83SMatthew G Knepley dms[j] = subdm[j]; 366bafc1b83SMatthew G Knepley } 367bafc1b83SMatthew G Knepley } 3687b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3697b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 370e7c4fc90SDmitry Karpeev if (dms) { 3718b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 372bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3737287d2a3SDmitry Karpeev const char *prefix; 3747287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr); 3757287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr); 3767b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3777b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3787287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr); 379e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 3802fa5ba8aSJed Brown } 3817b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3828b8307b2SJed Brown } 38366ffff09SJed Brown } else { 384521d7252SBarry Smith if (jac->bs <= 0) { 385704ba839SBarry Smith if (pc->pmat) { 386521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 387704ba839SBarry Smith } else { 388704ba839SBarry Smith jac->bs = 1; 389704ba839SBarry Smith } 390521d7252SBarry Smith } 391d32f9abdSBarry Smith 3927287d2a3SDmitry Karpeev ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr); 3936ce1633cSBarry Smith if (stokes) { 3946ce1633cSBarry Smith IS zerodiags,rest; 3956ce1633cSBarry Smith PetscInt nmin,nmax; 3966ce1633cSBarry Smith 3976ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3986ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3996ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 4007287d2a3SDmitry Karpeev if (jac->reset) { 4017287d2a3SDmitry Karpeev jac->head->is = rest; 4027287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 4037287d2a3SDmitry Karpeev } 4047287d2a3SDmitry Karpeev else { 4056ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 4066ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 4077287d2a3SDmitry Karpeev } 408fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 409fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 4106ce1633cSBarry Smith } else { 411*f23aa3ddSBarry Smith if (jac->reset) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 4127287d2a3SDmitry Karpeev if (!fieldsplit_default) { 413d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 414d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 4156c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 4166c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 417d32f9abdSBarry Smith } 4187287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 419d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 420db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 4216c924f48SJed Brown char splitname[8]; 4228caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 4235d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 42479416396SBarry Smith } 4255d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 426521d7252SBarry Smith } 42766ffff09SJed Brown } 4286ce1633cSBarry Smith } 429edf189efSBarry Smith } else if (jac->nsplits == 1) { 430edf189efSBarry Smith if (ilink->is) { 431edf189efSBarry Smith IS is2; 432edf189efSBarry Smith PetscInt nmin,nmax; 433edf189efSBarry Smith 434edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 435edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 436db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 437fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 4387b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 439edf189efSBarry Smith } 440d0af7cd3SBarry Smith 441d0af7cd3SBarry Smith 4427b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 44369a612a9SBarry Smith PetscFunctionReturn(0); 44469a612a9SBarry Smith } 44569a612a9SBarry Smith 446514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 447514bf10dSMatthew G Knepley 44869a612a9SBarry Smith #undef __FUNCT__ 44969a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 45069a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 45169a612a9SBarry Smith { 45269a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 45369a612a9SBarry Smith PetscErrorCode ierr; 4545a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4552c9966d7SBarry Smith PetscInt i,nsplit; 45669a612a9SBarry Smith MatStructure flag = pc->flag; 4575f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 45869a612a9SBarry Smith 45969a612a9SBarry Smith PetscFunctionBegin; 46069a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 46197bbdb24SBarry Smith nsplit = jac->nsplits; 4625a9f2f41SSatish Balay ilink = jac->head; 46397bbdb24SBarry Smith 46497bbdb24SBarry Smith /* get the matrices for each split */ 465704ba839SBarry Smith if (!jac->issetup) { 4661b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 46797bbdb24SBarry Smith 468704ba839SBarry Smith jac->issetup = PETSC_TRUE; 469704ba839SBarry Smith 4705d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4712c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4722c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4732c9966d7SBarry Smith } 47451f519a2SBarry Smith bs = jac->bs; 47597bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4761b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4771b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4781b9fc7fcSBarry Smith if (jac->defaultsplit) { 479704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4805f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 481704ba839SBarry Smith } else if (!ilink->is) { 482ccb205f8SBarry Smith if (ilink->nfields > 1) { 4835f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4845a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4855f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4861b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4871b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4881b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4895f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 49097bbdb24SBarry Smith } 49197bbdb24SBarry Smith } 492d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4935f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 494ccb205f8SBarry Smith } else { 495704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4965f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 497ccb205f8SBarry Smith } 4983e197d65SBarry Smith } 499704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 5005f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 5015f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 502704ba839SBarry Smith ilink = ilink->next; 5031b9fc7fcSBarry Smith } 5041b9fc7fcSBarry Smith } 5051b9fc7fcSBarry Smith 506704ba839SBarry Smith ilink = jac->head; 50797bbdb24SBarry Smith if (!jac->pmat) { 508bdddcaaaSMatthew G Knepley Vec xtmp; 509bdddcaaaSMatthew G Knepley 510bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 511cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 512bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 513cf502942SBarry Smith for (i=0; i<nsplit; i++) { 514bdddcaaaSMatthew G Knepley MatNullSpace sp; 515bdddcaaaSMatthew G Knepley 516a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 517a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 518a3df900dSMatthew G Knepley if (jac->pmat[i]) { 519a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 520a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 521a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 522a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 523a3df900dSMatthew G Knepley } 524a3df900dSMatthew G Knepley } else { 5255f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 526a3df900dSMatthew G Knepley } 527bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 528bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 529443836d0SMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = PETSC_NULL; 530bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 531bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 532ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 533bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr); 534bafc1b83SMatthew G Knepley if (sp) { 535bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 536bafc1b83SMatthew G Knepley } 537ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 538ed1f0337SMatthew G Knepley if (sp) { 539ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 540ed1f0337SMatthew G Knepley } 541704ba839SBarry Smith ilink = ilink->next; 542cf502942SBarry Smith } 543bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 54497bbdb24SBarry Smith } else { 545cf502942SBarry Smith for (i=0; i<nsplit; i++) { 546a3df900dSMatthew G Knepley Mat pmat; 547a3df900dSMatthew G Knepley 548a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 549a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 550a3df900dSMatthew G Knepley if (!pmat) { 5515f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 552a3df900dSMatthew G Knepley } 553704ba839SBarry Smith ilink = ilink->next; 554cf502942SBarry Smith } 55597bbdb24SBarry Smith } 556519d70e2SJed Brown if (jac->realdiagonal) { 557519d70e2SJed Brown ilink = jac->head; 558519d70e2SJed Brown if (!jac->mat) { 559519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 560519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5615f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 562519d70e2SJed Brown ilink = ilink->next; 563519d70e2SJed Brown } 564519d70e2SJed Brown } else { 565519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5665f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 567519d70e2SJed Brown ilink = ilink->next; 568519d70e2SJed Brown } 569519d70e2SJed Brown } 570519d70e2SJed Brown } else { 571519d70e2SJed Brown jac->mat = jac->pmat; 572519d70e2SJed Brown } 57397bbdb24SBarry Smith 5746c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 57568dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 57668dd23aaSBarry Smith ilink = jac->head; 57768dd23aaSBarry Smith if (!jac->Afield) { 57868dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 57968dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5804aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58168dd23aaSBarry Smith ilink = ilink->next; 58268dd23aaSBarry Smith } 58368dd23aaSBarry Smith } else { 58468dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5854aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58668dd23aaSBarry Smith ilink = ilink->next; 58768dd23aaSBarry Smith } 58868dd23aaSBarry Smith } 58968dd23aaSBarry Smith } 59068dd23aaSBarry Smith 5913b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5923b224e63SBarry Smith IS ccis; 5934aa3045dSJed Brown PetscInt rstart,rend; 594093c86ffSJed Brown char lscname[256]; 595093c86ffSJed Brown PetscObject LSC_L; 596e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 59768dd23aaSBarry Smith 598e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 599e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 600e6cab6aaSJed Brown 6013b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 6023b224e63SBarry Smith if (jac->schur) { 603a3314f2cSMatthew G Knepley KSP kspA = jac->head->ksp, kspInner = PETSC_NULL, kspUpper = jac->kspupper; 604443836d0SMatthew G Knepley 605fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 6063b224e63SBarry Smith ilink = jac->head; 60749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6084aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 609fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6103b224e63SBarry Smith ilink = ilink->next; 61149bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6124aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 613fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 614a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 615443836d0SMatthew G Knepley if (kspA != kspInner) { 616443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 617443836d0SMatthew G Knepley } 618443836d0SMatthew G Knepley if (kspUpper != kspA) { 619443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 620443836d0SMatthew G Knepley } 621084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 6223b224e63SBarry Smith } else { 6231cee3971SBarry Smith KSP ksp; 624bafc1b83SMatthew G Knepley const char *Dprefix; 625bafc1b83SMatthew G Knepley char schurprefix[256]; 626514bf10dSMatthew G Knepley char schurtestoption[256]; 627bdddcaaaSMatthew G Knepley MatNullSpace sp; 628514bf10dSMatthew G Knepley PetscBool flg; 6293b224e63SBarry Smith 630a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 6313b224e63SBarry Smith ilink = jac->head; 63249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6334aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 634fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6353b224e63SBarry Smith ilink = ilink->next; 63649bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6374aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 638fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 63920252d06SBarry Smith 64020252d06SBarry Smith /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */ 64120252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 64220252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 64320252d06SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr); 64420252d06SBarry Smith ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 64520252d06SBarry Smith /* Indent this deeper to emphasize the "inner" nature of this solver. */ 64620252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr); 64720252d06SBarry Smith ierr = KSPSetOptionsPrefix(ksp, schurprefix);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); 656514bf10dSMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 657514bf10dSMatthew G Knepley if (flg) { 658514bf10dSMatthew G Knepley DM dmInner; 659514bf10dSMatthew G Knepley 660514bf10dSMatthew G Knepley /* Set DM for new solver */ 661bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 662bafc1b83SMatthew G Knepley ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr); 6637287d2a3SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 664443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 665443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 666514bf10dSMatthew G Knepley } else { 667514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr); 668514bf10dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 669bafc1b83SMatthew G Knepley } 67020b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 6713b224e63SBarry Smith 672443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 673443836d0SMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 674443836d0SMatthew G Knepley if (flg) { 675443836d0SMatthew G Knepley DM dmInner; 676443836d0SMatthew G Knepley 677443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 678443836d0SMatthew G Knepley ierr = KSPCreate(((PetscObject) pc)->comm, &jac->kspupper);CHKERRQ(ierr); 679443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 680443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 681443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 682443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 683443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 684443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 685443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 686443836d0SMatthew G Knepley } else { 687443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 688443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 689443836d0SMatthew G Knepley } 690443836d0SMatthew G Knepley 6913b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 6929005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 69320252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 694084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 695084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 6967233a360SDmitry Karpeev PC pcschur; 6977233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 6987233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 699084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 700e69d4d44SBarry Smith } 7017287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr); 7027287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix);CHKERRQ(ierr); 7033b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 70420b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 70520b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 7063b224e63SBarry Smith } 707093c86ffSJed Brown 708093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 7098caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 710093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 711093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 712093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 7138caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 714093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 715093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 716093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 7173b224e63SBarry Smith } else { 71868bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 71997bbdb24SBarry Smith i = 0; 7205a9f2f41SSatish Balay ilink = jac->head; 7215a9f2f41SSatish Balay while (ilink) { 722519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 7233b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 724c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 72597bbdb24SBarry Smith i++; 7265a9f2f41SSatish Balay ilink = ilink->next; 7270971522cSBarry Smith } 7283b224e63SBarry Smith } 7293b224e63SBarry Smith 730c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 7310971522cSBarry Smith PetscFunctionReturn(0); 7320971522cSBarry Smith } 7330971522cSBarry Smith 7345a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 735ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 736ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 7375a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 738ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 739ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 74079416396SBarry Smith 7410971522cSBarry Smith #undef __FUNCT__ 7423b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7433b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7443b224e63SBarry Smith { 7453b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7463b224e63SBarry Smith PetscErrorCode ierr; 7473b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 748443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7493b224e63SBarry Smith 7503b224e63SBarry Smith PetscFunctionBegin; 751c5d2311dSJed Brown switch (jac->schurfactorization) { 752c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 753a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 754c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 755c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 756c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 757443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 758c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 759c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 760c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 761c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 762c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 763c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 764c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 765c5d2311dSJed Brown break; 766c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 767a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 768c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 769c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 770443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 771c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 772c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 773c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 774c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 775c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 776c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 777c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 778c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 779c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 780c5d2311dSJed Brown break; 781c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 782a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 783c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 784c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 785c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 786c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 787c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 788c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 789c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 790c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 791443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 792c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 793c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 794c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 795c5d2311dSJed Brown break; 796c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 797a04f6461SBarry 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 */ 7983b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7993b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 800443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 8013b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 8023b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 8033b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8043b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8053b224e63SBarry Smith 8063b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 8073b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8083b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8093b224e63SBarry Smith 810443836d0SMatthew G Knepley if (kspUpper == kspA) { 8113b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 8123b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 813443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 814443836d0SMatthew G Knepley } else { 815443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 816443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 817443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 818443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 819443836d0SMatthew G Knepley } 8203b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8213b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 822c5d2311dSJed Brown } 8233b224e63SBarry Smith PetscFunctionReturn(0); 8243b224e63SBarry Smith } 8253b224e63SBarry Smith 8263b224e63SBarry Smith #undef __FUNCT__ 8270971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 8280971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 8290971522cSBarry Smith { 8300971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8310971522cSBarry Smith PetscErrorCode ierr; 8325a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 833939b8a20SBarry Smith PetscInt cnt,bs; 8340971522cSBarry Smith 8350971522cSBarry Smith PetscFunctionBegin; 8364442daceSBarry Smith x->map->bs = jac->bs; 8374442daceSBarry Smith y->map->bs = jac->bs; 83851f519a2SBarry Smith CHKMEMQ; 83979416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8401b9fc7fcSBarry Smith if (jac->defaultsplit) { 841939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 842939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 843939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 844939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 8450971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8465a9f2f41SSatish Balay while (ilink) { 8475a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8485a9f2f41SSatish Balay ilink = ilink->next; 8490971522cSBarry Smith } 8500971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8511b9fc7fcSBarry Smith } else { 852efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8535a9f2f41SSatish Balay while (ilink) { 8545a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8555a9f2f41SSatish Balay ilink = ilink->next; 8561b9fc7fcSBarry Smith } 8571b9fc7fcSBarry Smith } 85816913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 85979416396SBarry Smith if (!jac->w1) { 86079416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 86179416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 86279416396SBarry Smith } 863efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8645a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8653e197d65SBarry Smith cnt = 1; 8665a9f2f41SSatish Balay while (ilink->next) { 8675a9f2f41SSatish Balay ilink = ilink->next; 8683e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8693e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8703e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8713e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8723e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8733e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8743e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8753e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8763e197d65SBarry Smith } 87751f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 87811755939SBarry Smith cnt -= 2; 87951f519a2SBarry Smith while (ilink->previous) { 88051f519a2SBarry Smith ilink = ilink->previous; 88111755939SBarry Smith /* compute the residual only over the part of the vector needed */ 88211755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 88311755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 88411755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88511755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88611755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 88711755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88811755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88951f519a2SBarry Smith } 89011755939SBarry Smith } 89165e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 89251f519a2SBarry Smith CHKMEMQ; 8930971522cSBarry Smith PetscFunctionReturn(0); 8940971522cSBarry Smith } 8950971522cSBarry Smith 896421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 897ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 898ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 899421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 900ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 901ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 902421e10b8SBarry Smith 903421e10b8SBarry Smith #undef __FUNCT__ 9048c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 905421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 906421e10b8SBarry Smith { 907421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 908421e10b8SBarry Smith PetscErrorCode ierr; 909421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 910939b8a20SBarry Smith PetscInt bs; 911421e10b8SBarry Smith 912421e10b8SBarry Smith PetscFunctionBegin; 913421e10b8SBarry Smith CHKMEMQ; 914421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 915421e10b8SBarry Smith if (jac->defaultsplit) { 916939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 917939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 918939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 919939b8a20SBarry Smith if (jac->bs > 0 && bs != jac->bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs); 920421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 921421e10b8SBarry Smith while (ilink) { 922421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 923421e10b8SBarry Smith ilink = ilink->next; 924421e10b8SBarry Smith } 925421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 926421e10b8SBarry Smith } else { 927421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 928421e10b8SBarry Smith while (ilink) { 929421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 930421e10b8SBarry Smith ilink = ilink->next; 931421e10b8SBarry Smith } 932421e10b8SBarry Smith } 933421e10b8SBarry Smith } else { 934421e10b8SBarry Smith if (!jac->w1) { 935421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 936421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 937421e10b8SBarry Smith } 938421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 939421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 940421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 941421e10b8SBarry Smith while (ilink->next) { 942421e10b8SBarry Smith ilink = ilink->next; 9439989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 944421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 945421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 946421e10b8SBarry Smith } 947421e10b8SBarry Smith while (ilink->previous) { 948421e10b8SBarry Smith ilink = ilink->previous; 9499989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 950421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 951421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 952421e10b8SBarry Smith } 953421e10b8SBarry Smith } else { 954421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 955421e10b8SBarry Smith ilink = ilink->next; 956421e10b8SBarry Smith } 957421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 958421e10b8SBarry Smith while (ilink->previous) { 959421e10b8SBarry Smith ilink = ilink->previous; 9609989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 961421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 962421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 963421e10b8SBarry Smith } 964421e10b8SBarry Smith } 965421e10b8SBarry Smith } 966421e10b8SBarry Smith CHKMEMQ; 967421e10b8SBarry Smith PetscFunctionReturn(0); 968421e10b8SBarry Smith } 969421e10b8SBarry Smith 9700971522cSBarry Smith #undef __FUNCT__ 971574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 972574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9730971522cSBarry Smith { 9740971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9750971522cSBarry Smith PetscErrorCode ierr; 9765a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9770971522cSBarry Smith 9780971522cSBarry Smith PetscFunctionBegin; 9795a9f2f41SSatish Balay while (ilink) { 980574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 981fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 982fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 983443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 984fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 985fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 986b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9875a9f2f41SSatish Balay next = ilink->next; 9885a9f2f41SSatish Balay ilink = next; 9890971522cSBarry Smith } 99005b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 991574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 992574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 993574deadeSBarry Smith } else if (jac->mat) { 994574deadeSBarry Smith jac->mat = PETSC_NULL; 995574deadeSBarry Smith } 99697bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 99768dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 9986bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 9996bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 10006bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 10016bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 10026bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 1003d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 10046bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 10056bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 100663ec74ffSBarry Smith jac->reset = PETSC_TRUE; 1007574deadeSBarry Smith PetscFunctionReturn(0); 1008574deadeSBarry Smith } 1009574deadeSBarry Smith 1010574deadeSBarry Smith #undef __FUNCT__ 1011574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1012574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1013574deadeSBarry Smith { 1014574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1015574deadeSBarry Smith PetscErrorCode ierr; 1016574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 1017574deadeSBarry Smith 1018574deadeSBarry Smith PetscFunctionBegin; 1019574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 1020574deadeSBarry Smith while (ilink) { 10216bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 1022574deadeSBarry Smith next = ilink->next; 1023574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 1024574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 10255d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 1026574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 1027574deadeSBarry Smith ilink = next; 1028574deadeSBarry Smith } 1029574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1030c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 10317b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 10327b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 10337b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 10347b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 10357b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 1036ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 1037c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 10380971522cSBarry Smith PetscFunctionReturn(0); 10390971522cSBarry Smith } 10400971522cSBarry Smith 10410971522cSBarry Smith #undef __FUNCT__ 10420971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10430971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10440971522cSBarry Smith { 10451b9fc7fcSBarry Smith PetscErrorCode ierr; 10466c924f48SJed Brown PetscInt bs; 1047bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10489dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10493b224e63SBarry Smith PCCompositeType ctype; 1050e7c4fc90SDmitry Karpeev DM ddm; 1051e7c4fc90SDmitry Karpeev char ddm_name[1025]; 10521b9fc7fcSBarry Smith 10530971522cSBarry Smith PetscFunctionBegin; 10541b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 1055e7c4fc90SDmitry Karpeev if (pc->dm) { 1056e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 1057e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024);CHKERRQ(ierr); 1058731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg);CHKERRQ(ierr); 1059e7c4fc90SDmitry Karpeev if (flg) { 106016621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm);CHKERRQ(ierr); 1061*f23aa3ddSBarry Smith if (!ddm) SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 106216621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 1063e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm);CHKERRQ(ierr); 1064e7c4fc90SDmitry Karpeev } 1065e7c4fc90SDmitry Karpeev } 1066acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 106751f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 106851f519a2SBarry Smith if (flg) { 106951f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 107051f519a2SBarry Smith } 1071704ba839SBarry Smith 1072435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 1073c0adfefeSBarry Smith if (stokes) { 1074c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1075c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1076c0adfefeSBarry Smith } 1077c0adfefeSBarry Smith 10783b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10793b224e63SBarry Smith if (flg) { 10803b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10813b224e63SBarry Smith } 1082c30613efSMatthew Knepley /* Only setup fields once */ 1083c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1084d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1085d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10866c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10876c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1088d32f9abdSBarry Smith } 1089c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1090c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1091c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 1092c9c6ffaaSJed Brown 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,PETSC_NULL);CHKERRQ(ierr); 1093084e4875SJed Brown ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,PETSC_NULL);CHKERRQ(ierr); 1094c5d2311dSJed Brown } 10951b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10960971522cSBarry Smith PetscFunctionReturn(0); 10970971522cSBarry Smith } 10980971522cSBarry Smith 10990971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 11000971522cSBarry Smith 11010971522cSBarry Smith EXTERN_C_BEGIN 11020971522cSBarry Smith #undef __FUNCT__ 11030971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 11045d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11050971522cSBarry Smith { 110697bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11070971522cSBarry Smith PetscErrorCode ierr; 11085a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 110969a612a9SBarry Smith char prefix[128]; 11105d4c12cdSJungho Lee PetscInt i; 11110971522cSBarry Smith 11120971522cSBarry Smith PetscFunctionBegin; 11136c924f48SJed Brown if (jac->splitdefined) { 11146c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11156c924f48SJed Brown PetscFunctionReturn(0); 11166c924f48SJed Brown } 111751f519a2SBarry Smith for (i=0; i<n; i++) { 1118e32f2f54SBarry 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); 1119e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 112051f519a2SBarry Smith } 1121704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1122a04f6461SBarry Smith if (splitname) { 1123db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1124a04f6461SBarry Smith } else { 1125a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1126a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1127a04f6461SBarry Smith } 1128704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 11295a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 11305d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 11315d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 11325a9f2f41SSatish Balay ilink->nfields = n; 11335a9f2f41SSatish Balay ilink->next = PETSC_NULL; 11347adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 113520252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11365a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11379005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 113869a612a9SBarry Smith 11398caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 11405a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11410971522cSBarry Smith 11420971522cSBarry Smith if (!next) { 11435a9f2f41SSatish Balay jac->head = ilink; 114451f519a2SBarry Smith ilink->previous = PETSC_NULL; 11450971522cSBarry Smith } else { 11460971522cSBarry Smith while (next->next) { 11470971522cSBarry Smith next = next->next; 11480971522cSBarry Smith } 11495a9f2f41SSatish Balay next->next = ilink; 115051f519a2SBarry Smith ilink->previous = next; 11510971522cSBarry Smith } 11520971522cSBarry Smith jac->nsplits++; 11530971522cSBarry Smith PetscFunctionReturn(0); 11540971522cSBarry Smith } 11550971522cSBarry Smith EXTERN_C_END 11560971522cSBarry Smith 1157e69d4d44SBarry Smith EXTERN_C_BEGIN 1158e69d4d44SBarry Smith #undef __FUNCT__ 1159e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11607087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1161e69d4d44SBarry Smith { 1162e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1163e69d4d44SBarry Smith PetscErrorCode ierr; 1164e69d4d44SBarry Smith 1165e69d4d44SBarry Smith PetscFunctionBegin; 1166519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1167e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1168e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 116913e0d083SBarry Smith if (n) *n = jac->nsplits; 1170e69d4d44SBarry Smith PetscFunctionReturn(0); 1171e69d4d44SBarry Smith } 1172e69d4d44SBarry Smith EXTERN_C_END 11730971522cSBarry Smith 11740971522cSBarry Smith EXTERN_C_BEGIN 11750971522cSBarry Smith #undef __FUNCT__ 117669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11777087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11780971522cSBarry Smith { 11790971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11800971522cSBarry Smith PetscErrorCode ierr; 11810971522cSBarry Smith PetscInt cnt = 0; 11825a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11830971522cSBarry Smith 11840971522cSBarry Smith PetscFunctionBegin; 11855d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11865a9f2f41SSatish Balay while (ilink) { 11875a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11885a9f2f41SSatish Balay ilink = ilink->next; 11890971522cSBarry Smith } 11905d480477SMatthew 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); 119113e0d083SBarry Smith if (n) *n = jac->nsplits; 11920971522cSBarry Smith PetscFunctionReturn(0); 11930971522cSBarry Smith } 11940971522cSBarry Smith EXTERN_C_END 11950971522cSBarry Smith 1196704ba839SBarry Smith EXTERN_C_BEGIN 1197704ba839SBarry Smith #undef __FUNCT__ 1198704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11997087cfbeSBarry Smith 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 } 121116913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1212a04f6461SBarry Smith if (splitname) { 1213db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1214a04f6461SBarry Smith } else { 1215b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&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; 1224704ba839SBarry Smith ilink->next = PETSC_NULL; 1225704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&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; 1235704ba839SBarry Smith ilink->previous = PETSC_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 1245704ba839SBarry Smith PetscFunctionReturn(0); 1246704ba839SBarry Smith } 1247704ba839SBarry Smith EXTERN_C_END 1248704ba839SBarry Smith 12490971522cSBarry Smith #undef __FUNCT__ 12500971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12510971522cSBarry Smith /*@ 12520971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12530971522cSBarry Smith 1254ad4df100SBarry Smith Logically Collective on PC 12550971522cSBarry Smith 12560971522cSBarry Smith Input Parameters: 12570971522cSBarry Smith + pc - the preconditioner context 1258a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 12590971522cSBarry Smith . n - the number of fields in this split 1260db4c96c1SJed Brown - fields - the fields in this split 12610971522cSBarry Smith 12620971522cSBarry Smith Level: intermediate 12630971522cSBarry Smith 1264d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1265d32f9abdSBarry Smith 12667287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1267d32f9abdSBarry 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 1268d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1269d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1270d32f9abdSBarry Smith 1271db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1272db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1273db4c96c1SJed Brown 12745d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12755d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12765d4c12cdSJungho Lee available when this routine is called. 12775d4c12cdSJungho Lee 1278d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12790971522cSBarry Smith 12800971522cSBarry Smith @*/ 12815d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12820971522cSBarry Smith { 12834ac538c5SBarry Smith PetscErrorCode ierr; 12840971522cSBarry Smith 12850971522cSBarry Smith PetscFunctionBegin; 12860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1287db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1288db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1289db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12905d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12910971522cSBarry Smith PetscFunctionReturn(0); 12920971522cSBarry Smith } 12930971522cSBarry Smith 12940971522cSBarry Smith #undef __FUNCT__ 1295704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1296704ba839SBarry Smith /*@ 1297704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1298704ba839SBarry Smith 1299ad4df100SBarry Smith Logically Collective on PC 1300704ba839SBarry Smith 1301704ba839SBarry Smith Input Parameters: 1302704ba839SBarry Smith + pc - the preconditioner context 1303a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1304db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1305704ba839SBarry Smith 1306d32f9abdSBarry Smith 1307a6ffb8dbSJed Brown Notes: 1308a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1309a6ffb8dbSJed Brown 1310db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1311db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1312d32f9abdSBarry Smith 1313704ba839SBarry Smith Level: intermediate 1314704ba839SBarry Smith 1315704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1316704ba839SBarry Smith 1317704ba839SBarry Smith @*/ 13187087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1319704ba839SBarry Smith { 13204ac538c5SBarry Smith PetscErrorCode ierr; 1321704ba839SBarry Smith 1322704ba839SBarry Smith PetscFunctionBegin; 13230700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13247b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1325db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 13264ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1327704ba839SBarry Smith PetscFunctionReturn(0); 1328704ba839SBarry Smith } 1329704ba839SBarry Smith 1330704ba839SBarry Smith #undef __FUNCT__ 133157a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 133257a9adfeSMatthew G Knepley /*@ 133357a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 133457a9adfeSMatthew G Knepley 133557a9adfeSMatthew G Knepley Logically Collective on PC 133657a9adfeSMatthew G Knepley 133757a9adfeSMatthew G Knepley Input Parameters: 133857a9adfeSMatthew G Knepley + pc - the preconditioner context 133957a9adfeSMatthew G Knepley - splitname - name of this split 134057a9adfeSMatthew G Knepley 134157a9adfeSMatthew G Knepley Output Parameter: 134257a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 134357a9adfeSMatthew G Knepley 134457a9adfeSMatthew G Knepley Level: intermediate 134557a9adfeSMatthew G Knepley 134657a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 134757a9adfeSMatthew G Knepley 134857a9adfeSMatthew G Knepley @*/ 134957a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 135057a9adfeSMatthew G Knepley { 135157a9adfeSMatthew G Knepley PetscErrorCode ierr; 135257a9adfeSMatthew G Knepley 135357a9adfeSMatthew G Knepley PetscFunctionBegin; 135457a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 135557a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 135657a9adfeSMatthew G Knepley PetscValidPointer(is,3); 135757a9adfeSMatthew G Knepley { 135857a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 135957a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 136057a9adfeSMatthew G Knepley PetscBool found; 136157a9adfeSMatthew G Knepley 136257a9adfeSMatthew G Knepley *is = PETSC_NULL; 136357a9adfeSMatthew G Knepley while (ilink) { 136457a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 136557a9adfeSMatthew G Knepley if (found) { 136657a9adfeSMatthew G Knepley *is = ilink->is; 136757a9adfeSMatthew G Knepley break; 136857a9adfeSMatthew G Knepley } 136957a9adfeSMatthew G Knepley ilink = ilink->next; 137057a9adfeSMatthew G Knepley } 137157a9adfeSMatthew G Knepley } 137257a9adfeSMatthew G Knepley PetscFunctionReturn(0); 137357a9adfeSMatthew G Knepley } 137457a9adfeSMatthew G Knepley 137557a9adfeSMatthew G Knepley #undef __FUNCT__ 137651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 137751f519a2SBarry Smith /*@ 137851f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 137951f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 138051f519a2SBarry Smith 1381ad4df100SBarry Smith Logically Collective on PC 138251f519a2SBarry Smith 138351f519a2SBarry Smith Input Parameters: 138451f519a2SBarry Smith + pc - the preconditioner context 138551f519a2SBarry Smith - bs - the block size 138651f519a2SBarry Smith 138751f519a2SBarry Smith Level: intermediate 138851f519a2SBarry Smith 138951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 139051f519a2SBarry Smith 139151f519a2SBarry Smith @*/ 13927087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 139351f519a2SBarry Smith { 13944ac538c5SBarry Smith PetscErrorCode ierr; 139551f519a2SBarry Smith 139651f519a2SBarry Smith PetscFunctionBegin; 13970700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1398c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13994ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 140051f519a2SBarry Smith PetscFunctionReturn(0); 140151f519a2SBarry Smith } 140251f519a2SBarry Smith 140351f519a2SBarry Smith #undef __FUNCT__ 140469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 14050971522cSBarry Smith /*@C 140669a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 14070971522cSBarry Smith 140869a612a9SBarry Smith Collective on KSP 14090971522cSBarry Smith 14100971522cSBarry Smith Input Parameter: 14110971522cSBarry Smith . pc - the preconditioner context 14120971522cSBarry Smith 14130971522cSBarry Smith Output Parameters: 141413e0d083SBarry Smith + n - the number of splits 141569a612a9SBarry Smith - pc - the array of KSP contexts 14160971522cSBarry Smith 14170971522cSBarry Smith Note: 1418d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1419d32f9abdSBarry Smith (not the KSP just the array that contains them). 14200971522cSBarry Smith 142169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 14220971522cSBarry Smith 1423196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 1424196cc216SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the 1425196cc216SBarry Smith KSP array must be. 1426196cc216SBarry Smith 1427196cc216SBarry Smith 14280971522cSBarry Smith Level: advanced 14290971522cSBarry Smith 14300971522cSBarry Smith .seealso: PCFIELDSPLIT 14310971522cSBarry Smith @*/ 14327087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 14330971522cSBarry Smith { 14344ac538c5SBarry Smith PetscErrorCode ierr; 14350971522cSBarry Smith 14360971522cSBarry Smith PetscFunctionBegin; 14370700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 143813e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14394ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14400971522cSBarry Smith PetscFunctionReturn(0); 14410971522cSBarry Smith } 14420971522cSBarry Smith 1443e69d4d44SBarry Smith #undef __FUNCT__ 1444e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1445e69d4d44SBarry Smith /*@ 1446e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1447a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1448e69d4d44SBarry Smith 1449e69d4d44SBarry Smith Collective on PC 1450e69d4d44SBarry Smith 1451e69d4d44SBarry Smith Input Parameters: 1452e69d4d44SBarry Smith + pc - the preconditioner context 1453e87fab1bSBarry Smith . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default 1454084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1455084e4875SJed Brown 1456e69d4d44SBarry Smith Options Database: 1457e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> default is a11 1458e69d4d44SBarry Smith 1459fd1303e9SJungho Lee Notes: 1460fd1303e9SJungho Lee $ If ptype is 1461fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1462fd1303e9SJungho Lee $ to this function). 1463e87fab1bSBarry Smith $ a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1464fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1465fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1466fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1467fd1303e9SJungho Lee $ preconditioner 1468fd1303e9SJungho Lee 1469e87fab1bSBarry Smith When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense 1470fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1471fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1472fd1303e9SJungho Lee 1473fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1474fd1303e9SJungho Lee the name since it sets a proceedure to use. 1475fd1303e9SJungho Lee 1476e69d4d44SBarry Smith Level: intermediate 1477e69d4d44SBarry Smith 1478fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1479e69d4d44SBarry Smith 1480e69d4d44SBarry Smith @*/ 14817087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1482e69d4d44SBarry Smith { 14834ac538c5SBarry Smith PetscErrorCode ierr; 1484e69d4d44SBarry Smith 1485e69d4d44SBarry Smith PetscFunctionBegin; 14860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14874ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1488e69d4d44SBarry Smith PetscFunctionReturn(0); 1489e69d4d44SBarry Smith } 1490e69d4d44SBarry Smith 1491e69d4d44SBarry Smith EXTERN_C_BEGIN 1492e69d4d44SBarry Smith #undef __FUNCT__ 1493e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14947087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1495e69d4d44SBarry Smith { 1496e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1497084e4875SJed Brown PetscErrorCode ierr; 1498e69d4d44SBarry Smith 1499e69d4d44SBarry Smith PetscFunctionBegin; 1500084e4875SJed Brown jac->schurpre = ptype; 1501084e4875SJed Brown if (pre) { 15026bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1503084e4875SJed Brown jac->schur_user = pre; 1504084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1505084e4875SJed Brown } 1506e69d4d44SBarry Smith PetscFunctionReturn(0); 1507e69d4d44SBarry Smith } 1508e69d4d44SBarry Smith EXTERN_C_END 1509e69d4d44SBarry Smith 151030ad9308SMatthew Knepley #undef __FUNCT__ 1511c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1512ab1df9f5SJed Brown /*@ 1513c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1514ab1df9f5SJed Brown 1515ab1df9f5SJed Brown Collective on PC 1516ab1df9f5SJed Brown 1517ab1df9f5SJed Brown Input Parameters: 1518ab1df9f5SJed Brown + pc - the preconditioner context 1519c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1520ab1df9f5SJed Brown 1521ab1df9f5SJed Brown Options Database: 1522c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1523ab1df9f5SJed Brown 1524ab1df9f5SJed Brown 1525ab1df9f5SJed Brown Level: intermediate 1526ab1df9f5SJed Brown 1527ab1df9f5SJed Brown Notes: 1528ab1df9f5SJed Brown The FULL factorization is 1529ab1df9f5SJed Brown 1530ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1531ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1532ab1df9f5SJed Brown 15336be4592eSBarry 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, 15346be4592eSBarry 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). 1535ab1df9f5SJed Brown 15366be4592eSBarry 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 15376be4592eSBarry 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 15386be4592eSBarry 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, 15396be4592eSBarry 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. 1540ab1df9f5SJed Brown 15416be4592eSBarry 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 15426be4592eSBarry 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). 1543ab1df9f5SJed Brown 1544ab1df9f5SJed Brown References: 1545ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1546ab1df9f5SJed Brown 1547ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1548ab1df9f5SJed Brown 1549ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1550ab1df9f5SJed Brown @*/ 1551c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1552ab1df9f5SJed Brown { 1553ab1df9f5SJed Brown PetscErrorCode ierr; 1554ab1df9f5SJed Brown 1555ab1df9f5SJed Brown PetscFunctionBegin; 1556ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1557c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1558ab1df9f5SJed Brown PetscFunctionReturn(0); 1559ab1df9f5SJed Brown } 1560ab1df9f5SJed Brown 1561ab1df9f5SJed Brown #undef __FUNCT__ 1562c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1563c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1564ab1df9f5SJed Brown { 1565ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1566ab1df9f5SJed Brown 1567ab1df9f5SJed Brown PetscFunctionBegin; 1568ab1df9f5SJed Brown jac->schurfactorization = ftype; 1569ab1df9f5SJed Brown PetscFunctionReturn(0); 1570ab1df9f5SJed Brown } 1571ab1df9f5SJed Brown 1572ab1df9f5SJed Brown #undef __FUNCT__ 157330ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 157430ad9308SMatthew Knepley /*@C 15758c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 157630ad9308SMatthew Knepley 157730ad9308SMatthew Knepley Collective on KSP 157830ad9308SMatthew Knepley 157930ad9308SMatthew Knepley Input Parameter: 158030ad9308SMatthew Knepley . pc - the preconditioner context 158130ad9308SMatthew Knepley 158230ad9308SMatthew Knepley Output Parameters: 1583a04f6461SBarry Smith + A00 - the (0,0) block 1584a04f6461SBarry Smith . A01 - the (0,1) block 1585a04f6461SBarry Smith . A10 - the (1,0) block 1586a04f6461SBarry Smith - A11 - the (1,1) block 158730ad9308SMatthew Knepley 158830ad9308SMatthew Knepley Level: advanced 158930ad9308SMatthew Knepley 159030ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 159130ad9308SMatthew Knepley @*/ 1592a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 159330ad9308SMatthew Knepley { 159430ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 159530ad9308SMatthew Knepley 159630ad9308SMatthew Knepley PetscFunctionBegin; 15970700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1598c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1599a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1600a04f6461SBarry Smith if (A01) *A01 = jac->B; 1601a04f6461SBarry Smith if (A10) *A10 = jac->C; 1602a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 160330ad9308SMatthew Knepley PetscFunctionReturn(0); 160430ad9308SMatthew Knepley } 160530ad9308SMatthew Knepley 160679416396SBarry Smith EXTERN_C_BEGIN 160779416396SBarry Smith #undef __FUNCT__ 160879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 16097087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 161079416396SBarry Smith { 161179416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1612e69d4d44SBarry Smith PetscErrorCode ierr; 161379416396SBarry Smith 161479416396SBarry Smith PetscFunctionBegin; 161579416396SBarry Smith jac->type = type; 16163b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 16173b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 16183b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1619e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1620e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1621c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1622e69d4d44SBarry Smith 16233b224e63SBarry Smith } else { 16243b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 16253b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1626e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 16279e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1628c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 16293b224e63SBarry Smith } 163079416396SBarry Smith PetscFunctionReturn(0); 163179416396SBarry Smith } 163279416396SBarry Smith EXTERN_C_END 163379416396SBarry Smith 163451f519a2SBarry Smith EXTERN_C_BEGIN 163551f519a2SBarry Smith #undef __FUNCT__ 163651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16377087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 163851f519a2SBarry Smith { 163951f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 164051f519a2SBarry Smith 164151f519a2SBarry Smith PetscFunctionBegin; 164265e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 164365e19b50SBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 164451f519a2SBarry Smith jac->bs = bs; 164551f519a2SBarry Smith PetscFunctionReturn(0); 164651f519a2SBarry Smith } 164751f519a2SBarry Smith EXTERN_C_END 164851f519a2SBarry Smith 164979416396SBarry Smith #undef __FUNCT__ 165079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1651bc08b0f1SBarry Smith /*@ 165279416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 165379416396SBarry Smith 165479416396SBarry Smith Collective on PC 165579416396SBarry Smith 165679416396SBarry Smith Input Parameter: 165779416396SBarry Smith . pc - the preconditioner context 165881540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 165979416396SBarry Smith 166079416396SBarry Smith Options Database Key: 1661a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 166279416396SBarry Smith 1663b02e2d75SMatthew G Knepley Level: Intermediate 166479416396SBarry Smith 166579416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 166679416396SBarry Smith 166779416396SBarry Smith .seealso: PCCompositeSetType() 166879416396SBarry Smith 166979416396SBarry Smith @*/ 16707087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 167179416396SBarry Smith { 16724ac538c5SBarry Smith PetscErrorCode ierr; 167379416396SBarry Smith 167479416396SBarry Smith PetscFunctionBegin; 16750700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16764ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 167779416396SBarry Smith PetscFunctionReturn(0); 167879416396SBarry Smith } 167979416396SBarry Smith 1680b02e2d75SMatthew G Knepley #undef __FUNCT__ 1681b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1682b02e2d75SMatthew G Knepley /*@ 1683b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1684b02e2d75SMatthew G Knepley 1685b02e2d75SMatthew G Knepley Not collective 1686b02e2d75SMatthew G Knepley 1687b02e2d75SMatthew G Knepley Input Parameter: 1688b02e2d75SMatthew G Knepley . pc - the preconditioner context 1689b02e2d75SMatthew G Knepley 1690b02e2d75SMatthew G Knepley Output Parameter: 1691b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1692b02e2d75SMatthew G Knepley 1693b02e2d75SMatthew G Knepley Level: Intermediate 1694b02e2d75SMatthew G Knepley 1695b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1696b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1697b02e2d75SMatthew G Knepley @*/ 1698b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1699b02e2d75SMatthew G Knepley { 1700b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1701b02e2d75SMatthew G Knepley 1702b02e2d75SMatthew G Knepley PetscFunctionBegin; 1703b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1704b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1705b02e2d75SMatthew G Knepley *type = jac->type; 1706b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1707b02e2d75SMatthew G Knepley } 1708b02e2d75SMatthew G Knepley 17090971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 17100971522cSBarry Smith /*MC 1711a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1712a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 17130971522cSBarry Smith 1714edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1715edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 17160971522cSBarry Smith 1717a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 171869a612a9SBarry Smith and set the options directly on the resulting KSP object 17190971522cSBarry Smith 17200971522cSBarry Smith Level: intermediate 17210971522cSBarry Smith 172279416396SBarry Smith Options Database Keys: 172381540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 172481540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 172581540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 172681540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 17270f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 1728e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11 1729435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 173079416396SBarry Smith 17315d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 17325d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 17335d4c12cdSJungho Lee 1734c8a0d604SMatthew G Knepley Notes: 1735c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1736d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1737d32f9abdSBarry Smith 1738d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1739d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1740d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1741d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1742d32f9abdSBarry Smith 1743c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1744c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1745c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1746c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1747c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1748a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1749c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1750c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1751c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1752a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1753c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1754c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 17555668aaf4SBarry Smith diag gives 1756c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1757c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 17585668aaf4SBarry 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 1759c8a0d604SMatthew G Knepley $ ( A00 0 ) 1760c8a0d604SMatthew G Knepley $ ( A10 S ) 1761c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1762c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1763c8a0d604SMatthew G Knepley $ ( 0 S ) 1764c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1765e69d4d44SBarry Smith 1766edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1767edf189efSBarry Smith is used automatically for a second block. 1768edf189efSBarry Smith 1769ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1770ff218e97SBarry Smith Generally it should be used with the AIJ format. 1771ff218e97SBarry Smith 1772ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1773ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1774ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 17750716a85fSBarry Smith 1776a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 17770971522cSBarry Smith 17787e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1779e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 17800971522cSBarry Smith M*/ 17810971522cSBarry Smith 17820971522cSBarry Smith EXTERN_C_BEGIN 17830971522cSBarry Smith #undef __FUNCT__ 17840971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 17857087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 17860971522cSBarry Smith { 17870971522cSBarry Smith PetscErrorCode ierr; 17880971522cSBarry Smith PC_FieldSplit *jac; 17890971522cSBarry Smith 17900971522cSBarry Smith PetscFunctionBegin; 179138f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 17920971522cSBarry Smith jac->bs = -1; 17930971522cSBarry Smith jac->nsplits = 0; 17943e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1795e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1796c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 179751f519a2SBarry Smith 17980971522cSBarry Smith pc->data = (void*)jac; 17990971522cSBarry Smith 18000971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1801421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 18020971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1803574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 18040971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 18050971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 18060971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 18070971522cSBarry Smith pc->ops->applyrichardson = 0; 18080971522cSBarry Smith 180969a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 181069a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 18110971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 18120971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1813704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1814704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 181579416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 181679416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 181751f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 181851f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 18190971522cSBarry Smith PetscFunctionReturn(0); 18200971522cSBarry Smith } 18210971522cSBarry Smith EXTERN_C_END 18220971522cSBarry Smith 18230971522cSBarry Smith 1824a541d17aSBarry Smith 1825