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; 80*d9884438SBarry 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); 86*d9884438SBarry 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); 117*d9884438SBarry Smith } if (isdraw) { 118*d9884438SBarry Smith PetscDraw draw; 119*d9884438SBarry Smith PetscReal x,y,w,wd; 120*d9884438SBarry Smith 121*d9884438SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 122*d9884438SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 123*d9884438SBarry Smith w = 2*PetscMin(1.0 - x,x); 124*d9884438SBarry Smith wd = w/(jac->nsplits + 1); 125*d9884438SBarry Smith x = x - wd*(jac->nsplits-1)/2.0; 126*d9884438SBarry Smith for (i=0; i<jac->nsplits; i++) { 127*d9884438SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr); 128*d9884438SBarry Smith ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 129*d9884438SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 130*d9884438SBarry Smith x += wd; 131*d9884438SBarry Smith ilink = ilink->next; 132*d9884438SBarry 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; 143ace3abfcSBarry Smith PetscBool iascii; 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); 1493b224e63SBarry Smith if (iascii) { 1502c9966d7SBarry Smith if (jac->bs > 0) { 151c9c6ffaaSJed Brown ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1522c9966d7SBarry Smith } else { 1532c9966d7SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr); 1542c9966d7SBarry Smith } 155a3df900dSMatthew G Knepley if (jac->realdiagonal) { 156a3df900dSMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr); 157a3df900dSMatthew G Knepley } 1583e8b8b31SMatthew G Knepley switch(jac->schurpre) { 1593e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_SELF: 1603e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break; 161b02e2d75SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_DIAG: 1623e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr);break; 1633e8b8b31SMatthew G Knepley case PC_FIELDSPLIT_SCHUR_PRE_USER: 1643e8b8b31SMatthew G Knepley if (jac->schur_user) { 1653e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr); 1663e8b8b31SMatthew G Knepley } else { 1673e8b8b31SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the block diagonal part of A11\n");CHKERRQ(ierr); 1683e8b8b31SMatthew G Knepley } 1693e8b8b31SMatthew G Knepley break; 1703e8b8b31SMatthew G Knepley default: 1713e8b8b31SMatthew G Knepley SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre); 1723e8b8b31SMatthew G Knepley } 1733b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1743b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1753b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1763b224e63SBarry Smith if (ilink->fields) { 1773b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1783b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1793b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1803b224e63SBarry Smith if (j > 0) { 1813b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1823b224e63SBarry Smith } 1833b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1843b224e63SBarry Smith } 1853b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1863b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1873b224e63SBarry Smith } else { 1883b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1893b224e63SBarry Smith } 1903b224e63SBarry Smith ilink = ilink->next; 1913b224e63SBarry Smith } 192435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr); 1933b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 194443836d0SMatthew G Knepley ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr); 1953b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 196443836d0SMatthew G Knepley if (jac->kspupper != jac->head->ksp) { 197443836d0SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr); 198443836d0SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 199443836d0SMatthew G Knepley ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr); 200443836d0SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 201443836d0SMatthew G Knepley } 202435f959eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr); 2033b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 20412cae6f2SJed Brown if (jac->kspschur) { 2053b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 20612cae6f2SJed Brown } else { 20712cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 20812cae6f2SJed Brown } 2093b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2103b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2113b224e63SBarry Smith } else { 21265e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 2133b224e63SBarry Smith } 2143b224e63SBarry Smith PetscFunctionReturn(0); 2153b224e63SBarry Smith } 2163b224e63SBarry Smith 2173b224e63SBarry Smith #undef __FUNCT__ 2186c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 2196c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 2206c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 2216c924f48SJed Brown { 2226c924f48SJed Brown PetscErrorCode ierr; 2236c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2245d4c12cdSJungho Lee PetscInt i,nfields,*ifields,nfields_col,*ifields_col; 2255d4c12cdSJungho Lee PetscBool flg,flg_col; 2265d4c12cdSJungho Lee char optionname[128],splitname[8],optionname_col[128]; 2276c924f48SJed Brown 2286c924f48SJed Brown PetscFunctionBegin; 2296c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 2305d4c12cdSJungho Lee ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr); 2316c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 2328caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 2338caf3d72SBarry Smith ierr = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 2348caf3d72SBarry Smith ierr = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr); 2356c924f48SJed Brown nfields = jac->bs; 23629499fbbSJungho Lee nfields_col = jac->bs; 2376c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 2385d4c12cdSJungho Lee ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr); 2396c924f48SJed Brown if (!flg) break; 2405d4c12cdSJungho Lee else if (flg && !flg_col) { 2415d4c12cdSJungho Lee if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2425d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr); 2435d4c12cdSJungho Lee } 2445d4c12cdSJungho Lee else { 2455d4c12cdSJungho Lee if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 2465d4c12cdSJungho Lee if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match"); 2475d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr); 2485d4c12cdSJungho Lee } 2496c924f48SJed Brown } 2506c924f48SJed Brown if (i > 0) { 2516c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 2526c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 2536c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 2546c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 2556c924f48SJed Brown } 2566c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 2575d4c12cdSJungho Lee ierr = PetscFree(ifields_col);CHKERRQ(ierr); 2586c924f48SJed Brown PetscFunctionReturn(0); 2596c924f48SJed Brown } 2606c924f48SJed Brown 2616c924f48SJed Brown #undef __FUNCT__ 26269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 26369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 2640971522cSBarry Smith { 2650971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2660971522cSBarry Smith PetscErrorCode ierr; 2675a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 2687287d2a3SDmitry Karpeev PetscBool fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE; 2696c924f48SJed Brown PetscInt i; 2700971522cSBarry Smith 2710971522cSBarry Smith PetscFunctionBegin; 2727287d2a3SDmitry Karpeev /* 2737287d2a3SDmitry Karpeev Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset. 2747287d2a3SDmitry Karpeev Should probably be rewritten. 2757287d2a3SDmitry Karpeev */ 2767287d2a3SDmitry Karpeev if (!ilink || jac->reset) { 277d0af7cd3SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 278d0af7cd3SBarry Smith if (pc->dm && !stokes) { 279bafc1b83SMatthew G Knepley PetscInt numFields, f, i, j; 2800784a22cSJed Brown char **fieldNames; 2817b62db95SJungho Lee IS *fields; 282e7c4fc90SDmitry Karpeev DM *dms; 283bafc1b83SMatthew G Knepley DM subdm[128]; 284bafc1b83SMatthew G Knepley PetscBool flg; 285bafc1b83SMatthew G Knepley 28616621825SDmitry Karpeev ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr); 287bafc1b83SMatthew G Knepley /* Allow the user to prescribe the splits */ 288bafc1b83SMatthew G Knepley for (i = 0, flg = PETSC_TRUE; ; i++) { 289bafc1b83SMatthew G Knepley PetscInt ifields[128]; 290bafc1b83SMatthew G Knepley IS compField; 291bafc1b83SMatthew G Knepley char optionname[128], splitname[8]; 292bafc1b83SMatthew G Knepley PetscInt nfields = numFields; 293bafc1b83SMatthew G Knepley 2948caf3d72SBarry Smith ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr); 295bafc1b83SMatthew G Knepley ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr); 296bafc1b83SMatthew G Knepley if (!flg) break; 297bafc1b83SMatthew G Knepley if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields); 298bafc1b83SMatthew G Knepley ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr); 299bafc1b83SMatthew G Knepley if (nfields == 1) { 300bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr); 3014d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr); 3024d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 303bafc1b83SMatthew G Knepley } else { 3048caf3d72SBarry Smith ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr); 305bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr); 3064d2389d7SMatthew G Knepley /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr); 3074d2389d7SMatthew G Knepley ierr = ISView(compField, PETSC_NULL);CHKERRQ(ierr); */ 3087287d2a3SDmitry Karpeev } 309bafc1b83SMatthew G Knepley ierr = ISDestroy(&compField);CHKERRQ(ierr); 310bafc1b83SMatthew G Knepley for (j = 0; j < nfields; ++j) { 311bafc1b83SMatthew G Knepley f = ifields[j]; 3127b62db95SJungho Lee ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 3137b62db95SJungho Lee ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 3147b62db95SJungho Lee } 315bafc1b83SMatthew G Knepley } 316bafc1b83SMatthew G Knepley if (i == 0) { 317bafc1b83SMatthew G Knepley for (f = 0; f < numFields; ++f) { 318bafc1b83SMatthew G Knepley ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr); 319bafc1b83SMatthew G Knepley ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr); 320bafc1b83SMatthew G Knepley ierr = ISDestroy(&fields[f]);CHKERRQ(ierr); 321bafc1b83SMatthew G Knepley } 322bafc1b83SMatthew G Knepley } else { 323bafc1b83SMatthew G Knepley ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr); 324bafc1b83SMatthew G Knepley for (j = 0; j < i; ++j) { 325bafc1b83SMatthew G Knepley dms[j] = subdm[j]; 326bafc1b83SMatthew G Knepley } 327bafc1b83SMatthew G Knepley } 3287b62db95SJungho Lee ierr = PetscFree(fieldNames);CHKERRQ(ierr); 3297b62db95SJungho Lee ierr = PetscFree(fields);CHKERRQ(ierr); 330e7c4fc90SDmitry Karpeev if (dms) { 3318b8307b2SJed Brown ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr); 332bafc1b83SMatthew G Knepley for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) { 3337287d2a3SDmitry Karpeev const char *prefix; 3347287d2a3SDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix); CHKERRQ(ierr); 3357287d2a3SDmitry Karpeev ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix); CHKERRQ(ierr); 3367b62db95SJungho Lee ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr); 3377b62db95SJungho Lee ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr); 3387287d2a3SDmitry Karpeev ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0); CHKERRQ(ierr); 339e7c4fc90SDmitry Karpeev ierr = DMDestroy(&dms[i]); CHKERRQ(ierr); 3402fa5ba8aSJed Brown } 3417b62db95SJungho Lee ierr = PetscFree(dms);CHKERRQ(ierr); 3428b8307b2SJed Brown } 34366ffff09SJed Brown } else { 344521d7252SBarry Smith if (jac->bs <= 0) { 345704ba839SBarry Smith if (pc->pmat) { 346521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 347704ba839SBarry Smith } else { 348704ba839SBarry Smith jac->bs = 1; 349704ba839SBarry Smith } 350521d7252SBarry Smith } 351d32f9abdSBarry Smith 3527287d2a3SDmitry Karpeev ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,PETSC_NULL);CHKERRQ(ierr); 3536ce1633cSBarry Smith if (stokes) { 3546ce1633cSBarry Smith IS zerodiags,rest; 3556ce1633cSBarry Smith PetscInt nmin,nmax; 3566ce1633cSBarry Smith 3576ce1633cSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 3586ce1633cSBarry Smith ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr); 3596ce1633cSBarry Smith ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr); 3607287d2a3SDmitry Karpeev if (jac->reset) { 3617287d2a3SDmitry Karpeev jac->head->is = rest; 3627287d2a3SDmitry Karpeev jac->head->next->is = zerodiags; 3637287d2a3SDmitry Karpeev } 3647287d2a3SDmitry Karpeev else { 3656ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr); 3666ce1633cSBarry Smith ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr); 3677287d2a3SDmitry Karpeev } 368fcfd50ebSBarry Smith ierr = ISDestroy(&zerodiags);CHKERRQ(ierr); 369fcfd50ebSBarry Smith ierr = ISDestroy(&rest);CHKERRQ(ierr); 3706ce1633cSBarry Smith } else { 3717287d2a3SDmitry Karpeev if (jac->reset) 3727287d2a3SDmitry Karpeev SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used"); 3737287d2a3SDmitry Karpeev if (!fieldsplit_default) { 374d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 375d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 3766c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 3776c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 378d32f9abdSBarry Smith } 3797287d2a3SDmitry Karpeev if (fieldsplit_default || !jac->splitdefined) { 380d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 381db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 3826c924f48SJed Brown char splitname[8]; 3838caf3d72SBarry Smith ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr); 3845d4c12cdSJungho Lee ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr); 38579416396SBarry Smith } 3865d4c12cdSJungho Lee jac->defaultsplit = PETSC_TRUE; 387521d7252SBarry Smith } 38866ffff09SJed Brown } 3896ce1633cSBarry Smith } 390edf189efSBarry Smith } else if (jac->nsplits == 1) { 391edf189efSBarry Smith if (ilink->is) { 392edf189efSBarry Smith IS is2; 393edf189efSBarry Smith PetscInt nmin,nmax; 394edf189efSBarry Smith 395edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 396edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 397db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 398fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 3997b62db95SJungho Lee } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 400edf189efSBarry Smith } 401d0af7cd3SBarry Smith 402d0af7cd3SBarry Smith 4037b62db95SJungho Lee if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits); 40469a612a9SBarry Smith PetscFunctionReturn(0); 40569a612a9SBarry Smith } 40669a612a9SBarry Smith 407514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg); 408514bf10dSMatthew G Knepley 40969a612a9SBarry Smith #undef __FUNCT__ 41069a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 41169a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 41269a612a9SBarry Smith { 41369a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 41469a612a9SBarry Smith PetscErrorCode ierr; 4155a9f2f41SSatish Balay PC_FieldSplitLink ilink; 4162c9966d7SBarry Smith PetscInt i,nsplit; 41769a612a9SBarry Smith MatStructure flag = pc->flag; 4185f4ab4e1SJungho Lee PetscBool sorted, sorted_col; 41969a612a9SBarry Smith 42069a612a9SBarry Smith PetscFunctionBegin; 42169a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 42297bbdb24SBarry Smith nsplit = jac->nsplits; 4235a9f2f41SSatish Balay ilink = jac->head; 42497bbdb24SBarry Smith 42597bbdb24SBarry Smith /* get the matrices for each split */ 426704ba839SBarry Smith if (!jac->issetup) { 4271b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 42897bbdb24SBarry Smith 429704ba839SBarry Smith jac->issetup = PETSC_TRUE; 430704ba839SBarry Smith 4315d4c12cdSJungho Lee /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 4322c9966d7SBarry Smith if (jac->defaultsplit || !ilink->is) { 4332c9966d7SBarry Smith if (jac->bs <= 0) jac->bs = nsplit; 4342c9966d7SBarry Smith } 43551f519a2SBarry Smith bs = jac->bs; 43697bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 4371b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 4381b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 4391b9fc7fcSBarry Smith if (jac->defaultsplit) { 440704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 4415f4ab4e1SJungho Lee ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr); 442704ba839SBarry Smith } else if (!ilink->is) { 443ccb205f8SBarry Smith if (ilink->nfields > 1) { 4445f4ab4e1SJungho Lee PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col; 4455a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 4465f4ab4e1SJungho Lee ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr); 4471b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 4481b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 4491b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 4505f4ab4e1SJungho Lee jj[nfields*j + k] = rstart + bs*j + fields_col[k]; 45197bbdb24SBarry Smith } 45297bbdb24SBarry Smith } 453d67e408aSBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr); 4545f4ab4e1SJungho Lee ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr); 455ccb205f8SBarry Smith } else { 456704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 4575f4ab4e1SJungho Lee ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr); 458ccb205f8SBarry Smith } 4593e197d65SBarry Smith } 460704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 4615f4ab4e1SJungho Lee if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); } 4625f4ab4e1SJungho Lee if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 463704ba839SBarry Smith ilink = ilink->next; 4641b9fc7fcSBarry Smith } 4651b9fc7fcSBarry Smith } 4661b9fc7fcSBarry Smith 467704ba839SBarry Smith ilink = jac->head; 46897bbdb24SBarry Smith if (!jac->pmat) { 469bdddcaaaSMatthew G Knepley Vec xtmp; 470bdddcaaaSMatthew G Knepley 471bdddcaaaSMatthew G Knepley ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 472cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 473bdddcaaaSMatthew G Knepley ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 474cf502942SBarry Smith for (i=0; i<nsplit; i++) { 475bdddcaaaSMatthew G Knepley MatNullSpace sp; 476bdddcaaaSMatthew G Knepley 477a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 478a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &jac->pmat[i]);CHKERRQ(ierr); 479a3df900dSMatthew G Knepley if (jac->pmat[i]) { 480a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr); 481a3df900dSMatthew G Knepley if (jac->type == PC_COMPOSITE_SCHUR) { 482a3df900dSMatthew G Knepley jac->schur_user = jac->pmat[i]; 483a3df900dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr); 484a3df900dSMatthew G Knepley } 485a3df900dSMatthew G Knepley } else { 4865f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 487a3df900dSMatthew G Knepley } 488bdddcaaaSMatthew G Knepley /* create work vectors for each split */ 489bdddcaaaSMatthew G Knepley ierr = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr); 490443836d0SMatthew G Knepley ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = PETSC_NULL; 491bdddcaaaSMatthew G Knepley /* compute scatter contexts needed by multiplicative versions and non-default splits */ 492bdddcaaaSMatthew G Knepley ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 493ed1f0337SMatthew G Knepley /* Check for null space attached to IS */ 494bafc1b83SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject *) &sp);CHKERRQ(ierr); 495bafc1b83SMatthew G Knepley if (sp) { 496bafc1b83SMatthew G Knepley ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 497bafc1b83SMatthew G Knepley } 498ed1f0337SMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject *) &sp);CHKERRQ(ierr); 499ed1f0337SMatthew G Knepley if (sp) { 500ed1f0337SMatthew G Knepley ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr); 501ed1f0337SMatthew G Knepley } 502704ba839SBarry Smith ilink = ilink->next; 503cf502942SBarry Smith } 504bdddcaaaSMatthew G Knepley ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 50597bbdb24SBarry Smith } else { 506cf502942SBarry Smith for (i=0; i<nsplit; i++) { 507a3df900dSMatthew G Knepley Mat pmat; 508a3df900dSMatthew G Knepley 509a3df900dSMatthew G Knepley /* Check for preconditioning matrix attached to IS */ 510a3df900dSMatthew G Knepley ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject *) &pmat);CHKERRQ(ierr); 511a3df900dSMatthew G Knepley if (!pmat) { 5125f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 513a3df900dSMatthew G Knepley } 514704ba839SBarry Smith ilink = ilink->next; 515cf502942SBarry Smith } 51697bbdb24SBarry Smith } 517519d70e2SJed Brown if (jac->realdiagonal) { 518519d70e2SJed Brown ilink = jac->head; 519519d70e2SJed Brown if (!jac->mat) { 520519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 521519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5225f4ab4e1SJungho Lee ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 523519d70e2SJed Brown ilink = ilink->next; 524519d70e2SJed Brown } 525519d70e2SJed Brown } else { 526519d70e2SJed Brown for (i=0; i<nsplit; i++) { 5275f4ab4e1SJungho Lee if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);} 528519d70e2SJed Brown ilink = ilink->next; 529519d70e2SJed Brown } 530519d70e2SJed Brown } 531519d70e2SJed Brown } else { 532519d70e2SJed Brown jac->mat = jac->pmat; 533519d70e2SJed Brown } 53497bbdb24SBarry Smith 5356c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 53668dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 53768dd23aaSBarry Smith ilink = jac->head; 53868dd23aaSBarry Smith if (!jac->Afield) { 53968dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 54068dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5414aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 54268dd23aaSBarry Smith ilink = ilink->next; 54368dd23aaSBarry Smith } 54468dd23aaSBarry Smith } else { 54568dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 5464aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 54768dd23aaSBarry Smith ilink = ilink->next; 54868dd23aaSBarry Smith } 54968dd23aaSBarry Smith } 55068dd23aaSBarry Smith } 55168dd23aaSBarry Smith 5523b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 5533b224e63SBarry Smith IS ccis; 5544aa3045dSJed Brown PetscInt rstart,rend; 555093c86ffSJed Brown char lscname[256]; 556093c86ffSJed Brown PetscObject LSC_L; 557e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 55868dd23aaSBarry Smith 559e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 560e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 561e6cab6aaSJed Brown 5623b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 5633b224e63SBarry Smith if (jac->schur) { 564a3314f2cSMatthew G Knepley KSP kspA = jac->head->ksp, kspInner = PETSC_NULL, kspUpper = jac->kspupper; 565443836d0SMatthew G Knepley 566fb3147dbSMatthew G Knepley ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr); 5673b224e63SBarry Smith ilink = jac->head; 56849bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5694aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 570fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5713b224e63SBarry Smith ilink = ilink->next; 57249bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5734aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 574fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 575a3df900dSMatthew G Knepley ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr); 576443836d0SMatthew G Knepley if (kspA != kspInner) { 577443836d0SMatthew G Knepley ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 578443836d0SMatthew G Knepley } 579443836d0SMatthew G Knepley if (kspUpper != kspA) { 580443836d0SMatthew G Knepley ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr); 581443836d0SMatthew G Knepley } 582084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 5833b224e63SBarry Smith } else { 5841cee3971SBarry Smith KSP ksp; 585bafc1b83SMatthew G Knepley const char *Dprefix; 586bafc1b83SMatthew G Knepley char schurprefix[256]; 587514bf10dSMatthew G Knepley char schurtestoption[256]; 588bdddcaaaSMatthew G Knepley MatNullSpace sp; 589514bf10dSMatthew G Knepley PetscBool flg; 5903b224e63SBarry Smith 591a04f6461SBarry Smith /* extract the A01 and A10 matrices */ 5923b224e63SBarry Smith ilink = jac->head; 59349bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5944aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 595fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 5963b224e63SBarry Smith ilink = ilink->next; 59749bb4cd7SJungho Lee ierr = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr); 5984aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 599fcfd50ebSBarry Smith ierr = ISDestroy(&ccis);CHKERRQ(ierr); 60020252d06SBarry Smith 60120252d06SBarry Smith /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */ 60220252d06SBarry Smith ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr); 60320252d06SBarry Smith ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr); 60420252d06SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr); 60520252d06SBarry Smith ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 60620252d06SBarry Smith /* Indent this deeper to emphasize the "inner" nature of this solver. */ 60720252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr); 60820252d06SBarry Smith ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr); 60920252d06SBarry Smith ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr); 61020252d06SBarry Smith 611bdddcaaaSMatthew G Knepley ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr); 61220252d06SBarry Smith if (sp) { 61320252d06SBarry Smith ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr); 61420252d06SBarry Smith } 61520252d06SBarry Smith 61620252d06SBarry Smith ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr); 617514bf10dSMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 618514bf10dSMatthew G Knepley if (flg) { 619514bf10dSMatthew G Knepley DM dmInner; 620514bf10dSMatthew G Knepley 621514bf10dSMatthew G Knepley /* Set DM for new solver */ 622bafc1b83SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 623bafc1b83SMatthew G Knepley ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr); 6247287d2a3SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 625443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 626443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr); 627514bf10dSMatthew G Knepley } else { 628514bf10dSMatthew G Knepley ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr); 629514bf10dSMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 630bafc1b83SMatthew G Knepley } 63120b26d62SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 6323b224e63SBarry Smith 633443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr); 634443836d0SMatthew G Knepley ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, PETSC_NULL, &flg);CHKERRQ(ierr); 635443836d0SMatthew G Knepley if (flg) { 636443836d0SMatthew G Knepley DM dmInner; 637443836d0SMatthew G Knepley 638443836d0SMatthew G Knepley ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr); 639443836d0SMatthew G Knepley ierr = KSPCreate(((PetscObject) pc)->comm, &jac->kspupper);CHKERRQ(ierr); 640443836d0SMatthew G Knepley ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr); 641443836d0SMatthew G Knepley ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr); 642443836d0SMatthew G Knepley ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr); 643443836d0SMatthew G Knepley ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr); 644443836d0SMatthew G Knepley ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr); 645443836d0SMatthew G Knepley ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr); 646443836d0SMatthew G Knepley ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr); 647443836d0SMatthew G Knepley } else { 648443836d0SMatthew G Knepley jac->kspupper = jac->head->ksp; 649443836d0SMatthew G Knepley ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr); 650443836d0SMatthew G Knepley } 651443836d0SMatthew G Knepley 6523b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur); CHKERRQ(ierr); 6539005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr); 65420252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1); CHKERRQ(ierr); 655084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 656084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 6577233a360SDmitry Karpeev PC pcschur; 6587233a360SDmitry Karpeev ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr); 6597233a360SDmitry Karpeev ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr); 660084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 661e69d4d44SBarry Smith } 6627287d2a3SDmitry Karpeev ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix); CHKERRQ(ierr); 6637287d2a3SDmitry Karpeev ierr = KSPSetOptionsPrefix(jac->kspschur, Dprefix); CHKERRQ(ierr); 6643b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 66520b26d62SBarry Smith /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */ 66620b26d62SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 6673b224e63SBarry Smith } 668093c86ffSJed Brown 669093c86ffSJed Brown /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */ 6708caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr); 671093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 672093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 673093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);} 6748caf3d72SBarry Smith ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr); 675093c86ffSJed Brown ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr); 676093c86ffSJed Brown if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);} 677093c86ffSJed Brown if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);} 6783b224e63SBarry Smith } else { 67968bd789dSDmitry Karpeev /* set up the individual splits' PCs */ 68097bbdb24SBarry Smith i = 0; 6815a9f2f41SSatish Balay ilink = jac->head; 6825a9f2f41SSatish Balay while (ilink) { 683519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 6843b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 685c1570756SJed Brown if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);} 68697bbdb24SBarry Smith i++; 6875a9f2f41SSatish Balay ilink = ilink->next; 6880971522cSBarry Smith } 6893b224e63SBarry Smith } 6903b224e63SBarry Smith 691c1570756SJed Brown jac->suboptionsset = PETSC_TRUE; 6920971522cSBarry Smith PetscFunctionReturn(0); 6930971522cSBarry Smith } 6940971522cSBarry Smith 6955a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 696ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 697ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 6985a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 699ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 700ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 70179416396SBarry Smith 7020971522cSBarry Smith #undef __FUNCT__ 7033b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 7043b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 7053b224e63SBarry Smith { 7063b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7073b224e63SBarry Smith PetscErrorCode ierr; 7083b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 709443836d0SMatthew G Knepley KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper; 7103b224e63SBarry Smith 7113b224e63SBarry Smith PetscFunctionBegin; 712c5d2311dSJed Brown switch (jac->schurfactorization) { 713c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_DIAG: 714a04f6461SBarry Smith /* [A00 0; 0 -S], positive definite, suitable for MINRES */ 715c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 716c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 717c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 718443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 719c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 720c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 721c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 722c5d2311dSJed Brown ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr); 723c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 724c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 725c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 726c5d2311dSJed Brown break; 727c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_LOWER: 728a04f6461SBarry Smith /* [A00 0; A10 S], suitable for left preconditioning */ 729c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 730c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 731443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 732c5d2311dSJed Brown ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 733c5d2311dSJed Brown ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr); 734c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 735c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 736c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 737c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 738c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 739c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 740c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 741c5d2311dSJed Brown break; 742c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_UPPER: 743a04f6461SBarry Smith /* [A00 A01; 0 S], suitable for right preconditioning */ 744c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 745c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 746c5d2311dSJed Brown ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 747c5d2311dSJed Brown ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 748c5d2311dSJed Brown ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr); 749c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 750c5d2311dSJed Brown ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 751c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 752443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 753c5d2311dSJed Brown ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 754c5d2311dSJed Brown ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 755c5d2311dSJed Brown ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 756c5d2311dSJed Brown break; 757c9c6ffaaSJed Brown case PC_FIELDSPLIT_SCHUR_FACT_FULL: 758a04f6461SBarry 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 */ 7593b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7603b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 761443836d0SMatthew G Knepley ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 7623b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 7633b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 7643b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7653b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7663b224e63SBarry Smith 7673b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 7683b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7693b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7703b224e63SBarry Smith 771443836d0SMatthew G Knepley if (kspUpper == kspA) { 7723b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 7733b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 774443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 775443836d0SMatthew G Knepley } else { 776443836d0SMatthew G Knepley ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 777443836d0SMatthew G Knepley ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr); 778443836d0SMatthew G Knepley ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr); 779443836d0SMatthew G Knepley ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr); 780443836d0SMatthew G Knepley } 7813b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7823b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 783c5d2311dSJed Brown } 7843b224e63SBarry Smith PetscFunctionReturn(0); 7853b224e63SBarry Smith } 7863b224e63SBarry Smith 7873b224e63SBarry Smith #undef __FUNCT__ 7880971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 7890971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 7900971522cSBarry Smith { 7910971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7920971522cSBarry Smith PetscErrorCode ierr; 7935a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 794939b8a20SBarry Smith PetscInt cnt,bs; 7950971522cSBarry Smith 7960971522cSBarry Smith PetscFunctionBegin; 7974442daceSBarry Smith x->map->bs = jac->bs; 7984442daceSBarry Smith y->map->bs = jac->bs; 79951f519a2SBarry Smith CHKMEMQ; 80079416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 8011b9fc7fcSBarry Smith if (jac->defaultsplit) { 802939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 803939b8a20SBarry 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); 804939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 805939b8a20SBarry 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); 8060971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 8075a9f2f41SSatish Balay while (ilink) { 8085a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8095a9f2f41SSatish Balay ilink = ilink->next; 8100971522cSBarry Smith } 8110971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 8121b9fc7fcSBarry Smith } else { 813efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8145a9f2f41SSatish Balay while (ilink) { 8155a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8165a9f2f41SSatish Balay ilink = ilink->next; 8171b9fc7fcSBarry Smith } 8181b9fc7fcSBarry Smith } 81916913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 82079416396SBarry Smith if (!jac->w1) { 82179416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 82279416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 82379416396SBarry Smith } 824efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 8255a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 8263e197d65SBarry Smith cnt = 1; 8275a9f2f41SSatish Balay while (ilink->next) { 8285a9f2f41SSatish Balay ilink = ilink->next; 8293e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 8303e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 8313e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 8323e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8333e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8343e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 8353e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8363e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8373e197d65SBarry Smith } 83851f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 83911755939SBarry Smith cnt -= 2; 84051f519a2SBarry Smith while (ilink->previous) { 84151f519a2SBarry Smith ilink = ilink->previous; 84211755939SBarry Smith /* compute the residual only over the part of the vector needed */ 84311755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 84411755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 84511755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 84611755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 84711755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 84811755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 84911755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 85051f519a2SBarry Smith } 85111755939SBarry Smith } 85265e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 85351f519a2SBarry Smith CHKMEMQ; 8540971522cSBarry Smith PetscFunctionReturn(0); 8550971522cSBarry Smith } 8560971522cSBarry Smith 857421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 858ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 859ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 860421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 861ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 862ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 863421e10b8SBarry Smith 864421e10b8SBarry Smith #undef __FUNCT__ 8658c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit" 866421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 867421e10b8SBarry Smith { 868421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 869421e10b8SBarry Smith PetscErrorCode ierr; 870421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 871939b8a20SBarry Smith PetscInt bs; 872421e10b8SBarry Smith 873421e10b8SBarry Smith PetscFunctionBegin; 874421e10b8SBarry Smith CHKMEMQ; 875421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 876421e10b8SBarry Smith if (jac->defaultsplit) { 877939b8a20SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 878939b8a20SBarry 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); 879939b8a20SBarry Smith ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr); 880939b8a20SBarry 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); 881421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 882421e10b8SBarry Smith while (ilink) { 883421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 884421e10b8SBarry Smith ilink = ilink->next; 885421e10b8SBarry Smith } 886421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 887421e10b8SBarry Smith } else { 888421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 889421e10b8SBarry Smith while (ilink) { 890421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 891421e10b8SBarry Smith ilink = ilink->next; 892421e10b8SBarry Smith } 893421e10b8SBarry Smith } 894421e10b8SBarry Smith } else { 895421e10b8SBarry Smith if (!jac->w1) { 896421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 897421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 898421e10b8SBarry Smith } 899421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 900421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 901421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 902421e10b8SBarry Smith while (ilink->next) { 903421e10b8SBarry Smith ilink = ilink->next; 9049989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 905421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 906421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 907421e10b8SBarry Smith } 908421e10b8SBarry Smith while (ilink->previous) { 909421e10b8SBarry Smith ilink = ilink->previous; 9109989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 911421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 912421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 913421e10b8SBarry Smith } 914421e10b8SBarry Smith } else { 915421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 916421e10b8SBarry Smith ilink = ilink->next; 917421e10b8SBarry Smith } 918421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 919421e10b8SBarry Smith while (ilink->previous) { 920421e10b8SBarry Smith ilink = ilink->previous; 9219989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 922421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 923421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 924421e10b8SBarry Smith } 925421e10b8SBarry Smith } 926421e10b8SBarry Smith } 927421e10b8SBarry Smith CHKMEMQ; 928421e10b8SBarry Smith PetscFunctionReturn(0); 929421e10b8SBarry Smith } 930421e10b8SBarry Smith 9310971522cSBarry Smith #undef __FUNCT__ 932574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit" 933574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc) 9340971522cSBarry Smith { 9350971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 9360971522cSBarry Smith PetscErrorCode ierr; 9375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 9380971522cSBarry Smith 9390971522cSBarry Smith PetscFunctionBegin; 9405a9f2f41SSatish Balay while (ilink) { 941574deadeSBarry Smith ierr = KSPReset(ilink->ksp);CHKERRQ(ierr); 942fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->x);CHKERRQ(ierr); 943fcfd50ebSBarry Smith ierr = VecDestroy(&ilink->y);CHKERRQ(ierr); 944443836d0SMatthew G Knepley ierr = VecDestroy(&ilink->z);CHKERRQ(ierr); 945fcfd50ebSBarry Smith ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr); 946fcfd50ebSBarry Smith ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 947b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 9485a9f2f41SSatish Balay next = ilink->next; 9495a9f2f41SSatish Balay ilink = next; 9500971522cSBarry Smith } 95105b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 952574deadeSBarry Smith if (jac->mat && jac->mat != jac->pmat) { 953574deadeSBarry Smith ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr); 954574deadeSBarry Smith } else if (jac->mat) { 955574deadeSBarry Smith jac->mat = PETSC_NULL; 956574deadeSBarry Smith } 95797bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 95868dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 9596bf464f9SBarry Smith ierr = VecDestroy(&jac->w1);CHKERRQ(ierr); 9606bf464f9SBarry Smith ierr = VecDestroy(&jac->w2);CHKERRQ(ierr); 9616bf464f9SBarry Smith ierr = MatDestroy(&jac->schur);CHKERRQ(ierr); 9626bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 9636bf464f9SBarry Smith ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr); 964d78dad28SBarry Smith ierr = KSPDestroy(&jac->kspupper);CHKERRQ(ierr); 9656bf464f9SBarry Smith ierr = MatDestroy(&jac->B);CHKERRQ(ierr); 9666bf464f9SBarry Smith ierr = MatDestroy(&jac->C);CHKERRQ(ierr); 96763ec74ffSBarry Smith jac->reset = PETSC_TRUE; 968574deadeSBarry Smith PetscFunctionReturn(0); 969574deadeSBarry Smith } 970574deadeSBarry Smith 971574deadeSBarry Smith #undef __FUNCT__ 972574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 973574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 974574deadeSBarry Smith { 975574deadeSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 976574deadeSBarry Smith PetscErrorCode ierr; 977574deadeSBarry Smith PC_FieldSplitLink ilink = jac->head,next; 978574deadeSBarry Smith 979574deadeSBarry Smith PetscFunctionBegin; 980574deadeSBarry Smith ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr); 981574deadeSBarry Smith while (ilink) { 9826bf464f9SBarry Smith ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr); 983574deadeSBarry Smith next = ilink->next; 984574deadeSBarry Smith ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 985574deadeSBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 9865d4c12cdSJungho Lee ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr); 987574deadeSBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 988574deadeSBarry Smith ilink = next; 989574deadeSBarry Smith } 990574deadeSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 991c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 9927b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr); 9937b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr); 9947b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr); 9957b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr); 9967b62db95SJungho Lee ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr); 997ab1df9f5SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",PETSC_NULL);CHKERRQ(ierr); 998c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",PETSC_NULL);CHKERRQ(ierr); 9990971522cSBarry Smith PetscFunctionReturn(0); 10000971522cSBarry Smith } 10010971522cSBarry Smith 10020971522cSBarry Smith #undef __FUNCT__ 10030971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 10040971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 10050971522cSBarry Smith { 10061b9fc7fcSBarry Smith PetscErrorCode ierr; 10076c924f48SJed Brown PetscInt bs; 1008bc59fbc5SBarry Smith PetscBool flg,stokes = PETSC_FALSE; 10099dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10103b224e63SBarry Smith PCCompositeType ctype; 1011e7c4fc90SDmitry Karpeev DM ddm; 1012e7c4fc90SDmitry Karpeev char ddm_name[1025]; 10131b9fc7fcSBarry Smith 10140971522cSBarry Smith PetscFunctionBegin; 10151b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 1016e7c4fc90SDmitry Karpeev if (pc->dm) { 1017e7c4fc90SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 1018e7c4fc90SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 1019731c8d9eSDmitry Karpeev ierr = PetscOptionsString("-pc_fieldsplit_decomposition", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 1020e7c4fc90SDmitry Karpeev if (flg) { 102116621825SDmitry Karpeev ierr = DMCreateFieldDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr); 1022e7c4fc90SDmitry Karpeev if (!ddm) { 102316621825SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown field decomposition name %s", ddm_name); 1024e7c4fc90SDmitry Karpeev } 102516621825SDmitry Karpeev ierr = PetscInfo(pc,"Using field decomposition DM defined using options database\n");CHKERRQ(ierr); 1026e7c4fc90SDmitry Karpeev ierr = PCSetDM(pc,ddm); CHKERRQ(ierr); 1027e7c4fc90SDmitry Karpeev } 1028e7c4fc90SDmitry Karpeev } 1029acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 103051f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 103151f519a2SBarry Smith if (flg) { 103251f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 103351f519a2SBarry Smith } 1034704ba839SBarry Smith 1035435f959eSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr); 1036c0adfefeSBarry Smith if (stokes) { 1037c0adfefeSBarry Smith ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1038c0adfefeSBarry Smith jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF; 1039c0adfefeSBarry Smith } 1040c0adfefeSBarry Smith 10413b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 10423b224e63SBarry Smith if (flg) { 10433b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 10443b224e63SBarry Smith } 1045c30613efSMatthew Knepley /* Only setup fields once */ 1046c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 1047d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 1048d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 10496c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 10506c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 1051d32f9abdSBarry Smith } 1052c5d2311dSJed Brown if (jac->type == PC_COMPOSITE_SCHUR) { 1053c9c6ffaaSJed Brown ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr); 1054c9c6ffaaSJed Brown if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);} 1055c9c6ffaaSJed 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); 1056084e4875SJed 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); 1057c5d2311dSJed Brown } 10581b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 10590971522cSBarry Smith PetscFunctionReturn(0); 10600971522cSBarry Smith } 10610971522cSBarry Smith 10620971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 10630971522cSBarry Smith 10640971522cSBarry Smith EXTERN_C_BEGIN 10650971522cSBarry Smith #undef __FUNCT__ 10660971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 10675d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 10680971522cSBarry Smith { 106997bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10700971522cSBarry Smith PetscErrorCode ierr; 10715a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 107269a612a9SBarry Smith char prefix[128]; 10735d4c12cdSJungho Lee PetscInt i; 10740971522cSBarry Smith 10750971522cSBarry Smith PetscFunctionBegin; 10766c924f48SJed Brown if (jac->splitdefined) { 10776c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 10786c924f48SJed Brown PetscFunctionReturn(0); 10796c924f48SJed Brown } 108051f519a2SBarry Smith for (i=0; i<n; i++) { 1081e32f2f54SBarry 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); 1082e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 108351f519a2SBarry Smith } 1084704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1085a04f6461SBarry Smith if (splitname) { 1086db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1087a04f6461SBarry Smith } else { 1088a04f6461SBarry Smith ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1089a04f6461SBarry Smith ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr); 1090a04f6461SBarry Smith } 1091704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 10925a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 10935d4c12cdSJungho Lee ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr); 10945d4c12cdSJungho Lee ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr); 10955a9f2f41SSatish Balay ilink->nfields = n; 10965a9f2f41SSatish Balay ilink->next = PETSC_NULL; 10977adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 109820252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 10995a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11009005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 110169a612a9SBarry Smith 11028caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 11035a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 11040971522cSBarry Smith 11050971522cSBarry Smith if (!next) { 11065a9f2f41SSatish Balay jac->head = ilink; 110751f519a2SBarry Smith ilink->previous = PETSC_NULL; 11080971522cSBarry Smith } else { 11090971522cSBarry Smith while (next->next) { 11100971522cSBarry Smith next = next->next; 11110971522cSBarry Smith } 11125a9f2f41SSatish Balay next->next = ilink; 111351f519a2SBarry Smith ilink->previous = next; 11140971522cSBarry Smith } 11150971522cSBarry Smith jac->nsplits++; 11160971522cSBarry Smith PetscFunctionReturn(0); 11170971522cSBarry Smith } 11180971522cSBarry Smith EXTERN_C_END 11190971522cSBarry Smith 1120e69d4d44SBarry Smith EXTERN_C_BEGIN 1121e69d4d44SBarry Smith #undef __FUNCT__ 1122e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 11237087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 1124e69d4d44SBarry Smith { 1125e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1126e69d4d44SBarry Smith PetscErrorCode ierr; 1127e69d4d44SBarry Smith 1128e69d4d44SBarry Smith PetscFunctionBegin; 1129519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 1130e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 1131e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 113213e0d083SBarry Smith if (n) *n = jac->nsplits; 1133e69d4d44SBarry Smith PetscFunctionReturn(0); 1134e69d4d44SBarry Smith } 1135e69d4d44SBarry Smith EXTERN_C_END 11360971522cSBarry Smith 11370971522cSBarry Smith EXTERN_C_BEGIN 11380971522cSBarry Smith #undef __FUNCT__ 113969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 11407087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 11410971522cSBarry Smith { 11420971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11430971522cSBarry Smith PetscErrorCode ierr; 11440971522cSBarry Smith PetscInt cnt = 0; 11455a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 11460971522cSBarry Smith 11470971522cSBarry Smith PetscFunctionBegin; 11485d480477SMatthew G Knepley ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 11495a9f2f41SSatish Balay while (ilink) { 11505a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 11515a9f2f41SSatish Balay ilink = ilink->next; 11520971522cSBarry Smith } 11535d480477SMatthew 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); 115413e0d083SBarry Smith if (n) *n = jac->nsplits; 11550971522cSBarry Smith PetscFunctionReturn(0); 11560971522cSBarry Smith } 11570971522cSBarry Smith EXTERN_C_END 11580971522cSBarry Smith 1159704ba839SBarry Smith EXTERN_C_BEGIN 1160704ba839SBarry Smith #undef __FUNCT__ 1161704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 11627087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 1163704ba839SBarry Smith { 1164704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1165704ba839SBarry Smith PetscErrorCode ierr; 1166704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 1167704ba839SBarry Smith char prefix[128]; 1168704ba839SBarry Smith 1169704ba839SBarry Smith PetscFunctionBegin; 11706c924f48SJed Brown if (jac->splitdefined) { 11716c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 11726c924f48SJed Brown PetscFunctionReturn(0); 11736c924f48SJed Brown } 117416913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 1175a04f6461SBarry Smith if (splitname) { 1176db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 1177a04f6461SBarry Smith } else { 1178b5787286SJed Brown ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr); 1179b5787286SJed Brown ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr); 1180a04f6461SBarry Smith } 11811ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1182b5787286SJed Brown ierr = ISDestroy(&ilink->is);CHKERRQ(ierr); 1183b5787286SJed Brown ilink->is = is; 1184b5787286SJed Brown ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 1185b5787286SJed Brown ierr = ISDestroy(&ilink->is_col);CHKERRQ(ierr); 1186b5787286SJed Brown ilink->is_col = is; 1187704ba839SBarry Smith ilink->next = PETSC_NULL; 1188704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 118920252d06SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1190704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 11919005cf84SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr); 1192704ba839SBarry Smith 11938caf3d72SBarry Smith ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 1194704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 1195704ba839SBarry Smith 1196704ba839SBarry Smith if (!next) { 1197704ba839SBarry Smith jac->head = ilink; 1198704ba839SBarry Smith ilink->previous = PETSC_NULL; 1199704ba839SBarry Smith } else { 1200704ba839SBarry Smith while (next->next) { 1201704ba839SBarry Smith next = next->next; 1202704ba839SBarry Smith } 1203704ba839SBarry Smith next->next = ilink; 1204704ba839SBarry Smith ilink->previous = next; 1205704ba839SBarry Smith } 1206704ba839SBarry Smith jac->nsplits++; 1207704ba839SBarry Smith 1208704ba839SBarry Smith PetscFunctionReturn(0); 1209704ba839SBarry Smith } 1210704ba839SBarry Smith EXTERN_C_END 1211704ba839SBarry Smith 12120971522cSBarry Smith #undef __FUNCT__ 12130971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 12140971522cSBarry Smith /*@ 12150971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 12160971522cSBarry Smith 1217ad4df100SBarry Smith Logically Collective on PC 12180971522cSBarry Smith 12190971522cSBarry Smith Input Parameters: 12200971522cSBarry Smith + pc - the preconditioner context 1221a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 12220971522cSBarry Smith . n - the number of fields in this split 1223db4c96c1SJed Brown - fields - the fields in this split 12240971522cSBarry Smith 12250971522cSBarry Smith Level: intermediate 12260971522cSBarry Smith 1227d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 1228d32f9abdSBarry Smith 12297287d2a3SDmitry Karpeev The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block 1230d32f9abdSBarry 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 1231d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 1232d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 1233d32f9abdSBarry Smith 1234db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1235db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1236db4c96c1SJed Brown 12375d4c12cdSJungho Lee Developer Note: This routine does not actually create the IS representing the split, that is delayed 12385d4c12cdSJungho Lee until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be 12395d4c12cdSJungho Lee available when this routine is called. 12405d4c12cdSJungho Lee 1241d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 12420971522cSBarry Smith 12430971522cSBarry Smith @*/ 12445d4c12cdSJungho Lee PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col) 12450971522cSBarry Smith { 12464ac538c5SBarry Smith PetscErrorCode ierr; 12470971522cSBarry Smith 12480971522cSBarry Smith PetscFunctionBegin; 12490700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1250db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 1251db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 1252db4c96c1SJed Brown PetscValidIntPointer(fields,3); 12535d4c12cdSJungho Lee ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr); 12540971522cSBarry Smith PetscFunctionReturn(0); 12550971522cSBarry Smith } 12560971522cSBarry Smith 12570971522cSBarry Smith #undef __FUNCT__ 1258704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 1259704ba839SBarry Smith /*@ 1260704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 1261704ba839SBarry Smith 1262ad4df100SBarry Smith Logically Collective on PC 1263704ba839SBarry Smith 1264704ba839SBarry Smith Input Parameters: 1265704ba839SBarry Smith + pc - the preconditioner context 1266a04f6461SBarry Smith . splitname - name of this split, if PETSC_NULL the number of the split is used 1267db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 1268704ba839SBarry Smith 1269d32f9abdSBarry Smith 1270a6ffb8dbSJed Brown Notes: 1271a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 1272a6ffb8dbSJed Brown 1273db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 1274db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 1275d32f9abdSBarry Smith 1276704ba839SBarry Smith Level: intermediate 1277704ba839SBarry Smith 1278704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 1279704ba839SBarry Smith 1280704ba839SBarry Smith @*/ 12817087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 1282704ba839SBarry Smith { 12834ac538c5SBarry Smith PetscErrorCode ierr; 1284704ba839SBarry Smith 1285704ba839SBarry Smith PetscFunctionBegin; 12860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12877b62db95SJungho Lee if (splitname) PetscValidCharPointer(splitname,2); 1288db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 12894ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr); 1290704ba839SBarry Smith PetscFunctionReturn(0); 1291704ba839SBarry Smith } 1292704ba839SBarry Smith 1293704ba839SBarry Smith #undef __FUNCT__ 129457a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS" 129557a9adfeSMatthew G Knepley /*@ 129657a9adfeSMatthew G Knepley PCFieldSplitGetIS - Retrieves the elements for a field as an IS 129757a9adfeSMatthew G Knepley 129857a9adfeSMatthew G Knepley Logically Collective on PC 129957a9adfeSMatthew G Knepley 130057a9adfeSMatthew G Knepley Input Parameters: 130157a9adfeSMatthew G Knepley + pc - the preconditioner context 130257a9adfeSMatthew G Knepley - splitname - name of this split 130357a9adfeSMatthew G Knepley 130457a9adfeSMatthew G Knepley Output Parameter: 130557a9adfeSMatthew G Knepley - is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found 130657a9adfeSMatthew G Knepley 130757a9adfeSMatthew G Knepley Level: intermediate 130857a9adfeSMatthew G Knepley 130957a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS() 131057a9adfeSMatthew G Knepley 131157a9adfeSMatthew G Knepley @*/ 131257a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is) 131357a9adfeSMatthew G Knepley { 131457a9adfeSMatthew G Knepley PetscErrorCode ierr; 131557a9adfeSMatthew G Knepley 131657a9adfeSMatthew G Knepley PetscFunctionBegin; 131757a9adfeSMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 131857a9adfeSMatthew G Knepley PetscValidCharPointer(splitname,2); 131957a9adfeSMatthew G Knepley PetscValidPointer(is,3); 132057a9adfeSMatthew G Knepley { 132157a9adfeSMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 132257a9adfeSMatthew G Knepley PC_FieldSplitLink ilink = jac->head; 132357a9adfeSMatthew G Knepley PetscBool found; 132457a9adfeSMatthew G Knepley 132557a9adfeSMatthew G Knepley *is = PETSC_NULL; 132657a9adfeSMatthew G Knepley while(ilink) { 132757a9adfeSMatthew G Knepley ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr); 132857a9adfeSMatthew G Knepley if (found) { 132957a9adfeSMatthew G Knepley *is = ilink->is; 133057a9adfeSMatthew G Knepley break; 133157a9adfeSMatthew G Knepley } 133257a9adfeSMatthew G Knepley ilink = ilink->next; 133357a9adfeSMatthew G Knepley } 133457a9adfeSMatthew G Knepley } 133557a9adfeSMatthew G Knepley PetscFunctionReturn(0); 133657a9adfeSMatthew G Knepley } 133757a9adfeSMatthew G Knepley 133857a9adfeSMatthew G Knepley #undef __FUNCT__ 133951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 134051f519a2SBarry Smith /*@ 134151f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 134251f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 134351f519a2SBarry Smith 1344ad4df100SBarry Smith Logically Collective on PC 134551f519a2SBarry Smith 134651f519a2SBarry Smith Input Parameters: 134751f519a2SBarry Smith + pc - the preconditioner context 134851f519a2SBarry Smith - bs - the block size 134951f519a2SBarry Smith 135051f519a2SBarry Smith Level: intermediate 135151f519a2SBarry Smith 135251f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 135351f519a2SBarry Smith 135451f519a2SBarry Smith @*/ 13557087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 135651f519a2SBarry Smith { 13574ac538c5SBarry Smith PetscErrorCode ierr; 135851f519a2SBarry Smith 135951f519a2SBarry Smith PetscFunctionBegin; 13600700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1361c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,bs,2); 13624ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr); 136351f519a2SBarry Smith PetscFunctionReturn(0); 136451f519a2SBarry Smith } 136551f519a2SBarry Smith 136651f519a2SBarry Smith #undef __FUNCT__ 136769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 13680971522cSBarry Smith /*@C 136969a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 13700971522cSBarry Smith 137169a612a9SBarry Smith Collective on KSP 13720971522cSBarry Smith 13730971522cSBarry Smith Input Parameter: 13740971522cSBarry Smith . pc - the preconditioner context 13750971522cSBarry Smith 13760971522cSBarry Smith Output Parameters: 137713e0d083SBarry Smith + n - the number of splits 137869a612a9SBarry Smith - pc - the array of KSP contexts 13790971522cSBarry Smith 13800971522cSBarry Smith Note: 1381d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 1382d32f9abdSBarry Smith (not the KSP just the array that contains them). 13830971522cSBarry Smith 138469a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 13850971522cSBarry Smith 1386196cc216SBarry Smith Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs. 1387196cc216SBarry Smith You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_OBJECT,ierr) to determine how large the 1388196cc216SBarry Smith KSP array must be. 1389196cc216SBarry Smith 1390196cc216SBarry Smith 13910971522cSBarry Smith Level: advanced 13920971522cSBarry Smith 13930971522cSBarry Smith .seealso: PCFIELDSPLIT 13940971522cSBarry Smith @*/ 13957087cfbeSBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 13960971522cSBarry Smith { 13974ac538c5SBarry Smith PetscErrorCode ierr; 13980971522cSBarry Smith 13990971522cSBarry Smith PetscFunctionBegin; 14000700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 140113e0d083SBarry Smith if (n) PetscValidIntPointer(n,2); 14024ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr); 14030971522cSBarry Smith PetscFunctionReturn(0); 14040971522cSBarry Smith } 14050971522cSBarry Smith 1406e69d4d44SBarry Smith #undef __FUNCT__ 1407e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 1408e69d4d44SBarry Smith /*@ 1409e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 1410a04f6461SBarry Smith A11 matrix. Otherwise no preconditioner is used. 1411e69d4d44SBarry Smith 1412e69d4d44SBarry Smith Collective on PC 1413e69d4d44SBarry Smith 1414e69d4d44SBarry Smith Input Parameters: 1415e69d4d44SBarry Smith + pc - the preconditioner context 1416fd1303e9SJungho Lee . ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default 1417084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1418084e4875SJed Brown 1419e69d4d44SBarry Smith Options Database: 1420084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1421e69d4d44SBarry Smith 1422fd1303e9SJungho Lee Notes: 1423fd1303e9SJungho Lee $ If ptype is 1424fd1303e9SJungho Lee $ user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument 1425fd1303e9SJungho Lee $ to this function). 1426fd1303e9SJungho Lee $ diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original 1427fd1303e9SJungho Lee $ matrix associated with the Schur complement (i.e. A11) 1428fd1303e9SJungho Lee $ self the preconditioner for the Schur complement is generated from the Schur complement matrix itself: 1429fd1303e9SJungho Lee $ The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC 1430fd1303e9SJungho Lee $ preconditioner 1431fd1303e9SJungho Lee 1432fd1303e9SJungho Lee When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense 1433fd1303e9SJungho Lee with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and 1434fd1303e9SJungho Lee -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement. 1435fd1303e9SJungho Lee 1436fd1303e9SJungho Lee Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in 1437fd1303e9SJungho Lee the name since it sets a proceedure to use. 1438fd1303e9SJungho Lee 1439e69d4d44SBarry Smith Level: intermediate 1440e69d4d44SBarry Smith 1441fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC 1442e69d4d44SBarry Smith 1443e69d4d44SBarry Smith @*/ 14447087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1445e69d4d44SBarry Smith { 14464ac538c5SBarry Smith PetscErrorCode ierr; 1447e69d4d44SBarry Smith 1448e69d4d44SBarry Smith PetscFunctionBegin; 14490700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14504ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr); 1451e69d4d44SBarry Smith PetscFunctionReturn(0); 1452e69d4d44SBarry Smith } 1453e69d4d44SBarry Smith 1454e69d4d44SBarry Smith EXTERN_C_BEGIN 1455e69d4d44SBarry Smith #undef __FUNCT__ 1456e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 14577087cfbeSBarry Smith PetscErrorCode PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1458e69d4d44SBarry Smith { 1459e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1460084e4875SJed Brown PetscErrorCode ierr; 1461e69d4d44SBarry Smith 1462e69d4d44SBarry Smith PetscFunctionBegin; 1463084e4875SJed Brown jac->schurpre = ptype; 1464084e4875SJed Brown if (pre) { 14656bf464f9SBarry Smith ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr); 1466084e4875SJed Brown jac->schur_user = pre; 1467084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1468084e4875SJed Brown } 1469e69d4d44SBarry Smith PetscFunctionReturn(0); 1470e69d4d44SBarry Smith } 1471e69d4d44SBarry Smith EXTERN_C_END 1472e69d4d44SBarry Smith 147330ad9308SMatthew Knepley #undef __FUNCT__ 1474c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType" 1475ab1df9f5SJed Brown /*@ 1476c9c6ffaaSJed Brown PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain 1477ab1df9f5SJed Brown 1478ab1df9f5SJed Brown Collective on PC 1479ab1df9f5SJed Brown 1480ab1df9f5SJed Brown Input Parameters: 1481ab1df9f5SJed Brown + pc - the preconditioner context 1482c9c6ffaaSJed Brown - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default 1483ab1df9f5SJed Brown 1484ab1df9f5SJed Brown Options Database: 1485c9c6ffaaSJed Brown . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full 1486ab1df9f5SJed Brown 1487ab1df9f5SJed Brown 1488ab1df9f5SJed Brown Level: intermediate 1489ab1df9f5SJed Brown 1490ab1df9f5SJed Brown Notes: 1491ab1df9f5SJed Brown The FULL factorization is 1492ab1df9f5SJed Brown 1493ab1df9f5SJed Brown $ (A B) = (1 0) (A 0) (1 Ainv*B) 1494ab1df9f5SJed Brown $ (C D) (C*Ainv 1) (0 S) (0 1 ) 1495ab1df9f5SJed Brown 14966be4592eSBarry 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, 14976be4592eSBarry 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). 1498ab1df9f5SJed Brown 14996be4592eSBarry 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 15006be4592eSBarry 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 15016be4592eSBarry 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, 15026be4592eSBarry 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. 1503ab1df9f5SJed Brown 15046be4592eSBarry 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 15056be4592eSBarry 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). 1506ab1df9f5SJed Brown 1507ab1df9f5SJed Brown References: 1508ab1df9f5SJed Brown Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972. 1509ab1df9f5SJed Brown 1510ab1df9f5SJed Brown Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051. 1511ab1df9f5SJed Brown 1512ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1513ab1df9f5SJed Brown @*/ 1514c9c6ffaaSJed Brown PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype) 1515ab1df9f5SJed Brown { 1516ab1df9f5SJed Brown PetscErrorCode ierr; 1517ab1df9f5SJed Brown 1518ab1df9f5SJed Brown PetscFunctionBegin; 1519ab1df9f5SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1520c9c6ffaaSJed Brown ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr); 1521ab1df9f5SJed Brown PetscFunctionReturn(0); 1522ab1df9f5SJed Brown } 1523ab1df9f5SJed Brown 1524ab1df9f5SJed Brown #undef __FUNCT__ 1525c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit" 1526c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype) 1527ab1df9f5SJed Brown { 1528ab1df9f5SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1529ab1df9f5SJed Brown 1530ab1df9f5SJed Brown PetscFunctionBegin; 1531ab1df9f5SJed Brown jac->schurfactorization = ftype; 1532ab1df9f5SJed Brown PetscFunctionReturn(0); 1533ab1df9f5SJed Brown } 1534ab1df9f5SJed Brown 1535ab1df9f5SJed Brown #undef __FUNCT__ 153630ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 153730ad9308SMatthew Knepley /*@C 15388c03b21aSDmitry Karpeev PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement 153930ad9308SMatthew Knepley 154030ad9308SMatthew Knepley Collective on KSP 154130ad9308SMatthew Knepley 154230ad9308SMatthew Knepley Input Parameter: 154330ad9308SMatthew Knepley . pc - the preconditioner context 154430ad9308SMatthew Knepley 154530ad9308SMatthew Knepley Output Parameters: 1546a04f6461SBarry Smith + A00 - the (0,0) block 1547a04f6461SBarry Smith . A01 - the (0,1) block 1548a04f6461SBarry Smith . A10 - the (1,0) block 1549a04f6461SBarry Smith - A11 - the (1,1) block 155030ad9308SMatthew Knepley 155130ad9308SMatthew Knepley Level: advanced 155230ad9308SMatthew Knepley 155330ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 155430ad9308SMatthew Knepley @*/ 1555a04f6461SBarry Smith PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11) 155630ad9308SMatthew Knepley { 155730ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 155830ad9308SMatthew Knepley 155930ad9308SMatthew Knepley PetscFunctionBegin; 15600700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1561c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 1562a04f6461SBarry Smith if (A00) *A00 = jac->pmat[0]; 1563a04f6461SBarry Smith if (A01) *A01 = jac->B; 1564a04f6461SBarry Smith if (A10) *A10 = jac->C; 1565a04f6461SBarry Smith if (A11) *A11 = jac->pmat[1]; 156630ad9308SMatthew Knepley PetscFunctionReturn(0); 156730ad9308SMatthew Knepley } 156830ad9308SMatthew Knepley 156979416396SBarry Smith EXTERN_C_BEGIN 157079416396SBarry Smith #undef __FUNCT__ 157179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 15727087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 157379416396SBarry Smith { 157479416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1575e69d4d44SBarry Smith PetscErrorCode ierr; 157679416396SBarry Smith 157779416396SBarry Smith PetscFunctionBegin; 157879416396SBarry Smith jac->type = type; 15793b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 15803b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 15813b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1582e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1583e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1584c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr); 1585e69d4d44SBarry Smith 15863b224e63SBarry Smith } else { 15873b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 15883b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1589e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 15909e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 1591c9c6ffaaSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr); 15923b224e63SBarry Smith } 159379416396SBarry Smith PetscFunctionReturn(0); 159479416396SBarry Smith } 159579416396SBarry Smith EXTERN_C_END 159679416396SBarry Smith 159751f519a2SBarry Smith EXTERN_C_BEGIN 159851f519a2SBarry Smith #undef __FUNCT__ 159951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 16007087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 160151f519a2SBarry Smith { 160251f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 160351f519a2SBarry Smith 160451f519a2SBarry Smith PetscFunctionBegin; 160565e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 160665e19b50SBarry 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); 160751f519a2SBarry Smith jac->bs = bs; 160851f519a2SBarry Smith PetscFunctionReturn(0); 160951f519a2SBarry Smith } 161051f519a2SBarry Smith EXTERN_C_END 161151f519a2SBarry Smith 161279416396SBarry Smith #undef __FUNCT__ 161379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1614bc08b0f1SBarry Smith /*@ 161579416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 161679416396SBarry Smith 161779416396SBarry Smith Collective on PC 161879416396SBarry Smith 161979416396SBarry Smith Input Parameter: 162079416396SBarry Smith . pc - the preconditioner context 162181540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 162279416396SBarry Smith 162379416396SBarry Smith Options Database Key: 1624a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 162579416396SBarry Smith 1626b02e2d75SMatthew G Knepley Level: Intermediate 162779416396SBarry Smith 162879416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 162979416396SBarry Smith 163079416396SBarry Smith .seealso: PCCompositeSetType() 163179416396SBarry Smith 163279416396SBarry Smith @*/ 16337087cfbeSBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type) 163479416396SBarry Smith { 16354ac538c5SBarry Smith PetscErrorCode ierr; 163679416396SBarry Smith 163779416396SBarry Smith PetscFunctionBegin; 16380700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16394ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 164079416396SBarry Smith PetscFunctionReturn(0); 164179416396SBarry Smith } 164279416396SBarry Smith 1643b02e2d75SMatthew G Knepley #undef __FUNCT__ 1644b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType" 1645b02e2d75SMatthew G Knepley /*@ 1646b02e2d75SMatthew G Knepley PCFieldSplitGetType - Gets the type of fieldsplit preconditioner. 1647b02e2d75SMatthew G Knepley 1648b02e2d75SMatthew G Knepley Not collective 1649b02e2d75SMatthew G Knepley 1650b02e2d75SMatthew G Knepley Input Parameter: 1651b02e2d75SMatthew G Knepley . pc - the preconditioner context 1652b02e2d75SMatthew G Knepley 1653b02e2d75SMatthew G Knepley Output Parameter: 1654b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 1655b02e2d75SMatthew G Knepley 1656b02e2d75SMatthew G Knepley Level: Intermediate 1657b02e2d75SMatthew G Knepley 1658b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative 1659b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType() 1660b02e2d75SMatthew G Knepley @*/ 1661b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type) 1662b02e2d75SMatthew G Knepley { 1663b02e2d75SMatthew G Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 1664b02e2d75SMatthew G Knepley 1665b02e2d75SMatthew G Knepley PetscFunctionBegin; 1666b02e2d75SMatthew G Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1667b02e2d75SMatthew G Knepley PetscValidIntPointer(type,2); 1668b02e2d75SMatthew G Knepley *type = jac->type; 1669b02e2d75SMatthew G Knepley PetscFunctionReturn(0); 1670b02e2d75SMatthew G Knepley } 1671b02e2d75SMatthew G Knepley 16720971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 16730971522cSBarry Smith /*MC 1674a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 1675a04f6461SBarry Smith fields or groups of fields. See the users manual section "Solving Block Matrices" for more details. 16760971522cSBarry Smith 1677edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1678edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 16790971522cSBarry Smith 1680a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 168169a612a9SBarry Smith and set the options directly on the resulting KSP object 16820971522cSBarry Smith 16830971522cSBarry Smith Level: intermediate 16840971522cSBarry Smith 168579416396SBarry Smith Options Database Keys: 168681540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 168781540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 168881540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 168981540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 16900f188ba9SJed Brown . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting 16910f188ba9SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> - default is diag 1692435f959eSBarry Smith . -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver 169379416396SBarry Smith 16945d4c12cdSJungho Lee - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 16955d4c12cdSJungho Lee for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 16965d4c12cdSJungho Lee 1697c8a0d604SMatthew G Knepley Notes: 1698c8a0d604SMatthew G Knepley Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1699d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1700d32f9abdSBarry Smith 1701d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1702d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1703d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1704d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1705d32f9abdSBarry Smith 1706c8a0d604SMatthew G Knepley $ For the Schur complement preconditioner if J = ( A00 A01 ) 1707c8a0d604SMatthew G Knepley $ ( A10 A11 ) 1708c8a0d604SMatthew G Knepley $ the preconditioner using full factorization is 1709c8a0d604SMatthew G Knepley $ ( I -A10 ksp(A00) ) ( inv(A00) 0 ) ( I 0 ) 1710c8a0d604SMatthew G Knepley $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I ) 1711a04f6461SBarry Smith where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1712c8a0d604SMatthew G Knepley ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given 1713c8a0d604SMatthew G Knepley in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0, 1714c8a0d604SMatthew G Knepley it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1715a04f6461SBarry Smith A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1716c8a0d604SMatthew G Knepley option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The 1717c9c6ffaaSJed Brown factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above, 17185668aaf4SBarry Smith diag gives 1719c8a0d604SMatthew G Knepley $ ( inv(A00) 0 ) 1720c8a0d604SMatthew G Knepley $ ( 0 -ksp(S) ) 17215668aaf4SBarry 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 1722c8a0d604SMatthew G Knepley $ ( A00 0 ) 1723c8a0d604SMatthew G Knepley $ ( A10 S ) 1724c8a0d604SMatthew G Knepley where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of 1725c8a0d604SMatthew G Knepley $ ( A00 A01 ) 1726c8a0d604SMatthew G Knepley $ ( 0 S ) 1727c8a0d604SMatthew G Knepley where again the inverses of A00 and S are applied using KSPs. 1728e69d4d44SBarry Smith 1729edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1730edf189efSBarry Smith is used automatically for a second block. 1731edf189efSBarry Smith 1732ff218e97SBarry Smith The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. 1733ff218e97SBarry Smith Generally it should be used with the AIJ format. 1734ff218e97SBarry Smith 1735ff218e97SBarry Smith The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see, 1736ff218e97SBarry Smith for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT 1737ff218e97SBarry Smith inside a smoother resulting in "Distributive Smoothers". 17380716a85fSBarry Smith 1739a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 17400971522cSBarry Smith 17417e8cb189SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC, 1742e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 17430971522cSBarry Smith M*/ 17440971522cSBarry Smith 17450971522cSBarry Smith EXTERN_C_BEGIN 17460971522cSBarry Smith #undef __FUNCT__ 17470971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 17487087cfbeSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 17490971522cSBarry Smith { 17500971522cSBarry Smith PetscErrorCode ierr; 17510971522cSBarry Smith PC_FieldSplit *jac; 17520971522cSBarry Smith 17530971522cSBarry Smith PetscFunctionBegin; 175438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 17550971522cSBarry Smith jac->bs = -1; 17560971522cSBarry Smith jac->nsplits = 0; 17573e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1758e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 1759c9c6ffaaSJed Brown jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL; 176051f519a2SBarry Smith 17610971522cSBarry Smith pc->data = (void*)jac; 17620971522cSBarry Smith 17630971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1764421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 17650971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 1766574deadeSBarry Smith pc->ops->reset = PCReset_FieldSplit; 17670971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 17680971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 17690971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 17700971522cSBarry Smith pc->ops->applyrichardson = 0; 17710971522cSBarry Smith 177269a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 177369a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 17740971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 17750971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1776704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1777704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 177879416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 177979416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 178051f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 178151f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 17820971522cSBarry Smith PetscFunctionReturn(0); 17830971522cSBarry Smith } 17840971522cSBarry Smith EXTERN_C_END 17850971522cSBarry Smith 17860971522cSBarry Smith 1787a541d17aSBarry Smith 1788