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 13f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","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; 66084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: 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; 143*4996c5bdSBarry 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); 149*4996c5bdSBarry 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; 162b02e2d75SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_DIAG: 1633e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of 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 { 1683e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of 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); 212*4996c5bdSBarry Smith } else if (isdraw) { 213*4996c5bdSBarry Smith PetscDraw draw; 214*4996c5bdSBarry Smith PetscReal x,y,w,wd,h; 215*4996c5bdSBarry Smith PetscInt cnt = 2; 216*4996c5bdSBarry Smith char str[32]; 217*4996c5bdSBarry Smith 218*4996c5bdSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 219*4996c5bdSBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 220*4996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 221*4996c5bdSBarry Smith ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,PETSC_NULL,&h);CHKERRQ(ierr); 222*4996c5bdSBarry Smith y -= h; 223*4996c5bdSBarry Smith if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) { 224*4996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_DIAG]);CHKERRQ(ierr); 2253b224e63SBarry Smith } else { 226*4996c5bdSBarry Smith ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr); 227*4996c5bdSBarry Smith } 228*4996c5bdSBarry Smith ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,PETSC_NULL,&h);CHKERRQ(ierr); 229*4996c5bdSBarry Smith y -= h; 230*4996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) cnt++; 231*4996c5bdSBarry Smith w = 2*PetscMin(1.0 - x,x); 232*4996c5bdSBarry Smith wd = w/(cnt + 1); 233*4996c5bdSBarry Smith x = x - wd*(cnt-1)/2.0; 234*4996c5bdSBarry Smith 235*4996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 236*4996c5bdSBarry Smith ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 237*4996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 238*4996c5bdSBarry Smith if (jac->kspupper != jac->head->ksp) { 239*4996c5bdSBarry Smith x += wd; 240*4996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 241*4996c5bdSBarry Smith ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 242*4996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 243*4996c5bdSBarry Smith } 244*4996c5bdSBarry Smith x += wd; 245*4996c5bdSBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 246*4996c5bdSBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 247*4996c5bdSBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 2483b224e63SBarry Smith } 2493b224e63SBarry Smith PetscFunctionReturn(0); 2503b224e63SBarry Smith } 2513b224e63SBarry Smith 2523b224e63SBarry Smith #undef __FUNCT__ 2536c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2546c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2556c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2566c924f48SJed Brown { 2576c924f48SJed Brown PetscErrorCode ierr; 2586c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2595d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2605d4c12cdSJungho Lee PetscBool flg,flg_col; 2615d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2626c924f48SJed Brown 2636c924f48SJed Brown PetscFunctionBegin; 2646c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2655d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2666c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2678caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2688caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2698caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2706c924f48SJed Brown nfields = jac->bs; 27129499fbbSJungho Lee nfields_col = jac->bs; 2726c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2735d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2746c924f48SJed Brown if (!flg) break; 2755d4c12cdSJungho Lee else if (flg && !flg_col) { 2765d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2775d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2785d4c12cdSJungho Lee } 2795d4c12cdSJungho Lee else { 2805d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2815d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2825d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2835d4c12cdSJungho Lee } 2846c924f48SJed Brown } 2856c924f48SJed Brown if (i > 0) { 2866c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2876c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2886c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2896c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2906c924f48SJed Brown } 2916c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2925d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2936c924f48SJed Brown PetscFunctionReturn(0); 2946c924f48SJed Brown } 2956c924f48SJed Brown 2966c924f48SJed Brown #undef __FUNCT__ 29769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 29869a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2990971522cSBarry Smith { 3000971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3010971522cSBarry Smith PetscErrorCode ierr; 3025a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3037287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 3046c924f48SJed Brown PetscInt i; 3050971522cSBarry Smith 3060971522cSBarry Smith PetscFunctionBegin; 3077287d2a3SDmitry Karpeev /* 3087287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 3097287d2a3SDmitry Karpeev Should probably be rewritten. 3107287d2a3SDmitry Karpeev */ 3117287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 312d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 313d0af7cd3SBarry Smith if (pc->dm && !stokes) { 314bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 3150784a22cSJed Brown char **fieldNames; 3167b62db95SJungho Lee IS *fields; 317e7c4fc90SDmitry Karpeev DM *dms; 318bafc1b83SMatthew G Knepley DM subdm[128]; 319bafc1b83SMatthew G Knepley PetscBool flg; 320bafc1b83SMatthew G Knepley 32116621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 322bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 323bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE; ; i++) { 324bafc1b83SMatthew G Knepley PetscInt ifields[128]; 325bafc1b83SMatthew G Knepley IS compField; 326bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 327bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 328bafc1b83SMatthew G Knepley 3298caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 330bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 331bafc1b83SMatthew G Knepley if (!flg) break; 332bafc1b83SMatthew G Knepley if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 333bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 334bafc1b83SMatthew G Knepley if (nfields == 1) { 335bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 3364d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 3374d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 338bafc1b83SMatthew G Knepley } else { 3398caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 340bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 3414d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr); 3424d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 3437287d2a3SDmitry Karpeev } 344bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 345bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 346bafc1b83SMatthew G Knepley f = ifields[j]; 3477b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 3487b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 3497b62db95SJungho Lee } 350bafc1b83SMatthew G Knepley } 351bafc1b83SMatthew G Knepley if (i == 0) { 352bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 353bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 354bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 355bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 356bafc1b83SMatthew G Knepley } 357bafc1b83SMatthew G Knepley } else { 358bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 359bafc1b83SMatthew G Knepley for (j = 0; j < i; ++j) { 360bafc1b83SMatthew G Knepley dms[j] = subdm[j]; 361bafc1b83SMatthew G Knepley } 362bafc1b83SMatthew G Knepley } 3637b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3647b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 365e7c4fc90SDmitry Karpeev if (dms) { 3668b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 367bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3687287d2a3SDmitry Karpeev const char *prefix; 3697287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr); 3707287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix); CHKERRQ(ierr); 3717b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3727b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3737287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr); 374e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 3752fa5ba8aSJed Brown } 3767b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3778b8307b2SJed Brown } 37866ffff09SJed Brown } else { 379521d7252SBarry Smith if (jac->bs <= 0) { 380704ba839SBarry Smith if (pc->pmat) { 381521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 382704ba839SBarry Smith } else { 383704ba839SBarry Smith jac->bs = 1; 384704ba839SBarry Smith } 385521d7252SBarry Smith } 386d32f9abdSBarry Smith 3877287d2a3SDmitry Karpeev ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr); 3886ce1633cSBarry Smith if (stokes) { 3896ce1633cSBarry Smith IS zerodiags,rest; 3906ce1633cSBarry Smith PetscInt nmin,nmax; 3916ce1633cSBarry Smith 3926ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3936ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3946ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 3957287d2a3SDmitry Karpeev if (jac->reset) { 3967287d2a3SDmitry Karpeev jac->head->is = rest; 3977287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 3987287d2a3SDmitry Karpeev } 3997287d2a3SDmitry Karpeev else { 4006ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 4016ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 4027287d2a3SDmitry Karpeev } 403fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 404fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 4056ce1633cSBarry Smith } else { 4067287d2a3SDmitry Karpeev if (jac->reset) 4077287d2a3SDmitry Karpeev SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 4087287d2a3SDmitry Karpeev if (!fieldsplit_default) { 409d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 410d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 4116c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 4126c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 413d32f9abdSBarry Smith } 4147287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 415d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 416db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 4176c924f48SJed Brown char splitname[8]; 4188caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 4195d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 42079416396SBarry Smith } 4215d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 422521d7252SBarry Smith } 42366ffff09SJed Brown } 4246ce1633cSBarry Smith } 425edf189efSBarry Smith } else if (jac->nsplits == 1) { 426edf189efSBarry Smith if (ilink->is) { 427edf189efSBarry Smith IS is2; 428edf189efSBarry Smith PetscInt nmin,nmax; 429edf189efSBarry Smith 430edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 431edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 432db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 433fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 4347b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 435edf189efSBarry Smith } 436d0af7cd3SBarry Smith 437d0af7cd3SBarry Smith 4387b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 43969a612a9SBarry Smith PetscFunctionReturn(0); 44069a612a9SBarry Smith } 44169a612a9SBarry Smith 442514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 443514bf10dSMatthew G Knepley 44469a612a9SBarry Smith #undef __FUNCT__ 44569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 44669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 44769a612a9SBarry Smith { 44869a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 44969a612a9SBarry Smith PetscErrorCode ierr; 4505a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4512c9966d7SBarry Smith PetscInt i,nsplit; 45269a612a9SBarry Smith MatStructure flag = pc->flag; 4535f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 45469a612a9SBarry Smith 45569a612a9SBarry Smith PetscFunctionBegin; 45669a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 45797bbdb24SBarry Smith nsplit = jac->nsplits; 4585a9f2f41SSatish Balay ilink = jac->head; 45997bbdb24SBarry Smith 46097bbdb24SBarry Smith /* get the matrices for each split */ 461704ba839SBarry Smith if (!jac->issetup) { 4621b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 46397bbdb24SBarry Smith 464704ba839SBarry Smith jac->issetup = PETSC_TRUE; 465704ba839SBarry Smith 4665d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4672c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4682c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4692c9966d7SBarry Smith } 47051f519a2SBarry Smith bs = jac->bs; 47197bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4721b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4731b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4741b9fc7fcSBarry Smith if (jac->defaultsplit) { 475704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4765f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 477704ba839SBarry Smith } else if (!ilink->is) { 478ccb205f8SBarry Smith if (ilink->nfields > 1) { 4795f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4805a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4815f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4821b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4831b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4841b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4855f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 48697bbdb24SBarry Smith } 48797bbdb24SBarry Smith } 488d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4895f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 490ccb205f8SBarry Smith } else { 491704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4925f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 493ccb205f8SBarry Smith } 4943e197d65SBarry Smith } 495704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4965f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4975f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 498704ba839SBarry Smith ilink = ilink->next; 4991b9fc7fcSBarry Smith } 5001b9fc7fcSBarry Smith } 5011b9fc7fcSBarry Smith 502704ba839SBarry Smith ilink = jac->head; 50397bbdb24SBarry Smith if (!jac->pmat) { 504bdddcaaaSMatthew G Knepley Vec xtmp; 505bdddcaaaSMatthew G Knepley 506bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 507cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 508bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 509cf502942SBarry Smith for (i=0; i<nsplit; i++) { 510bdddcaaaSMatthew G Knepley MatNullSpace sp; 511bdddcaaaSMatthew G Knepley 512a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 513a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 514a3df900dSMatthew G Knepley if (jac->pmat[i]) { 515a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 516a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 517a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 518a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 519a3df900dSMatthew G Knepley } 520a3df900dSMatthew G Knepley } else { 5215f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 522a3df900dSMatthew G Knepley } 523bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 524bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 525443836d0SMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = PETSC_NULL; 526bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 527bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 528ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 529bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr); 530bafc1b83SMatthew G Knepley if (sp) { 531bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 532bafc1b83SMatthew G Knepley } 533ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 534ed1f0337SMatthew G Knepley if (sp) { 535ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 536ed1f0337SMatthew G Knepley } 537704ba839SBarry Smith ilink = ilink->next; 538cf502942SBarry Smith } 539bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 54097bbdb24SBarry Smith } else { 541cf502942SBarry Smith for (i=0; i<nsplit; i++) { 542a3df900dSMatthew G Knepley Mat pmat; 543a3df900dSMatthew G Knepley 544a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 545a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 546a3df900dSMatthew G Knepley if (!pmat) { 5475f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 548a3df900dSMatthew G Knepley } 549704ba839SBarry Smith ilink = ilink->next; 550cf502942SBarry Smith } 55197bbdb24SBarry Smith } 552519d70e2SJed Brown if (jac->realdiagonal) { 553519d70e2SJed Brown ilink = jac->head; 554519d70e2SJed Brown if (!jac->mat) { 555519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 556519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5575f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 558519d70e2SJed Brown ilink = ilink->next; 559519d70e2SJed Brown } 560519d70e2SJed Brown } else { 561519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5625f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 563519d70e2SJed Brown ilink = ilink->next; 564519d70e2SJed Brown } 565519d70e2SJed Brown } 566519d70e2SJed Brown } else { 567519d70e2SJed Brown jac->mat = jac->pmat; 568519d70e2SJed Brown } 56997bbdb24SBarry Smith 5706c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 57168dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 57268dd23aaSBarry Smith ilink = jac->head; 57368dd23aaSBarry Smith if (!jac->Afield) { 57468dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 57568dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5764aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 57768dd23aaSBarry Smith ilink = ilink->next; 57868dd23aaSBarry Smith } 57968dd23aaSBarry Smith } else { 58068dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5814aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 58268dd23aaSBarry Smith ilink = ilink->next; 58368dd23aaSBarry Smith } 58468dd23aaSBarry Smith } 58568dd23aaSBarry Smith } 58668dd23aaSBarry Smith 5873b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5883b224e63SBarry Smith IS ccis; 5894aa3045dSJed Brown PetscInt rstart,rend; 590093c86ffSJed Brown char lscname[256]; 591093c86ffSJed Brown PetscObject LSC_L; 592e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 59368dd23aaSBarry Smith 594e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 595e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 596e6cab6aaSJed Brown 5973b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5983b224e63SBarry Smith if (jac->schur) { 599a3314f2cSMatthew G Knepley KSP kspA = jac->head->ksp, kspInner = PETSC_NULL, kspUpper = jac->kspupper; 600443836d0SMatthew G Knepley 601fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 6023b224e63SBarry Smith ilink = jac->head; 60349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6044aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 605fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6063b224e63SBarry Smith ilink = ilink->next; 60749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6084aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 609fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 610a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 611443836d0SMatthew G Knepley if (kspA != kspInner) { 612443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 613443836d0SMatthew G Knepley } 614443836d0SMatthew G Knepley if (kspUpper != kspA) { 615443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 616443836d0SMatthew G Knepley } 617084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 6183b224e63SBarry Smith } else { 6191cee3971SBarry Smith KSP ksp; 620bafc1b83SMatthew G Knepley const char *Dprefix; 621bafc1b83SMatthew G Knepley char schurprefix[256]; 622514bf10dSMatthew G Knepley char schurtestoption[256]; 623bdddcaaaSMatthew G Knepley MatNullSpace sp; 624514bf10dSMatthew G Knepley PetscBool flg; 6253b224e63SBarry Smith 626a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 6273b224e63SBarry Smith ilink = jac->head; 62849bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6294aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 630fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 6313b224e63SBarry Smith ilink = ilink->next; 63249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 6334aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 634fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 63520252d06SBarry Smith 63620252d06SBarry Smith /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */ 63720252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 63820252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 63920252d06SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr); 64020252d06SBarry Smith ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 64120252d06SBarry Smith /* Indent this deeper to emphasize the "inner" nature of this solver. */ 64220252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr); 64320252d06SBarry Smith ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr); 64420252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 64520252d06SBarry Smith 646bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 64720252d06SBarry Smith if (sp) { 64820252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 64920252d06SBarry Smith } 65020252d06SBarry Smith 65120252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 652514bf10dSMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 653514bf10dSMatthew G Knepley if (flg) { 654514bf10dSMatthew G Knepley DM dmInner; 655514bf10dSMatthew G Knepley 656514bf10dSMatthew G Knepley /* Set DM for new solver */ 657bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 658bafc1b83SMatthew G Knepley ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr); 6597287d2a3SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 660443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 661443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 662514bf10dSMatthew G Knepley } else { 663514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr); 664514bf10dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 665bafc1b83SMatthew G Knepley } 66620b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 6673b224e63SBarry Smith 668443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 669443836d0SMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 670443836d0SMatthew G Knepley if (flg) { 671443836d0SMatthew G Knepley DM dmInner; 672443836d0SMatthew G Knepley 673443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 674443836d0SMatthew G Knepley ierr = KSPCreate(((PetscObject) pc)->comm, &jac->kspupper);CHKERRQ(ierr); 675443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 676443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 677443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 678443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 679443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 680443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 681443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 682443836d0SMatthew G Knepley } else { 683443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 684443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 685443836d0SMatthew G Knepley } 686443836d0SMatthew G Knepley 6873b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur); CHKERRQ(ierr); 6889005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 68920252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1); CHKERRQ(ierr); 690084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 691084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 6927233a360SDmitry Karpeev PC pcschur; 6937233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 6947233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 695084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 696e69d4d44SBarry Smith } 6977287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr); 6987287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix); CHKERRQ(ierr); 6993b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 70020b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 70120b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 7023b224e63SBarry Smith } 703093c86ffSJed Brown 704093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 7058caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 706093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 707093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 708093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 7098caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 710093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 711093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 712093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 7133b224e63SBarry Smith } else { 71468bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 71597bbdb24SBarry Smith i = 0; 7165a9f2f41SSatish Balay ilink = jac->head; 7175a9f2f41SSatish Balay while (ilink) { 718519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 7193b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 720c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 72197bbdb24SBarry Smith i++; 7225a9f2f41SSatish Balay ilink = ilink->next; 7230971522cSBarry Smith } 7243b224e63SBarry Smith } 7253b224e63SBarry Smith 726c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 7270971522cSBarry Smith PetscFunctionReturn(0); 7280971522cSBarry Smith } 7290971522cSBarry Smith 7305a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 731ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 732ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 7335a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 734ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 735ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 73679416396SBarry Smith 7370971522cSBarry Smith #undef __FUNCT__ 7383b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7393b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7403b224e63SBarry Smith { 7413b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7423b224e63SBarry Smith PetscErrorCode ierr; 7433b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 744443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7453b224e63SBarry Smith 7463b224e63SBarry Smith PetscFunctionBegin; 747c5d2311dSJed Brown switch (jac->schurfactorization) { 748c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 749a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 750c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 751c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 752c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 753443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 754c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 755c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 756c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 757c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 758c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 759c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 760c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 761c5d2311dSJed Brown break; 762c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 763a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 764c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 765c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 766443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 767c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 768c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 769c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 770c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 771c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 772c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 773c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 774c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 775c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 776c5d2311dSJed Brown break; 777c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 778a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 779c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 780c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 781c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 782c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 783c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 784c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 785c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 786c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 787443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 788c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 789c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 790c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 791c5d2311dSJed Brown break; 792c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 793a04f6461SBarry 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 */ 7943b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7953b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 796443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 7973b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 7983b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 7993b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8003b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8013b224e63SBarry Smith 8023b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 8033b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8043b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8053b224e63SBarry Smith 806443836d0SMatthew G Knepley if (kspUpper == kspA) { 8073b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 8083b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 809443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 810443836d0SMatthew G Knepley } else { 811443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 812443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 813443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 814443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 815443836d0SMatthew G Knepley } 8163b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8173b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 818c5d2311dSJed Brown } 8193b224e63SBarry Smith PetscFunctionReturn(0); 8203b224e63SBarry Smith } 8213b224e63SBarry Smith 8223b224e63SBarry Smith #undef __FUNCT__ 8230971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 8240971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 8250971522cSBarry Smith { 8260971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 8270971522cSBarry Smith PetscErrorCode ierr; 8285a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 829939b8a20SBarry Smith PetscInt cnt,bs; 8300971522cSBarry Smith 8310971522cSBarry Smith PetscFunctionBegin; 8324442daceSBarry Smith x->map->bs = jac->bs; 8334442daceSBarry Smith y->map->bs = jac->bs; 83451f519a2SBarry Smith CHKMEMQ; 83579416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8361b9fc7fcSBarry Smith if (jac->defaultsplit) { 837939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 838939b8a20SBarry 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); 839939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 840939b8a20SBarry 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); 8410971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8425a9f2f41SSatish Balay while (ilink) { 8435a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8445a9f2f41SSatish Balay ilink = ilink->next; 8450971522cSBarry Smith } 8460971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8471b9fc7fcSBarry Smith } else { 848efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8495a9f2f41SSatish Balay while (ilink) { 8505a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8515a9f2f41SSatish Balay ilink = ilink->next; 8521b9fc7fcSBarry Smith } 8531b9fc7fcSBarry Smith } 85416913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 85579416396SBarry Smith if (!jac->w1) { 85679416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 85779416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 85879416396SBarry Smith } 859efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8605a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8613e197d65SBarry Smith cnt = 1; 8625a9f2f41SSatish Balay while (ilink->next) { 8635a9f2f41SSatish Balay ilink = ilink->next; 8643e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8653e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8663e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8673e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8683e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8693e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8703e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8713e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8723e197d65SBarry Smith } 87351f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 87411755939SBarry Smith cnt -= 2; 87551f519a2SBarry Smith while (ilink->previous) { 87651f519a2SBarry Smith ilink = ilink->previous; 87711755939SBarry Smith /* compute the residual only over the part of the vector needed */ 87811755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 87911755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 88011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88211755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 88311755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88411755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88551f519a2SBarry Smith } 88611755939SBarry Smith } 88765e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 88851f519a2SBarry Smith CHKMEMQ; 8890971522cSBarry Smith PetscFunctionReturn(0); 8900971522cSBarry Smith } 8910971522cSBarry Smith 892421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 893ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 894ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 895421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 896ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 897ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 898421e10b8SBarry Smith 899421e10b8SBarry Smith #undef __FUNCT__ 9008c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 901421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 902421e10b8SBarry Smith { 903421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 904421e10b8SBarry Smith PetscErrorCode ierr; 905421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 906939b8a20SBarry Smith PetscInt bs; 907421e10b8SBarry Smith 908421e10b8SBarry Smith PetscFunctionBegin; 909421e10b8SBarry Smith CHKMEMQ; 910421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 911421e10b8SBarry Smith if (jac->defaultsplit) { 912939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 913939b8a20SBarry 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); 914939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 915939b8a20SBarry 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); 916421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 917421e10b8SBarry Smith while (ilink) { 918421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 919421e10b8SBarry Smith ilink = ilink->next; 920421e10b8SBarry Smith } 921421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 922421e10b8SBarry Smith } else { 923421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 924421e10b8SBarry Smith while (ilink) { 925421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 926421e10b8SBarry Smith ilink = ilink->next; 927421e10b8SBarry Smith } 928421e10b8SBarry Smith } 929421e10b8SBarry Smith } else { 930421e10b8SBarry Smith if (!jac->w1) { 931421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 932421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 933421e10b8SBarry Smith } 934421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 935421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 936421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 937421e10b8SBarry Smith while (ilink->next) { 938421e10b8SBarry Smith ilink = ilink->next; 9399989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 940421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 941421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 942421e10b8SBarry Smith } 943421e10b8SBarry Smith while (ilink->previous) { 944421e10b8SBarry Smith ilink = ilink->previous; 9459989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 946421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 947421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 948421e10b8SBarry Smith } 949421e10b8SBarry Smith } else { 950421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 951421e10b8SBarry Smith ilink = ilink->next; 952421e10b8SBarry Smith } 953421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 954421e10b8SBarry Smith while (ilink->previous) { 955421e10b8SBarry Smith ilink = ilink->previous; 9569989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 957421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 958421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 959421e10b8SBarry Smith } 960421e10b8SBarry Smith } 961421e10b8SBarry Smith } 962421e10b8SBarry Smith CHKMEMQ; 963421e10b8SBarry Smith PetscFunctionReturn(0); 964421e10b8SBarry Smith } 965421e10b8SBarry Smith 9660971522cSBarry Smith #undef __FUNCT__ 967574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 968574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9690971522cSBarry Smith { 9700971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9710971522cSBarry Smith PetscErrorCode ierr; 9725a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9730971522cSBarry Smith 9740971522cSBarry Smith PetscFunctionBegin; 9755a9f2f41SSatish Balay while (ilink) { 976574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 977fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 978fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 979443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 980fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 981fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 982b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9835a9f2f41SSatish Balay next = ilink->next; 9845a9f2f41SSatish Balay ilink = next; 9850971522cSBarry Smith } 98605b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 987574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 988574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 989574deadeSBarry Smith } else if (jac->mat) { 990574deadeSBarry Smith jac->mat = PETSC_NULL; 991574deadeSBarry Smith } 99297bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 99368dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 9946bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 9956bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 9966bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 9976bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 9986bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 999d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 10006bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 10016bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 100263ec74ffSBarry Smith jac->reset = PETSC_TRUE; 1003574deadeSBarry Smith PetscFunctionReturn(0); 1004574deadeSBarry Smith } 1005574deadeSBarry Smith 1006574deadeSBarry Smith #undef __FUNCT__ 1007574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1008574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1009574deadeSBarry Smith { 1010574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1011574deadeSBarry Smith PetscErrorCode ierr; 1012574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 1013574deadeSBarry Smith 1014574deadeSBarry Smith PetscFunctionBegin; 1015574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 1016574deadeSBarry Smith while (ilink) { 10176bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 1018574deadeSBarry Smith next = ilink->next; 1019574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 1020574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 10215d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 1022574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 1023574deadeSBarry Smith ilink = next; 1024574deadeSBarry Smith } 1025574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 1026c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 10277b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 10287b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 10297b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 10307b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 10317b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 1032ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 1033c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 10340971522cSBarry Smith PetscFunctionReturn(0); 10350971522cSBarry Smith } 10360971522cSBarry Smith 10370971522cSBarry Smith #undef __FUNCT__ 10380971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10390971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10400971522cSBarry Smith { 10411b9fc7fcSBarry Smith PetscErrorCode ierr; 10426c924f48SJed Brown PetscInt bs; 1043bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10449dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10453b224e63SBarry Smith PCCompositeType ctype; 1046e7c4fc90SDmitry Karpeev DM ddm; 1047e7c4fc90SDmitry Karpeev char ddm_name[1025]; 10481b9fc7fcSBarry Smith 10490971522cSBarry Smith PetscFunctionBegin; 10501b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 1051e7c4fc90SDmitry Karpeev if (pc->dm) { 1052e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 1053e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 1054731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 1055e7c4fc90SDmitry Karpeev if (flg) { 105616621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 1057e7c4fc90SDmitry Karpeev if (!ddm) { 105816621825SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 1059e7c4fc90SDmitry Karpeev } 106016621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 1061e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 1062e7c4fc90SDmitry Karpeev } 1063e7c4fc90SDmitry Karpeev } 1064acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 106551f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 106651f519a2SBarry Smith if (flg) { 106751f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 106851f519a2SBarry Smith } 1069704ba839SBarry Smith 1070435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 1071c0adfefeSBarry Smith if (stokes) { 1072c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1073c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1074c0adfefeSBarry Smith } 1075c0adfefeSBarry Smith 10763b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10773b224e63SBarry Smith if (flg) { 10783b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10793b224e63SBarry Smith } 1080c30613efSMatthew Knepley /* Only setup fields once */ 1081c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1082d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1083d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10846c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10856c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1086d32f9abdSBarry Smith } 1087c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1088c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1089c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 1090c9c6ffaaSJed 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); 1091084e4875SJed 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); 1092c5d2311dSJed Brown } 10931b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10940971522cSBarry Smith PetscFunctionReturn(0); 10950971522cSBarry Smith } 10960971522cSBarry Smith 10970971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10980971522cSBarry Smith 10990971522cSBarry Smith EXTERN_C_BEGIN 11000971522cSBarry Smith #undef __FUNCT__ 11010971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 11025d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 11030971522cSBarry Smith { 110497bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11050971522cSBarry Smith PetscErrorCode ierr; 11065a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 110769a612a9SBarry Smith char prefix[128]; 11085d4c12cdSJungho Lee PetscInt i; 11090971522cSBarry Smith 11100971522cSBarry Smith PetscFunctionBegin; 11116c924f48SJed Brown if (jac->splitdefined) { 11126c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11136c924f48SJed Brown PetscFunctionReturn(0); 11146c924f48SJed Brown } 111551f519a2SBarry Smith for (i=0; i<n; i++) { 1116e32f2f54SBarry 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); 1117e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 111851f519a2SBarry Smith } 1119704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1120a04f6461SBarry Smith if (splitname) { 1121db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1122a04f6461SBarry Smith } else { 1123a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1124a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1125a04f6461SBarry Smith } 1126704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 11275a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 11285d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 11295d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 11305a9f2f41SSatish Balay ilink->nfields = n; 11315a9f2f41SSatish Balay ilink->next = PETSC_NULL; 11327adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 113320252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 11345a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11359005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 113669a612a9SBarry Smith 11378caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 11385a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11390971522cSBarry Smith 11400971522cSBarry Smith if (!next) { 11415a9f2f41SSatish Balay jac->head = ilink; 114251f519a2SBarry Smith ilink->previous = PETSC_NULL; 11430971522cSBarry Smith } else { 11440971522cSBarry Smith while (next->next) { 11450971522cSBarry Smith next = next->next; 11460971522cSBarry Smith } 11475a9f2f41SSatish Balay next->next = ilink; 114851f519a2SBarry Smith ilink->previous = next; 11490971522cSBarry Smith } 11500971522cSBarry Smith jac->nsplits++; 11510971522cSBarry Smith PetscFunctionReturn(0); 11520971522cSBarry Smith } 11530971522cSBarry Smith EXTERN_C_END 11540971522cSBarry Smith 1155e69d4d44SBarry Smith EXTERN_C_BEGIN 1156e69d4d44SBarry Smith #undef __FUNCT__ 1157e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11587087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1159e69d4d44SBarry Smith { 1160e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1161e69d4d44SBarry Smith PetscErrorCode ierr; 1162e69d4d44SBarry Smith 1163e69d4d44SBarry Smith PetscFunctionBegin; 1164519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1165e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1166e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 116713e0d083SBarry Smith if (n) *n = jac->nsplits; 1168e69d4d44SBarry Smith PetscFunctionReturn(0); 1169e69d4d44SBarry Smith } 1170e69d4d44SBarry Smith EXTERN_C_END 11710971522cSBarry Smith 11720971522cSBarry Smith EXTERN_C_BEGIN 11730971522cSBarry Smith #undef __FUNCT__ 117469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11757087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11760971522cSBarry Smith { 11770971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11780971522cSBarry Smith PetscErrorCode ierr; 11790971522cSBarry Smith PetscInt cnt = 0; 11805a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11810971522cSBarry Smith 11820971522cSBarry Smith PetscFunctionBegin; 11835d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11845a9f2f41SSatish Balay while (ilink) { 11855a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11865a9f2f41SSatish Balay ilink = ilink->next; 11870971522cSBarry Smith } 11885d480477SMatthew 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); 118913e0d083SBarry Smith if (n) *n = jac->nsplits; 11900971522cSBarry Smith PetscFunctionReturn(0); 11910971522cSBarry Smith } 11920971522cSBarry Smith EXTERN_C_END 11930971522cSBarry Smith 1194704ba839SBarry Smith EXTERN_C_BEGIN 1195704ba839SBarry Smith #undef __FUNCT__ 1196704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11977087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1198704ba839SBarry Smith { 1199704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1200704ba839SBarry Smith PetscErrorCode ierr; 1201704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1202704ba839SBarry Smith char prefix[128]; 1203704ba839SBarry Smith 1204704ba839SBarry Smith PetscFunctionBegin; 12056c924f48SJed Brown if (jac->splitdefined) { 12066c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 12076c924f48SJed Brown PetscFunctionReturn(0); 12086c924f48SJed Brown } 120916913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1210a04f6461SBarry Smith if (splitname) { 1211db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1212a04f6461SBarry Smith } else { 1213b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1214b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1215a04f6461SBarry Smith } 12161ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1217b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1218b5787286SJed Brown ilink->is = is; 1219b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1220b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1221b5787286SJed Brown ilink->is_col = is; 1222704ba839SBarry Smith ilink->next = PETSC_NULL; 1223704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 122420252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1225704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 12269005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1227704ba839SBarry Smith 12288caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1229704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1230704ba839SBarry Smith 1231704ba839SBarry Smith if (!next) { 1232704ba839SBarry Smith jac->head = ilink; 1233704ba839SBarry Smith ilink->previous = PETSC_NULL; 1234704ba839SBarry Smith } else { 1235704ba839SBarry Smith while (next->next) { 1236704ba839SBarry Smith next = next->next; 1237704ba839SBarry Smith } 1238704ba839SBarry Smith next->next = ilink; 1239704ba839SBarry Smith ilink->previous = next; 1240704ba839SBarry Smith } 1241704ba839SBarry Smith jac->nsplits++; 1242704ba839SBarry Smith 1243704ba839SBarry Smith PetscFunctionReturn(0); 1244704ba839SBarry Smith } 1245704ba839SBarry Smith EXTERN_C_END 1246704ba839SBarry Smith 12470971522cSBarry Smith #undef __FUNCT__ 12480971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12490971522cSBarry Smith /*@ 12500971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12510971522cSBarry Smith 1252ad4df100SBarry Smith Logically Collective on PC 12530971522cSBarry Smith 12540971522cSBarry Smith Input Parameters: 12550971522cSBarry Smith + pc - the preconditioner context 1256a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 12570971522cSBarry Smith . n - the number of fields in this split 1258db4c96c1SJed Brown - fields - the fields in this split 12590971522cSBarry Smith 12600971522cSBarry Smith Level: intermediate 12610971522cSBarry Smith 1262d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1263d32f9abdSBarry Smith 12647287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1265d32f9abdSBarry Smith size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean 1266d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1267d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1268d32f9abdSBarry Smith 1269db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1270db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1271db4c96c1SJed Brown 12725d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12735d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12745d4c12cdSJungho Lee available when this routine is called. 12755d4c12cdSJungho Lee 1276d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12770971522cSBarry Smith 12780971522cSBarry Smith @*/ 12795d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12800971522cSBarry Smith { 12814ac538c5SBarry Smith PetscErrorCode ierr; 12820971522cSBarry Smith 12830971522cSBarry Smith PetscFunctionBegin; 12840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1285db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1286db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1287db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12885d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12890971522cSBarry Smith PetscFunctionReturn(0); 12900971522cSBarry Smith } 12910971522cSBarry Smith 12920971522cSBarry Smith #undef __FUNCT__ 1293704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1294704ba839SBarry Smith /*@ 1295704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1296704ba839SBarry Smith 1297ad4df100SBarry Smith Logically Collective on PC 1298704ba839SBarry Smith 1299704ba839SBarry Smith Input Parameters: 1300704ba839SBarry Smith + pc - the preconditioner context 1301a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1302db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1303704ba839SBarry Smith 1304d32f9abdSBarry Smith 1305a6ffb8dbSJed Brown Notes: 1306a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1307a6ffb8dbSJed Brown 1308db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1309db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1310d32f9abdSBarry Smith 1311704ba839SBarry Smith Level: intermediate 1312704ba839SBarry Smith 1313704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1314704ba839SBarry Smith 1315704ba839SBarry Smith @*/ 13167087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1317704ba839SBarry Smith { 13184ac538c5SBarry Smith PetscErrorCode ierr; 1319704ba839SBarry Smith 1320704ba839SBarry Smith PetscFunctionBegin; 13210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13227b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1323db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 13244ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1325704ba839SBarry Smith PetscFunctionReturn(0); 1326704ba839SBarry Smith } 1327704ba839SBarry Smith 1328704ba839SBarry Smith #undef __FUNCT__ 132957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 133057a9adfeSMatthew G Knepley /*@ 133157a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 133257a9adfeSMatthew G Knepley 133357a9adfeSMatthew G Knepley Logically Collective on PC 133457a9adfeSMatthew G Knepley 133557a9adfeSMatthew G Knepley Input Parameters: 133657a9adfeSMatthew G Knepley + pc - the preconditioner context 133757a9adfeSMatthew G Knepley - splitname - name of this split 133857a9adfeSMatthew G Knepley 133957a9adfeSMatthew G Knepley Output Parameter: 134057a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 134157a9adfeSMatthew G Knepley 134257a9adfeSMatthew G Knepley Level: intermediate 134357a9adfeSMatthew G Knepley 134457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 134557a9adfeSMatthew G Knepley 134657a9adfeSMatthew G Knepley @*/ 134757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 134857a9adfeSMatthew G Knepley { 134957a9adfeSMatthew G Knepley PetscErrorCode ierr; 135057a9adfeSMatthew G Knepley 135157a9adfeSMatthew G Knepley PetscFunctionBegin; 135257a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 135357a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 135457a9adfeSMatthew G Knepley PetscValidPointer(is,3); 135557a9adfeSMatthew G Knepley { 135657a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 135757a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 135857a9adfeSMatthew G Knepley PetscBool found; 135957a9adfeSMatthew G Knepley 136057a9adfeSMatthew G Knepley *is = PETSC_NULL; 136157a9adfeSMatthew G Knepley while(ilink) { 136257a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 136357a9adfeSMatthew G Knepley if (found) { 136457a9adfeSMatthew G Knepley *is = ilink->is; 136557a9adfeSMatthew G Knepley break; 136657a9adfeSMatthew G Knepley } 136757a9adfeSMatthew G Knepley ilink = ilink->next; 136857a9adfeSMatthew G Knepley } 136957a9adfeSMatthew G Knepley } 137057a9adfeSMatthew G Knepley PetscFunctionReturn(0); 137157a9adfeSMatthew G Knepley } 137257a9adfeSMatthew G Knepley 137357a9adfeSMatthew G Knepley #undef __FUNCT__ 137451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 137551f519a2SBarry Smith /*@ 137651f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 137751f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 137851f519a2SBarry Smith 1379ad4df100SBarry Smith Logically Collective on PC 138051f519a2SBarry Smith 138151f519a2SBarry Smith Input Parameters: 138251f519a2SBarry Smith + pc - the preconditioner context 138351f519a2SBarry Smith - bs - the block size 138451f519a2SBarry Smith 138551f519a2SBarry Smith Level: intermediate 138651f519a2SBarry Smith 138751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 138851f519a2SBarry Smith 138951f519a2SBarry Smith @*/ 13907087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 139151f519a2SBarry Smith { 13924ac538c5SBarry Smith PetscErrorCode ierr; 139351f519a2SBarry Smith 139451f519a2SBarry Smith PetscFunctionBegin; 13950700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1396c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13974ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 139851f519a2SBarry Smith PetscFunctionReturn(0); 139951f519a2SBarry Smith } 140051f519a2SBarry Smith 140151f519a2SBarry Smith #undef __FUNCT__ 140269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 14030971522cSBarry Smith /*@C 140469a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 14050971522cSBarry Smith 140669a612a9SBarry Smith Collective on KSP 14070971522cSBarry Smith 14080971522cSBarry Smith Input Parameter: 14090971522cSBarry Smith . pc - the preconditioner context 14100971522cSBarry Smith 14110971522cSBarry Smith Output Parameters: 141213e0d083SBarry Smith + n - the number of splits 141369a612a9SBarry Smith - pc - the array of KSP contexts 14140971522cSBarry Smith 14150971522cSBarry Smith Note: 1416d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1417d32f9abdSBarry Smith (not the KSP just the array that contains them). 14180971522cSBarry Smith 141969a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 14200971522cSBarry Smith 1421196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 1422196cc216SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the 1423196cc216SBarry Smith KSP array must be. 1424196cc216SBarry Smith 1425196cc216SBarry Smith 14260971522cSBarry Smith Level: advanced 14270971522cSBarry Smith 14280971522cSBarry Smith .seealso: PCFIELDSPLIT 14290971522cSBarry Smith @*/ 14307087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 14310971522cSBarry Smith { 14324ac538c5SBarry Smith PetscErrorCode ierr; 14330971522cSBarry Smith 14340971522cSBarry Smith PetscFunctionBegin; 14350700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 143613e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14374ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14380971522cSBarry Smith PetscFunctionReturn(0); 14390971522cSBarry Smith } 14400971522cSBarry Smith 1441e69d4d44SBarry Smith #undef __FUNCT__ 1442e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1443e69d4d44SBarry Smith /*@ 1444e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1445a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1446e69d4d44SBarry Smith 1447e69d4d44SBarry Smith Collective on PC 1448e69d4d44SBarry Smith 1449e69d4d44SBarry Smith Input Parameters: 1450e69d4d44SBarry Smith + pc - the preconditioner context 1451fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1452084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1453084e4875SJed Brown 1454e69d4d44SBarry Smith Options Database: 1455084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1456e69d4d44SBarry Smith 1457fd1303e9SJungho Lee Notes: 1458fd1303e9SJungho Lee $ If ptype is 1459fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1460fd1303e9SJungho Lee $ to this function). 1461fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1462fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1463fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1464fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1465fd1303e9SJungho Lee $ preconditioner 1466fd1303e9SJungho Lee 1467fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1468fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1469fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1470fd1303e9SJungho Lee 1471fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1472fd1303e9SJungho Lee the name since it sets a proceedure to use. 1473fd1303e9SJungho Lee 1474e69d4d44SBarry Smith Level: intermediate 1475e69d4d44SBarry Smith 1476fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1477e69d4d44SBarry Smith 1478e69d4d44SBarry Smith @*/ 14797087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1480e69d4d44SBarry Smith { 14814ac538c5SBarry Smith PetscErrorCode ierr; 1482e69d4d44SBarry Smith 1483e69d4d44SBarry Smith PetscFunctionBegin; 14840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14854ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1486e69d4d44SBarry Smith PetscFunctionReturn(0); 1487e69d4d44SBarry Smith } 1488e69d4d44SBarry Smith 1489e69d4d44SBarry Smith EXTERN_C_BEGIN 1490e69d4d44SBarry Smith #undef __FUNCT__ 1491e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14927087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1493e69d4d44SBarry Smith { 1494e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1495084e4875SJed Brown PetscErrorCode ierr; 1496e69d4d44SBarry Smith 1497e69d4d44SBarry Smith PetscFunctionBegin; 1498084e4875SJed Brown jac->schurpre = ptype; 1499084e4875SJed Brown if (pre) { 15006bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1501084e4875SJed Brown jac->schur_user = pre; 1502084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1503084e4875SJed Brown } 1504e69d4d44SBarry Smith PetscFunctionReturn(0); 1505e69d4d44SBarry Smith } 1506e69d4d44SBarry Smith EXTERN_C_END 1507e69d4d44SBarry Smith 150830ad9308SMatthew Knepley #undef __FUNCT__ 1509c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1510ab1df9f5SJed Brown /*@ 1511c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1512ab1df9f5SJed Brown 1513ab1df9f5SJed Brown Collective on PC 1514ab1df9f5SJed Brown 1515ab1df9f5SJed Brown Input Parameters: 1516ab1df9f5SJed Brown + pc - the preconditioner context 1517c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1518ab1df9f5SJed Brown 1519ab1df9f5SJed Brown Options Database: 1520c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1521ab1df9f5SJed Brown 1522ab1df9f5SJed Brown 1523ab1df9f5SJed Brown Level: intermediate 1524ab1df9f5SJed Brown 1525ab1df9f5SJed Brown Notes: 1526ab1df9f5SJed Brown The FULL factorization is 1527ab1df9f5SJed Brown 1528ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1529ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1530ab1df9f5SJed Brown 15316be4592eSBarry 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, 15326be4592eSBarry 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). 1533ab1df9f5SJed Brown 15346be4592eSBarry 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 15356be4592eSBarry 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 15366be4592eSBarry 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, 15376be4592eSBarry 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. 1538ab1df9f5SJed Brown 15396be4592eSBarry 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 15406be4592eSBarry 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). 1541ab1df9f5SJed Brown 1542ab1df9f5SJed Brown References: 1543ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1544ab1df9f5SJed Brown 1545ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1546ab1df9f5SJed Brown 1547ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1548ab1df9f5SJed Brown @*/ 1549c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1550ab1df9f5SJed Brown { 1551ab1df9f5SJed Brown PetscErrorCode ierr; 1552ab1df9f5SJed Brown 1553ab1df9f5SJed Brown PetscFunctionBegin; 1554ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1555c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1556ab1df9f5SJed Brown PetscFunctionReturn(0); 1557ab1df9f5SJed Brown } 1558ab1df9f5SJed Brown 1559ab1df9f5SJed Brown #undef __FUNCT__ 1560c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1561c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1562ab1df9f5SJed Brown { 1563ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1564ab1df9f5SJed Brown 1565ab1df9f5SJed Brown PetscFunctionBegin; 1566ab1df9f5SJed Brown jac->schurfactorization = ftype; 1567ab1df9f5SJed Brown PetscFunctionReturn(0); 1568ab1df9f5SJed Brown } 1569ab1df9f5SJed Brown 1570ab1df9f5SJed Brown #undef __FUNCT__ 157130ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 157230ad9308SMatthew Knepley /*@C 15738c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 157430ad9308SMatthew Knepley 157530ad9308SMatthew Knepley Collective on KSP 157630ad9308SMatthew Knepley 157730ad9308SMatthew Knepley Input Parameter: 157830ad9308SMatthew Knepley . pc - the preconditioner context 157930ad9308SMatthew Knepley 158030ad9308SMatthew Knepley Output Parameters: 1581a04f6461SBarry Smith + A00 - the (0,0) block 1582a04f6461SBarry Smith . A01 - the (0,1) block 1583a04f6461SBarry Smith . A10 - the (1,0) block 1584a04f6461SBarry Smith - A11 - the (1,1) block 158530ad9308SMatthew Knepley 158630ad9308SMatthew Knepley Level: advanced 158730ad9308SMatthew Knepley 158830ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 158930ad9308SMatthew Knepley @*/ 1590a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 159130ad9308SMatthew Knepley { 159230ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 159330ad9308SMatthew Knepley 159430ad9308SMatthew Knepley PetscFunctionBegin; 15950700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1596c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1597a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1598a04f6461SBarry Smith if (A01) *A01 = jac->B; 1599a04f6461SBarry Smith if (A10) *A10 = jac->C; 1600a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 160130ad9308SMatthew Knepley PetscFunctionReturn(0); 160230ad9308SMatthew Knepley } 160330ad9308SMatthew Knepley 160479416396SBarry Smith EXTERN_C_BEGIN 160579416396SBarry Smith #undef __FUNCT__ 160679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 16077087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 160879416396SBarry Smith { 160979416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1610e69d4d44SBarry Smith PetscErrorCode ierr; 161179416396SBarry Smith 161279416396SBarry Smith PetscFunctionBegin; 161379416396SBarry Smith jac->type = type; 16143b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 16153b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 16163b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1617e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1618e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1619c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1620e69d4d44SBarry Smith 16213b224e63SBarry Smith } else { 16223b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 16233b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1624e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 16259e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1626c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 16273b224e63SBarry Smith } 162879416396SBarry Smith PetscFunctionReturn(0); 162979416396SBarry Smith } 163079416396SBarry Smith EXTERN_C_END 163179416396SBarry Smith 163251f519a2SBarry Smith EXTERN_C_BEGIN 163351f519a2SBarry Smith #undef __FUNCT__ 163451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16357087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 163651f519a2SBarry Smith { 163751f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 163851f519a2SBarry Smith 163951f519a2SBarry Smith PetscFunctionBegin; 164065e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 164165e19b50SBarry 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); 164251f519a2SBarry Smith jac->bs = bs; 164351f519a2SBarry Smith PetscFunctionReturn(0); 164451f519a2SBarry Smith } 164551f519a2SBarry Smith EXTERN_C_END 164651f519a2SBarry Smith 164779416396SBarry Smith #undef __FUNCT__ 164879416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1649bc08b0f1SBarry Smith /*@ 165079416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 165179416396SBarry Smith 165279416396SBarry Smith Collective on PC 165379416396SBarry Smith 165479416396SBarry Smith Input Parameter: 165579416396SBarry Smith . pc - the preconditioner context 165681540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 165779416396SBarry Smith 165879416396SBarry Smith Options Database Key: 1659a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 166079416396SBarry Smith 1661b02e2d75SMatthew G Knepley Level: Intermediate 166279416396SBarry Smith 166379416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 166479416396SBarry Smith 166579416396SBarry Smith .seealso: PCCompositeSetType() 166679416396SBarry Smith 166779416396SBarry Smith @*/ 16687087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 166979416396SBarry Smith { 16704ac538c5SBarry Smith PetscErrorCode ierr; 167179416396SBarry Smith 167279416396SBarry Smith PetscFunctionBegin; 16730700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16744ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 167579416396SBarry Smith PetscFunctionReturn(0); 167679416396SBarry Smith } 167779416396SBarry Smith 1678b02e2d75SMatthew G Knepley #undef __FUNCT__ 1679b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1680b02e2d75SMatthew G Knepley /*@ 1681b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1682b02e2d75SMatthew G Knepley 1683b02e2d75SMatthew G Knepley Not collective 1684b02e2d75SMatthew G Knepley 1685b02e2d75SMatthew G Knepley Input Parameter: 1686b02e2d75SMatthew G Knepley . pc - the preconditioner context 1687b02e2d75SMatthew G Knepley 1688b02e2d75SMatthew G Knepley Output Parameter: 1689b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1690b02e2d75SMatthew G Knepley 1691b02e2d75SMatthew G Knepley Level: Intermediate 1692b02e2d75SMatthew G Knepley 1693b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1694b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1695b02e2d75SMatthew G Knepley @*/ 1696b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1697b02e2d75SMatthew G Knepley { 1698b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1699b02e2d75SMatthew G Knepley 1700b02e2d75SMatthew G Knepley PetscFunctionBegin; 1701b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1702b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1703b02e2d75SMatthew G Knepley *type = jac->type; 1704b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1705b02e2d75SMatthew G Knepley } 1706b02e2d75SMatthew G Knepley 17070971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 17080971522cSBarry Smith /*MC 1709a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1710a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 17110971522cSBarry Smith 1712edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1713edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 17140971522cSBarry Smith 1715a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 171669a612a9SBarry Smith and set the options directly on the resulting KSP object 17170971522cSBarry Smith 17180971522cSBarry Smith Level: intermediate 17190971522cSBarry Smith 172079416396SBarry Smith Options Database Keys: 172181540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 172281540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 172381540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 172481540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 17250f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 17260f188ba9SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag 1727435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 172879416396SBarry Smith 17295d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 17305d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 17315d4c12cdSJungho Lee 1732c8a0d604SMatthew G Knepley Notes: 1733c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1734d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1735d32f9abdSBarry Smith 1736d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1737d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1738d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1739d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1740d32f9abdSBarry Smith 1741c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1742c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1743c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1744c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1745c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1746a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1747c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1748c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1749c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1750a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1751c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1752c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 17535668aaf4SBarry Smith diag gives 1754c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1755c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 17565668aaf4SBarry 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 1757c8a0d604SMatthew G Knepley $ ( A00 0 ) 1758c8a0d604SMatthew G Knepley $ ( A10 S ) 1759c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1760c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1761c8a0d604SMatthew G Knepley $ ( 0 S ) 1762c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1763e69d4d44SBarry Smith 1764edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1765edf189efSBarry Smith is used automatically for a second block. 1766edf189efSBarry Smith 1767ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1768ff218e97SBarry Smith Generally it should be used with the AIJ format. 1769ff218e97SBarry Smith 1770ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1771ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1772ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 17730716a85fSBarry Smith 1774a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 17750971522cSBarry Smith 17767e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1777e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 17780971522cSBarry Smith M*/ 17790971522cSBarry Smith 17800971522cSBarry Smith EXTERN_C_BEGIN 17810971522cSBarry Smith #undef __FUNCT__ 17820971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 17837087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 17840971522cSBarry Smith { 17850971522cSBarry Smith PetscErrorCode ierr; 17860971522cSBarry Smith PC_FieldSplit *jac; 17870971522cSBarry Smith 17880971522cSBarry Smith PetscFunctionBegin; 178938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 17900971522cSBarry Smith jac->bs = -1; 17910971522cSBarry Smith jac->nsplits = 0; 17923e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1793e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1794c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 179551f519a2SBarry Smith 17960971522cSBarry Smith pc->data = (void*)jac; 17970971522cSBarry Smith 17980971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1799421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 18000971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1801574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 18020971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 18030971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 18040971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 18050971522cSBarry Smith pc->ops->applyrichardson = 0; 18060971522cSBarry Smith 180769a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 180869a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 18090971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 18100971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1811704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1812704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 181379416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 181479416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 181551f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 181651f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 18170971522cSBarry Smith PetscFunctionReturn(0); 18180971522cSBarry Smith } 18190971522cSBarry Smith EXTERN_C_END 18200971522cSBarry Smith 18210971522cSBarry Smith 1822a541d17aSBarry Smith 1823