1dba47a55SKris Buschelman 2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 33c48a1e8SJed Brown #include <petscdmcomposite.h> /*I "petscdmcomposite.h" I*/ 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 13*e87fab1bSBarry 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; 66*e87fab1bSBarry 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; 162*e87fab1bSBarry Smith case PC_FIELDSPLIT_SCHUR_PRE_A11: 163*e87fab1bSBarry 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 { 168*e87fab1bSBarry 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) { 228*e87fab1bSBarry 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 { 359bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 360bafc1b83SMatthew G Knepley for (j = 0; j < i; ++j) { 361bafc1b83SMatthew G Knepley dms[j] = subdm[j]; 362bafc1b83SMatthew G Knepley } 363bafc1b83SMatthew G Knepley } 3647b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3657b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 366e7c4fc90SDmitry Karpeev if (dms) { 3678b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 368bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3697287d2a3SDmitry Karpeev const char *prefix; 3707287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr); 3717287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix); CHKERRQ(ierr); 3727b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3737b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3747287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr); 375e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 3762fa5ba8aSJed Brown } 3777b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3788b8307b2SJed Brown } 37966ffff09SJed Brown } else { 380521d7252SBarry Smith if (jac->bs <= 0) { 381704ba839SBarry Smith if (pc->pmat) { 382521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 383704ba839SBarry Smith } else { 384704ba839SBarry Smith jac->bs = 1; 385704ba839SBarry Smith } 386521d7252SBarry Smith } 387d32f9abdSBarry Smith 3887287d2a3SDmitry Karpeev ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr); 3896ce1633cSBarry Smith if (stokes) { 3906ce1633cSBarry Smith IS zerodiags,rest; 3916ce1633cSBarry Smith PetscInt nmin,nmax; 3926ce1633cSBarry Smith 3936ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3946ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3956ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 3967287d2a3SDmitry Karpeev if (jac->reset) { 3977287d2a3SDmitry Karpeev jac->head->is = rest; 3987287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 3997287d2a3SDmitry Karpeev } 4007287d2a3SDmitry Karpeev else { 4016ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 4026ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 4037287d2a3SDmitry Karpeev } 404fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 405fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 4066ce1633cSBarry Smith } else { 4077287d2a3SDmitry Karpeev if (jac->reset) 4087287d2a3SDmitry Karpeev SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 4097287d2a3SDmitry Karpeev if (!fieldsplit_default) { 410d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 411d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 4126c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 4136c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 414d32f9abdSBarry Smith } 4157287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 416d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 417db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 4186c924f48SJed Brown char splitname[8]; 4198caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 4205d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 42179416396SBarry Smith } 4225d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 423521d7252SBarry Smith } 42466ffff09SJed Brown } 4256ce1633cSBarry Smith } 426edf189efSBarry Smith } else if (jac->nsplits == 1) { 427edf189efSBarry Smith if (ilink->is) { 428edf189efSBarry Smith IS is2; 429edf189efSBarry Smith PetscInt nmin,nmax; 430edf189efSBarry Smith 431edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 432edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 433db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 434fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 4357b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 436edf189efSBarry Smith } 437d0af7cd3SBarry Smith 438d0af7cd3SBarry Smith 4397b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 44069a612a9SBarry Smith PetscFunctionReturn(0); 44169a612a9SBarry Smith } 44269a612a9SBarry Smith 443514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 444514bf10dSMatthew G Knepley 44569a612a9SBarry Smith #undef __FUNCT__ 44669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 44769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 44869a612a9SBarry Smith { 44969a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 45069a612a9SBarry Smith PetscErrorCode ierr; 4515a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4522c9966d7SBarry Smith PetscInt i,nsplit; 45369a612a9SBarry Smith MatStructure flag = pc->flag; 4545f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 45569a612a9SBarry Smith 45669a612a9SBarry Smith PetscFunctionBegin; 45769a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 45897bbdb24SBarry Smith nsplit = jac->nsplits; 4595a9f2f41SSatish Balay ilink = jac->head; 46097bbdb24SBarry Smith 46197bbdb24SBarry Smith /* get the matrices for each split */ 462704ba839SBarry Smith if (!jac->issetup) { 4631b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 46497bbdb24SBarry Smith 465704ba839SBarry Smith jac->issetup = PETSC_TRUE; 466704ba839SBarry Smith 4675d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4682c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4692c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4702c9966d7SBarry Smith } 47151f519a2SBarry Smith bs = jac->bs; 47297bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4731b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4741b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4751b9fc7fcSBarry Smith if (jac->defaultsplit) { 476704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4775f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 478704ba839SBarry Smith } else if (!ilink->is) { 479ccb205f8SBarry Smith if (ilink->nfields > 1) { 4805f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4815a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4825f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4831b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4841b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4851b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4865f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 48797bbdb24SBarry Smith } 48897bbdb24SBarry Smith } 489d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4905f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 491ccb205f8SBarry Smith } else { 492704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4935f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 494ccb205f8SBarry Smith } 4953e197d65SBarry Smith } 496704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4975f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4985f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 499704ba839SBarry Smith ilink = ilink->next; 5001b9fc7fcSBarry Smith } 5011b9fc7fcSBarry Smith } 5021b9fc7fcSBarry Smith 503704ba839SBarry Smith ilink = jac->head; 50497bbdb24SBarry Smith if (!jac->pmat) { 505bdddcaaaSMatthew G Knepley Vec xtmp; 506bdddcaaaSMatthew G Knepley 507bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 508cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 509bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 510cf502942SBarry Smith for (i=0; i<nsplit; i++) { 511bdddcaaaSMatthew G Knepley MatNullSpace sp; 512bdddcaaaSMatthew G Knepley 513a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 514a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 515a3df900dSMatthew G Knepley if (jac->pmat[i]) { 516a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 517a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 518a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 519a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 520a3df900dSMatthew G Knepley } 521a3df900dSMatthew G Knepley } else { 5225f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 523a3df900dSMatthew G Knepley } 524bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 525bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 526443836d0SMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = PETSC_NULL; 527bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 528bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 529ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 530bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr); 531bafc1b83SMatthew G Knepley if (sp) { 532bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 533bafc1b83SMatthew G Knepley } 534ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 535ed1f0337SMatthew G Knepley if (sp) { 536ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 537ed1f0337SMatthew G Knepley } 538704ba839SBarry Smith ilink = ilink->next; 539cf502942SBarry Smith } 540bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 54197bbdb24SBarry Smith } else { 542cf502942SBarry Smith for (i=0; i<nsplit; i++) { 543a3df900dSMatthew G Knepley Mat pmat; 544a3df900dSMatthew G Knepley 545a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 546a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 547a3df900dSMatthew G Knepley if (!pmat) { 5485f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 549a3df900dSMatthew G Knepley } 550704ba839SBarry Smith ilink = ilink->next; 551cf502942SBarry Smith } 55297bbdb24SBarry Smith } 553519d70e2SJed Brown if (jac->realdiagonal) { 554519d70e2SJed Brown ilink = jac->head; 555519d70e2SJed Brown if (!jac->mat) { 556519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 557519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5585f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 559519d70e2SJed Brown ilink = ilink->next; 560519d70e2SJed Brown } 561519d70e2SJed Brown } else { 562519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5635f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 564519d70e2SJed Brown ilink = ilink->next; 565519d70e2SJed Brown } 566519d70e2SJed Brown } 567519d70e2SJed Brown } else { 568519d70e2SJed Brown jac->mat = jac->pmat; 569519d70e2SJed Brown } 57097bbdb24SBarry Smith 5716c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 57268dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 57368dd23aaSBarry Smith ilink = jac->head; 57468dd23aaSBarry Smith if (!jac->Afield) { 57568dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 57668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5774aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 57868dd23aaSBarry Smith ilink = ilink->next; 57968dd23aaSBarry Smith } 58068dd23aaSBarry Smith } else { 58168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5824aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58368dd23aaSBarry Smith ilink = ilink->next; 58468dd23aaSBarry Smith } 58568dd23aaSBarry Smith } 58668dd23aaSBarry Smith } 58768dd23aaSBarry Smith 5883b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5893b224e63SBarry Smith IS ccis; 5904aa3045dSJed Brown PetscInt rstart,rend; 591093c86ffSJed Brown char lscname[256]; 592093c86ffSJed Brown PetscObject LSC_L; 593e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 59468dd23aaSBarry Smith 595e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 596e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 597e6cab6aaSJed Brown 5983b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5993b224e63SBarry Smith if (jac->schur) { 600a3314f2cSMatthew G Knepley KSP kspA = jac->head->ksp, kspInner = PETSC_NULL, kspUpper = jac->kspupper; 601443836d0SMatthew G Knepley 602fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 6033b224e63SBarry Smith ilink = jac->head; 60449bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6054aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 606fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6073b224e63SBarry Smith ilink = ilink->next; 60849bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6094aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 610fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 611a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 612443836d0SMatthew G Knepley if (kspA != kspInner) { 613443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 614443836d0SMatthew G Knepley } 615443836d0SMatthew G Knepley if (kspUpper != kspA) { 616443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 617443836d0SMatthew G Knepley } 618084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 6193b224e63SBarry Smith } else { 6201cee3971SBarry Smith KSP ksp; 621bafc1b83SMatthew G Knepley const char *Dprefix; 622bafc1b83SMatthew G Knepley char schurprefix[256]; 623514bf10dSMatthew G Knepley char schurtestoption[256]; 624bdddcaaaSMatthew G Knepley MatNullSpace sp; 625514bf10dSMatthew G Knepley PetscBool flg; 6263b224e63SBarry Smith 627a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 6283b224e63SBarry Smith ilink = jac->head; 62949bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6304aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 631fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6323b224e63SBarry Smith ilink = ilink->next; 63349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6344aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 635fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 63620252d06SBarry Smith 63720252d06SBarry Smith /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */ 63820252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 63920252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 64020252d06SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr); 64120252d06SBarry Smith ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 64220252d06SBarry Smith /* Indent this deeper to emphasize the "inner" nature of this solver. */ 64320252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr); 64420252d06SBarry Smith ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr); 64520252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 64620252d06SBarry Smith 647bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 64820252d06SBarry Smith if (sp) { 64920252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 65020252d06SBarry Smith } 65120252d06SBarry Smith 65220252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 653514bf10dSMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 654514bf10dSMatthew G Knepley if (flg) { 655514bf10dSMatthew G Knepley DM dmInner; 656514bf10dSMatthew G Knepley 657514bf10dSMatthew G Knepley /* Set DM for new solver */ 658bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 659bafc1b83SMatthew G Knepley ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr); 6607287d2a3SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 661443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 662443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 663514bf10dSMatthew G Knepley } else { 664514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr); 665514bf10dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 666bafc1b83SMatthew G Knepley } 66720b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 6683b224e63SBarry Smith 669443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 670443836d0SMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 671443836d0SMatthew G Knepley if (flg) { 672443836d0SMatthew G Knepley DM dmInner; 673443836d0SMatthew G Knepley 674443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 675443836d0SMatthew G Knepley ierr = KSPCreate(((PetscObject) pc)->comm, &jac->kspupper);CHKERRQ(ierr); 676443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 677443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 678443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 679443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 680443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 681443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 682443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 683443836d0SMatthew G Knepley } else { 684443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 685443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 686443836d0SMatthew G Knepley } 687443836d0SMatthew G Knepley 6883b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur); CHKERRQ(ierr); 6899005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 69020252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1); CHKERRQ(ierr); 691084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 692084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 6937233a360SDmitry Karpeev PC pcschur; 6947233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 6957233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 696084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 697e69d4d44SBarry Smith } 6987287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr); 6997287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix); CHKERRQ(ierr); 7003b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 70120b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 70220b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 7033b224e63SBarry Smith } 704093c86ffSJed Brown 705093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 7068caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 707093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 708093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 709093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 7108caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 711093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 712093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 713093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 7143b224e63SBarry Smith } else { 71568bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 71697bbdb24SBarry Smith i = 0; 7175a9f2f41SSatish Balay ilink = jac->head; 7185a9f2f41SSatish Balay while (ilink) { 719519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 7203b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 721c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 72297bbdb24SBarry Smith i++; 7235a9f2f41SSatish Balay ilink = ilink->next; 7240971522cSBarry Smith } 7253b224e63SBarry Smith } 7263b224e63SBarry Smith 727c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 7280971522cSBarry Smith PetscFunctionReturn(0); 7290971522cSBarry Smith } 7300971522cSBarry Smith 7315a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 732ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 733ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 7345a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 735ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 736ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 73779416396SBarry Smith 7380971522cSBarry Smith #undef __FUNCT__ 7393b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7403b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7413b224e63SBarry Smith { 7423b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7433b224e63SBarry Smith PetscErrorCode ierr; 7443b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 745443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7463b224e63SBarry Smith 7473b224e63SBarry Smith PetscFunctionBegin; 748c5d2311dSJed Brown switch (jac->schurfactorization) { 749c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 750a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 751c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 752c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 753c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 754443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 755c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 756c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 757c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 758c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 759c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 760c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 761c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 762c5d2311dSJed Brown break; 763c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 764a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 765c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 766c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 767443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 768c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 769c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 770c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 771c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 772c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 773c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 774c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 775c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 776c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 777c5d2311dSJed Brown break; 778c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 779a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 780c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 781c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 782c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 783c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 784c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 785c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 786c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 787c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 788443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 789c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 790c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 791c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 792c5d2311dSJed Brown break; 793c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 794a04f6461SBarry 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 */ 7953b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7963b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 797443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 7983b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 7993b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 8003b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8013b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8023b224e63SBarry Smith 8033b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 8043b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8053b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8063b224e63SBarry Smith 807443836d0SMatthew G Knepley if (kspUpper == kspA) { 8083b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 8093b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 810443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 811443836d0SMatthew G Knepley } else { 812443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 813443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 814443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 815443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 816443836d0SMatthew G Knepley } 8173b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8183b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 819c5d2311dSJed Brown } 8203b224e63SBarry Smith PetscFunctionReturn(0); 8213b224e63SBarry Smith } 8223b224e63SBarry Smith 8233b224e63SBarry Smith #undef __FUNCT__ 8240971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 8250971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 8260971522cSBarry Smith { 8270971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8280971522cSBarry Smith PetscErrorCode ierr; 8295a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 830939b8a20SBarry Smith PetscInt cnt,bs; 8310971522cSBarry Smith 8320971522cSBarry Smith PetscFunctionBegin; 8334442daceSBarry Smith x->map->bs = jac->bs; 8344442daceSBarry Smith y->map->bs = jac->bs; 83551f519a2SBarry Smith CHKMEMQ; 83679416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8371b9fc7fcSBarry Smith if (jac->defaultsplit) { 838939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 839939b8a20SBarry 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); 840939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 841939b8a20SBarry 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); 8420971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8435a9f2f41SSatish Balay while (ilink) { 8445a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8455a9f2f41SSatish Balay ilink = ilink->next; 8460971522cSBarry Smith } 8470971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8481b9fc7fcSBarry Smith } else { 849efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8505a9f2f41SSatish Balay while (ilink) { 8515a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8525a9f2f41SSatish Balay ilink = ilink->next; 8531b9fc7fcSBarry Smith } 8541b9fc7fcSBarry Smith } 85516913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 85679416396SBarry Smith if (!jac->w1) { 85779416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 85879416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 85979416396SBarry Smith } 860efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8615a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8623e197d65SBarry Smith cnt = 1; 8635a9f2f41SSatish Balay while (ilink->next) { 8645a9f2f41SSatish Balay ilink = ilink->next; 8653e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8663e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8673e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8683e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8693e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8703e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8713e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8723e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8733e197d65SBarry Smith } 87451f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 87511755939SBarry Smith cnt -= 2; 87651f519a2SBarry Smith while (ilink->previous) { 87751f519a2SBarry Smith ilink = ilink->previous; 87811755939SBarry Smith /* compute the residual only over the part of the vector needed */ 87911755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 88011755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 88111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88311755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 88411755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88511755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88651f519a2SBarry Smith } 88711755939SBarry Smith } 88865e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 88951f519a2SBarry Smith CHKMEMQ; 8900971522cSBarry Smith PetscFunctionReturn(0); 8910971522cSBarry Smith } 8920971522cSBarry Smith 893421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 894ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 895ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 896421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 897ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 898ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 899421e10b8SBarry Smith 900421e10b8SBarry Smith #undef __FUNCT__ 9018c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 902421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 903421e10b8SBarry Smith { 904421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 905421e10b8SBarry Smith PetscErrorCode ierr; 906421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 907939b8a20SBarry Smith PetscInt bs; 908421e10b8SBarry Smith 909421e10b8SBarry Smith PetscFunctionBegin; 910421e10b8SBarry Smith CHKMEMQ; 911421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 912421e10b8SBarry Smith if (jac->defaultsplit) { 913939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 914939b8a20SBarry 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); 915939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 916939b8a20SBarry 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); 917421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 918421e10b8SBarry Smith while (ilink) { 919421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 920421e10b8SBarry Smith ilink = ilink->next; 921421e10b8SBarry Smith } 922421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 923421e10b8SBarry Smith } else { 924421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 925421e10b8SBarry Smith while (ilink) { 926421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 927421e10b8SBarry Smith ilink = ilink->next; 928421e10b8SBarry Smith } 929421e10b8SBarry Smith } 930421e10b8SBarry Smith } else { 931421e10b8SBarry Smith if (!jac->w1) { 932421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 933421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 934421e10b8SBarry Smith } 935421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 936421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 937421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 938421e10b8SBarry Smith while (ilink->next) { 939421e10b8SBarry Smith ilink = ilink->next; 9409989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 941421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 942421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 943421e10b8SBarry Smith } 944421e10b8SBarry Smith while (ilink->previous) { 945421e10b8SBarry Smith ilink = ilink->previous; 9469989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 947421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 948421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 949421e10b8SBarry Smith } 950421e10b8SBarry Smith } else { 951421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 952421e10b8SBarry Smith ilink = ilink->next; 953421e10b8SBarry Smith } 954421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 955421e10b8SBarry Smith while (ilink->previous) { 956421e10b8SBarry Smith ilink = ilink->previous; 9579989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 958421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 959421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 960421e10b8SBarry Smith } 961421e10b8SBarry Smith } 962421e10b8SBarry Smith } 963421e10b8SBarry Smith CHKMEMQ; 964421e10b8SBarry Smith PetscFunctionReturn(0); 965421e10b8SBarry Smith } 966421e10b8SBarry Smith 9670971522cSBarry Smith #undef __FUNCT__ 968574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 969574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9700971522cSBarry Smith { 9710971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9720971522cSBarry Smith PetscErrorCode ierr; 9735a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9740971522cSBarry Smith 9750971522cSBarry Smith PetscFunctionBegin; 9765a9f2f41SSatish Balay while (ilink) { 977574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 978fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 979fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 980443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 981fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 982fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 983b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9845a9f2f41SSatish Balay next = ilink->next; 9855a9f2f41SSatish Balay ilink = next; 9860971522cSBarry Smith } 98705b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 988574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 989574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 990574deadeSBarry Smith } else if (jac->mat) { 991574deadeSBarry Smith jac->mat = PETSC_NULL; 992574deadeSBarry Smith } 99397bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 99468dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 9956bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 9966bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 9976bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 9986bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 9996bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 1000d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 10016bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 10026bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 100363ec74ffSBarry Smith jac->reset = PETSC_TRUE; 1004574deadeSBarry Smith PetscFunctionReturn(0); 1005574deadeSBarry Smith } 1006574deadeSBarry Smith 1007574deadeSBarry Smith #undef __FUNCT__ 1008574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1009574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1010574deadeSBarry Smith { 1011574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1012574deadeSBarry Smith PetscErrorCode ierr; 1013574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 1014574deadeSBarry Smith 1015574deadeSBarry Smith PetscFunctionBegin; 1016574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 1017574deadeSBarry Smith while (ilink) { 10186bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 1019574deadeSBarry Smith next = ilink->next; 1020574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 1021574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 10225d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 1023574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 1024574deadeSBarry Smith ilink = next; 1025574deadeSBarry Smith } 1026574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1027c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 10287b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 10297b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 10307b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 10317b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 10327b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 1033ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 1034c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 10350971522cSBarry Smith PetscFunctionReturn(0); 10360971522cSBarry Smith } 10370971522cSBarry Smith 10380971522cSBarry Smith #undef __FUNCT__ 10390971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10400971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10410971522cSBarry Smith { 10421b9fc7fcSBarry Smith PetscErrorCode ierr; 10436c924f48SJed Brown PetscInt bs; 1044bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10459dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10463b224e63SBarry Smith PCCompositeType ctype; 1047e7c4fc90SDmitry Karpeev DM ddm; 1048e7c4fc90SDmitry Karpeev char ddm_name[1025]; 10491b9fc7fcSBarry Smith 10500971522cSBarry Smith PetscFunctionBegin; 10511b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 1052e7c4fc90SDmitry Karpeev if (pc->dm) { 1053e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 1054e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 1055731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 1056e7c4fc90SDmitry Karpeev if (flg) { 105716621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 1058e7c4fc90SDmitry Karpeev if (!ddm) { 105916621825SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 1060e7c4fc90SDmitry Karpeev } 106116621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 1062e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 1063e7c4fc90SDmitry Karpeev } 1064e7c4fc90SDmitry Karpeev } 1065acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 106651f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 106751f519a2SBarry Smith if (flg) { 106851f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 106951f519a2SBarry Smith } 1070704ba839SBarry Smith 1071435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 1072c0adfefeSBarry Smith if (stokes) { 1073c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1074c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1075c0adfefeSBarry Smith } 1076c0adfefeSBarry Smith 10773b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10783b224e63SBarry Smith if (flg) { 10793b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10803b224e63SBarry Smith } 1081c30613efSMatthew Knepley /* Only setup fields once */ 1082c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1083d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1084d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10856c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10866c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1087d32f9abdSBarry Smith } 1088c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1089c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1090c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 1091c9c6ffaaSJed 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); 1092084e4875SJed 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); 1093c5d2311dSJed Brown } 10941b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10950971522cSBarry Smith PetscFunctionReturn(0); 10960971522cSBarry Smith } 10970971522cSBarry Smith 10980971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10990971522cSBarry Smith 11000971522cSBarry Smith EXTERN_C_BEGIN 11010971522cSBarry Smith #undef __FUNCT__ 11020971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 11035d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11040971522cSBarry Smith { 110597bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11060971522cSBarry Smith PetscErrorCode ierr; 11075a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 110869a612a9SBarry Smith char prefix[128]; 11095d4c12cdSJungho Lee PetscInt i; 11100971522cSBarry Smith 11110971522cSBarry Smith PetscFunctionBegin; 11126c924f48SJed Brown if (jac->splitdefined) { 11136c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11146c924f48SJed Brown PetscFunctionReturn(0); 11156c924f48SJed Brown } 111651f519a2SBarry Smith for (i=0; i<n; i++) { 1117e32f2f54SBarry 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); 1118e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 111951f519a2SBarry Smith } 1120704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1121a04f6461SBarry Smith if (splitname) { 1122db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1123a04f6461SBarry Smith } else { 1124a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1125a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1126a04f6461SBarry Smith } 1127704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 11285a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 11295d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 11305d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 11315a9f2f41SSatish Balay ilink->nfields = n; 11325a9f2f41SSatish Balay ilink->next = PETSC_NULL; 11337adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 113420252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11355a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11369005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 113769a612a9SBarry Smith 11388caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 11395a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11400971522cSBarry Smith 11410971522cSBarry Smith if (!next) { 11425a9f2f41SSatish Balay jac->head = ilink; 114351f519a2SBarry Smith ilink->previous = PETSC_NULL; 11440971522cSBarry Smith } else { 11450971522cSBarry Smith while (next->next) { 11460971522cSBarry Smith next = next->next; 11470971522cSBarry Smith } 11485a9f2f41SSatish Balay next->next = ilink; 114951f519a2SBarry Smith ilink->previous = next; 11500971522cSBarry Smith } 11510971522cSBarry Smith jac->nsplits++; 11520971522cSBarry Smith PetscFunctionReturn(0); 11530971522cSBarry Smith } 11540971522cSBarry Smith EXTERN_C_END 11550971522cSBarry Smith 1156e69d4d44SBarry Smith EXTERN_C_BEGIN 1157e69d4d44SBarry Smith #undef __FUNCT__ 1158e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11597087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1160e69d4d44SBarry Smith { 1161e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1162e69d4d44SBarry Smith PetscErrorCode ierr; 1163e69d4d44SBarry Smith 1164e69d4d44SBarry Smith PetscFunctionBegin; 1165519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1166e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1167e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 116813e0d083SBarry Smith if (n) *n = jac->nsplits; 1169e69d4d44SBarry Smith PetscFunctionReturn(0); 1170e69d4d44SBarry Smith } 1171e69d4d44SBarry Smith EXTERN_C_END 11720971522cSBarry Smith 11730971522cSBarry Smith EXTERN_C_BEGIN 11740971522cSBarry Smith #undef __FUNCT__ 117569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11767087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11770971522cSBarry Smith { 11780971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11790971522cSBarry Smith PetscErrorCode ierr; 11800971522cSBarry Smith PetscInt cnt = 0; 11815a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11820971522cSBarry Smith 11830971522cSBarry Smith PetscFunctionBegin; 11845d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11855a9f2f41SSatish Balay while (ilink) { 11865a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11875a9f2f41SSatish Balay ilink = ilink->next; 11880971522cSBarry Smith } 11895d480477SMatthew 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); 119013e0d083SBarry Smith if (n) *n = jac->nsplits; 11910971522cSBarry Smith PetscFunctionReturn(0); 11920971522cSBarry Smith } 11930971522cSBarry Smith EXTERN_C_END 11940971522cSBarry Smith 1195704ba839SBarry Smith EXTERN_C_BEGIN 1196704ba839SBarry Smith #undef __FUNCT__ 1197704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11987087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1199704ba839SBarry Smith { 1200704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1201704ba839SBarry Smith PetscErrorCode ierr; 1202704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1203704ba839SBarry Smith char prefix[128]; 1204704ba839SBarry Smith 1205704ba839SBarry Smith PetscFunctionBegin; 12066c924f48SJed Brown if (jac->splitdefined) { 12076c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 12086c924f48SJed Brown PetscFunctionReturn(0); 12096c924f48SJed Brown } 121016913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1211a04f6461SBarry Smith if (splitname) { 1212db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1213a04f6461SBarry Smith } else { 1214b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1215b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1216a04f6461SBarry Smith } 12171ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1218b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1219b5787286SJed Brown ilink->is = is; 1220b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1221b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1222b5787286SJed Brown ilink->is_col = is; 1223704ba839SBarry Smith ilink->next = PETSC_NULL; 1224704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 122520252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1226704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 12279005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1228704ba839SBarry Smith 12298caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1230704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1231704ba839SBarry Smith 1232704ba839SBarry Smith if (!next) { 1233704ba839SBarry Smith jac->head = ilink; 1234704ba839SBarry Smith ilink->previous = PETSC_NULL; 1235704ba839SBarry Smith } else { 1236704ba839SBarry Smith while (next->next) { 1237704ba839SBarry Smith next = next->next; 1238704ba839SBarry Smith } 1239704ba839SBarry Smith next->next = ilink; 1240704ba839SBarry Smith ilink->previous = next; 1241704ba839SBarry Smith } 1242704ba839SBarry Smith jac->nsplits++; 1243704ba839SBarry Smith 1244704ba839SBarry Smith PetscFunctionReturn(0); 1245704ba839SBarry Smith } 1246704ba839SBarry Smith EXTERN_C_END 1247704ba839SBarry Smith 12480971522cSBarry Smith #undef __FUNCT__ 12490971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12500971522cSBarry Smith /*@ 12510971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12520971522cSBarry Smith 1253ad4df100SBarry Smith Logically Collective on PC 12540971522cSBarry Smith 12550971522cSBarry Smith Input Parameters: 12560971522cSBarry Smith + pc - the preconditioner context 1257a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 12580971522cSBarry Smith . n - the number of fields in this split 1259db4c96c1SJed Brown - fields - the fields in this split 12600971522cSBarry Smith 12610971522cSBarry Smith Level: intermediate 12620971522cSBarry Smith 1263d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1264d32f9abdSBarry Smith 12657287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1266d32f9abdSBarry 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 1267d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1268d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1269d32f9abdSBarry Smith 1270db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1271db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1272db4c96c1SJed Brown 12735d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12745d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12755d4c12cdSJungho Lee available when this routine is called. 12765d4c12cdSJungho Lee 1277d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12780971522cSBarry Smith 12790971522cSBarry Smith @*/ 12805d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12810971522cSBarry Smith { 12824ac538c5SBarry Smith PetscErrorCode ierr; 12830971522cSBarry Smith 12840971522cSBarry Smith PetscFunctionBegin; 12850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1286db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1287db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1288db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12895d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12900971522cSBarry Smith PetscFunctionReturn(0); 12910971522cSBarry Smith } 12920971522cSBarry Smith 12930971522cSBarry Smith #undef __FUNCT__ 1294704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1295704ba839SBarry Smith /*@ 1296704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1297704ba839SBarry Smith 1298ad4df100SBarry Smith Logically Collective on PC 1299704ba839SBarry Smith 1300704ba839SBarry Smith Input Parameters: 1301704ba839SBarry Smith + pc - the preconditioner context 1302a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1303db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1304704ba839SBarry Smith 1305d32f9abdSBarry Smith 1306a6ffb8dbSJed Brown Notes: 1307a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1308a6ffb8dbSJed Brown 1309db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1310db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1311d32f9abdSBarry Smith 1312704ba839SBarry Smith Level: intermediate 1313704ba839SBarry Smith 1314704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1315704ba839SBarry Smith 1316704ba839SBarry Smith @*/ 13177087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1318704ba839SBarry Smith { 13194ac538c5SBarry Smith PetscErrorCode ierr; 1320704ba839SBarry Smith 1321704ba839SBarry Smith PetscFunctionBegin; 13220700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13237b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1324db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 13254ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1326704ba839SBarry Smith PetscFunctionReturn(0); 1327704ba839SBarry Smith } 1328704ba839SBarry Smith 1329704ba839SBarry Smith #undef __FUNCT__ 133057a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 133157a9adfeSMatthew G Knepley /*@ 133257a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 133357a9adfeSMatthew G Knepley 133457a9adfeSMatthew G Knepley Logically Collective on PC 133557a9adfeSMatthew G Knepley 133657a9adfeSMatthew G Knepley Input Parameters: 133757a9adfeSMatthew G Knepley + pc - the preconditioner context 133857a9adfeSMatthew G Knepley - splitname - name of this split 133957a9adfeSMatthew G Knepley 134057a9adfeSMatthew G Knepley Output Parameter: 134157a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 134257a9adfeSMatthew G Knepley 134357a9adfeSMatthew G Knepley Level: intermediate 134457a9adfeSMatthew G Knepley 134557a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 134657a9adfeSMatthew G Knepley 134757a9adfeSMatthew G Knepley @*/ 134857a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 134957a9adfeSMatthew G Knepley { 135057a9adfeSMatthew G Knepley PetscErrorCode ierr; 135157a9adfeSMatthew G Knepley 135257a9adfeSMatthew G Knepley PetscFunctionBegin; 135357a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 135457a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 135557a9adfeSMatthew G Knepley PetscValidPointer(is,3); 135657a9adfeSMatthew G Knepley { 135757a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 135857a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 135957a9adfeSMatthew G Knepley PetscBool found; 136057a9adfeSMatthew G Knepley 136157a9adfeSMatthew G Knepley *is = PETSC_NULL; 136257a9adfeSMatthew G Knepley while(ilink) { 136357a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 136457a9adfeSMatthew G Knepley if (found) { 136557a9adfeSMatthew G Knepley *is = ilink->is; 136657a9adfeSMatthew G Knepley break; 136757a9adfeSMatthew G Knepley } 136857a9adfeSMatthew G Knepley ilink = ilink->next; 136957a9adfeSMatthew G Knepley } 137057a9adfeSMatthew G Knepley } 137157a9adfeSMatthew G Knepley PetscFunctionReturn(0); 137257a9adfeSMatthew G Knepley } 137357a9adfeSMatthew G Knepley 137457a9adfeSMatthew G Knepley #undef __FUNCT__ 137551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 137651f519a2SBarry Smith /*@ 137751f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 137851f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 137951f519a2SBarry Smith 1380ad4df100SBarry Smith Logically Collective on PC 138151f519a2SBarry Smith 138251f519a2SBarry Smith Input Parameters: 138351f519a2SBarry Smith + pc - the preconditioner context 138451f519a2SBarry Smith - bs - the block size 138551f519a2SBarry Smith 138651f519a2SBarry Smith Level: intermediate 138751f519a2SBarry Smith 138851f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 138951f519a2SBarry Smith 139051f519a2SBarry Smith @*/ 13917087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 139251f519a2SBarry Smith { 13934ac538c5SBarry Smith PetscErrorCode ierr; 139451f519a2SBarry Smith 139551f519a2SBarry Smith PetscFunctionBegin; 13960700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1397c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13984ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 139951f519a2SBarry Smith PetscFunctionReturn(0); 140051f519a2SBarry Smith } 140151f519a2SBarry Smith 140251f519a2SBarry Smith #undef __FUNCT__ 140369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 14040971522cSBarry Smith /*@C 140569a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 14060971522cSBarry Smith 140769a612a9SBarry Smith Collective on KSP 14080971522cSBarry Smith 14090971522cSBarry Smith Input Parameter: 14100971522cSBarry Smith . pc - the preconditioner context 14110971522cSBarry Smith 14120971522cSBarry Smith Output Parameters: 141313e0d083SBarry Smith + n - the number of splits 141469a612a9SBarry Smith - pc - the array of KSP contexts 14150971522cSBarry Smith 14160971522cSBarry Smith Note: 1417d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1418d32f9abdSBarry Smith (not the KSP just the array that contains them). 14190971522cSBarry Smith 142069a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 14210971522cSBarry Smith 1422196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 1423196cc216SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the 1424196cc216SBarry Smith KSP array must be. 1425196cc216SBarry Smith 1426196cc216SBarry Smith 14270971522cSBarry Smith Level: advanced 14280971522cSBarry Smith 14290971522cSBarry Smith .seealso: PCFIELDSPLIT 14300971522cSBarry Smith @*/ 14317087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 14320971522cSBarry Smith { 14334ac538c5SBarry Smith PetscErrorCode ierr; 14340971522cSBarry Smith 14350971522cSBarry Smith PetscFunctionBegin; 14360700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 143713e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14384ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14390971522cSBarry Smith PetscFunctionReturn(0); 14400971522cSBarry Smith } 14410971522cSBarry Smith 1442e69d4d44SBarry Smith #undef __FUNCT__ 1443e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1444e69d4d44SBarry Smith /*@ 1445e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1446a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1447e69d4d44SBarry Smith 1448e69d4d44SBarry Smith Collective on PC 1449e69d4d44SBarry Smith 1450e69d4d44SBarry Smith Input Parameters: 1451e69d4d44SBarry Smith + pc - the preconditioner context 1452*e87fab1bSBarry Smith . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default 1453084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1454084e4875SJed Brown 1455e69d4d44SBarry Smith Options Database: 1456*e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> default is a11 1457e69d4d44SBarry Smith 1458fd1303e9SJungho Lee Notes: 1459fd1303e9SJungho Lee $ If ptype is 1460fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1461fd1303e9SJungho Lee $ to this function). 1462*e87fab1bSBarry Smith $ a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1463fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1464fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1465fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1466fd1303e9SJungho Lee $ preconditioner 1467fd1303e9SJungho Lee 1468*e87fab1bSBarry Smith When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense 1469fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1470fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1471fd1303e9SJungho Lee 1472fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1473fd1303e9SJungho Lee the name since it sets a proceedure to use. 1474fd1303e9SJungho Lee 1475e69d4d44SBarry Smith Level: intermediate 1476e69d4d44SBarry Smith 1477fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1478e69d4d44SBarry Smith 1479e69d4d44SBarry Smith @*/ 14807087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1481e69d4d44SBarry Smith { 14824ac538c5SBarry Smith PetscErrorCode ierr; 1483e69d4d44SBarry Smith 1484e69d4d44SBarry Smith PetscFunctionBegin; 14850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14864ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1487e69d4d44SBarry Smith PetscFunctionReturn(0); 1488e69d4d44SBarry Smith } 1489e69d4d44SBarry Smith 1490e69d4d44SBarry Smith EXTERN_C_BEGIN 1491e69d4d44SBarry Smith #undef __FUNCT__ 1492e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14937087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1494e69d4d44SBarry Smith { 1495e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1496084e4875SJed Brown PetscErrorCode ierr; 1497e69d4d44SBarry Smith 1498e69d4d44SBarry Smith PetscFunctionBegin; 1499084e4875SJed Brown jac->schurpre = ptype; 1500084e4875SJed Brown if (pre) { 15016bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1502084e4875SJed Brown jac->schur_user = pre; 1503084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1504084e4875SJed Brown } 1505e69d4d44SBarry Smith PetscFunctionReturn(0); 1506e69d4d44SBarry Smith } 1507e69d4d44SBarry Smith EXTERN_C_END 1508e69d4d44SBarry Smith 150930ad9308SMatthew Knepley #undef __FUNCT__ 1510c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1511ab1df9f5SJed Brown /*@ 1512c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1513ab1df9f5SJed Brown 1514ab1df9f5SJed Brown Collective on PC 1515ab1df9f5SJed Brown 1516ab1df9f5SJed Brown Input Parameters: 1517ab1df9f5SJed Brown + pc - the preconditioner context 1518c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1519ab1df9f5SJed Brown 1520ab1df9f5SJed Brown Options Database: 1521c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1522ab1df9f5SJed Brown 1523ab1df9f5SJed Brown 1524ab1df9f5SJed Brown Level: intermediate 1525ab1df9f5SJed Brown 1526ab1df9f5SJed Brown Notes: 1527ab1df9f5SJed Brown The FULL factorization is 1528ab1df9f5SJed Brown 1529ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1530ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1531ab1df9f5SJed Brown 15326be4592eSBarry 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, 15336be4592eSBarry 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). 1534ab1df9f5SJed Brown 15356be4592eSBarry 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 15366be4592eSBarry 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 15376be4592eSBarry 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, 15386be4592eSBarry 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. 1539ab1df9f5SJed Brown 15406be4592eSBarry 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 15416be4592eSBarry 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). 1542ab1df9f5SJed Brown 1543ab1df9f5SJed Brown References: 1544ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1545ab1df9f5SJed Brown 1546ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1547ab1df9f5SJed Brown 1548ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1549ab1df9f5SJed Brown @*/ 1550c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1551ab1df9f5SJed Brown { 1552ab1df9f5SJed Brown PetscErrorCode ierr; 1553ab1df9f5SJed Brown 1554ab1df9f5SJed Brown PetscFunctionBegin; 1555ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1556c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1557ab1df9f5SJed Brown PetscFunctionReturn(0); 1558ab1df9f5SJed Brown } 1559ab1df9f5SJed Brown 1560ab1df9f5SJed Brown #undef __FUNCT__ 1561c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1562c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1563ab1df9f5SJed Brown { 1564ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1565ab1df9f5SJed Brown 1566ab1df9f5SJed Brown PetscFunctionBegin; 1567ab1df9f5SJed Brown jac->schurfactorization = ftype; 1568ab1df9f5SJed Brown PetscFunctionReturn(0); 1569ab1df9f5SJed Brown } 1570ab1df9f5SJed Brown 1571ab1df9f5SJed Brown #undef __FUNCT__ 157230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 157330ad9308SMatthew Knepley /*@C 15748c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 157530ad9308SMatthew Knepley 157630ad9308SMatthew Knepley Collective on KSP 157730ad9308SMatthew Knepley 157830ad9308SMatthew Knepley Input Parameter: 157930ad9308SMatthew Knepley . pc - the preconditioner context 158030ad9308SMatthew Knepley 158130ad9308SMatthew Knepley Output Parameters: 1582a04f6461SBarry Smith + A00 - the (0,0) block 1583a04f6461SBarry Smith . A01 - the (0,1) block 1584a04f6461SBarry Smith . A10 - the (1,0) block 1585a04f6461SBarry Smith - A11 - the (1,1) block 158630ad9308SMatthew Knepley 158730ad9308SMatthew Knepley Level: advanced 158830ad9308SMatthew Knepley 158930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 159030ad9308SMatthew Knepley @*/ 1591a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 159230ad9308SMatthew Knepley { 159330ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 159430ad9308SMatthew Knepley 159530ad9308SMatthew Knepley PetscFunctionBegin; 15960700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1597c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1598a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1599a04f6461SBarry Smith if (A01) *A01 = jac->B; 1600a04f6461SBarry Smith if (A10) *A10 = jac->C; 1601a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 160230ad9308SMatthew Knepley PetscFunctionReturn(0); 160330ad9308SMatthew Knepley } 160430ad9308SMatthew Knepley 160579416396SBarry Smith EXTERN_C_BEGIN 160679416396SBarry Smith #undef __FUNCT__ 160779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 16087087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 160979416396SBarry Smith { 161079416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1611e69d4d44SBarry Smith PetscErrorCode ierr; 161279416396SBarry Smith 161379416396SBarry Smith PetscFunctionBegin; 161479416396SBarry Smith jac->type = type; 16153b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 16163b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 16173b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1618e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1619e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1620c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1621e69d4d44SBarry Smith 16223b224e63SBarry Smith } else { 16233b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 16243b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1625e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 16269e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1627c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 16283b224e63SBarry Smith } 162979416396SBarry Smith PetscFunctionReturn(0); 163079416396SBarry Smith } 163179416396SBarry Smith EXTERN_C_END 163279416396SBarry Smith 163351f519a2SBarry Smith EXTERN_C_BEGIN 163451f519a2SBarry Smith #undef __FUNCT__ 163551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16367087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 163751f519a2SBarry Smith { 163851f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 163951f519a2SBarry Smith 164051f519a2SBarry Smith PetscFunctionBegin; 164165e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 164265e19b50SBarry 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); 164351f519a2SBarry Smith jac->bs = bs; 164451f519a2SBarry Smith PetscFunctionReturn(0); 164551f519a2SBarry Smith } 164651f519a2SBarry Smith EXTERN_C_END 164751f519a2SBarry Smith 164879416396SBarry Smith #undef __FUNCT__ 164979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1650bc08b0f1SBarry Smith /*@ 165179416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 165279416396SBarry Smith 165379416396SBarry Smith Collective on PC 165479416396SBarry Smith 165579416396SBarry Smith Input Parameter: 165679416396SBarry Smith . pc - the preconditioner context 165781540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 165879416396SBarry Smith 165979416396SBarry Smith Options Database Key: 1660a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 166179416396SBarry Smith 1662b02e2d75SMatthew G Knepley Level: Intermediate 166379416396SBarry Smith 166479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 166579416396SBarry Smith 166679416396SBarry Smith .seealso: PCCompositeSetType() 166779416396SBarry Smith 166879416396SBarry Smith @*/ 16697087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 167079416396SBarry Smith { 16714ac538c5SBarry Smith PetscErrorCode ierr; 167279416396SBarry Smith 167379416396SBarry Smith PetscFunctionBegin; 16740700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16754ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 167679416396SBarry Smith PetscFunctionReturn(0); 167779416396SBarry Smith } 167879416396SBarry Smith 1679b02e2d75SMatthew G Knepley #undef __FUNCT__ 1680b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1681b02e2d75SMatthew G Knepley /*@ 1682b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1683b02e2d75SMatthew G Knepley 1684b02e2d75SMatthew G Knepley Not collective 1685b02e2d75SMatthew G Knepley 1686b02e2d75SMatthew G Knepley Input Parameter: 1687b02e2d75SMatthew G Knepley . pc - the preconditioner context 1688b02e2d75SMatthew G Knepley 1689b02e2d75SMatthew G Knepley Output Parameter: 1690b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1691b02e2d75SMatthew G Knepley 1692b02e2d75SMatthew G Knepley Level: Intermediate 1693b02e2d75SMatthew G Knepley 1694b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1695b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1696b02e2d75SMatthew G Knepley @*/ 1697b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1698b02e2d75SMatthew G Knepley { 1699b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1700b02e2d75SMatthew G Knepley 1701b02e2d75SMatthew G Knepley PetscFunctionBegin; 1702b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1703b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1704b02e2d75SMatthew G Knepley *type = jac->type; 1705b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1706b02e2d75SMatthew G Knepley } 1707b02e2d75SMatthew G Knepley 17080971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 17090971522cSBarry Smith /*MC 1710a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1711a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 17120971522cSBarry Smith 1713edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1714edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 17150971522cSBarry Smith 1716a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 171769a612a9SBarry Smith and set the options directly on the resulting KSP object 17180971522cSBarry Smith 17190971522cSBarry Smith Level: intermediate 17200971522cSBarry Smith 172179416396SBarry Smith Options Database Keys: 172281540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 172381540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 172481540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 172581540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 17260f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 1727*e87fab1bSBarry Smith . -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11 1728435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 172979416396SBarry Smith 17305d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 17315d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 17325d4c12cdSJungho Lee 1733c8a0d604SMatthew G Knepley Notes: 1734c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1735d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1736d32f9abdSBarry Smith 1737d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1738d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1739d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1740d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1741d32f9abdSBarry Smith 1742c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1743c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1744c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1745c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1746c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1747a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1748c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1749c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1750c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1751a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1752c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1753c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 17545668aaf4SBarry Smith diag gives 1755c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1756c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 17575668aaf4SBarry 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 1758c8a0d604SMatthew G Knepley $ ( A00 0 ) 1759c8a0d604SMatthew G Knepley $ ( A10 S ) 1760c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1761c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1762c8a0d604SMatthew G Knepley $ ( 0 S ) 1763c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1764e69d4d44SBarry Smith 1765edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1766edf189efSBarry Smith is used automatically for a second block. 1767edf189efSBarry Smith 1768ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1769ff218e97SBarry Smith Generally it should be used with the AIJ format. 1770ff218e97SBarry Smith 1771ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1772ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1773ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 17740716a85fSBarry Smith 1775a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 17760971522cSBarry Smith 17777e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1778e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 17790971522cSBarry Smith M*/ 17800971522cSBarry Smith 17810971522cSBarry Smith EXTERN_C_BEGIN 17820971522cSBarry Smith #undef __FUNCT__ 17830971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 17847087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 17850971522cSBarry Smith { 17860971522cSBarry Smith PetscErrorCode ierr; 17870971522cSBarry Smith PC_FieldSplit *jac; 17880971522cSBarry Smith 17890971522cSBarry Smith PetscFunctionBegin; 179038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 17910971522cSBarry Smith jac->bs = -1; 17920971522cSBarry Smith jac->nsplits = 0; 17933e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1794e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1795c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 179651f519a2SBarry Smith 17970971522cSBarry Smith pc->data = (void*)jac; 17980971522cSBarry Smith 17990971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1800421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 18010971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1802574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 18030971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 18040971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 18050971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 18060971522cSBarry Smith pc->ops->applyrichardson = 0; 18070971522cSBarry Smith 180869a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 180969a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 18100971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 18110971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1812704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1813704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 181479416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 181579416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 181651f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 181751f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 18180971522cSBarry Smith PetscFunctionReturn(0); 18190971522cSBarry Smith } 18200971522cSBarry Smith EXTERN_C_END 18210971522cSBarry Smith 18220971522cSBarry Smith 1823a541d17aSBarry Smith 1824