xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 1720468b2959f38b91bb4b2d5f82b6e6f9336da9)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    Solvers
4a0d3c3abSStefano Zampini    - Add support for cholesky for coarse solver (similar to local solvers)
5eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
6eb97c9d2SStefano Zampini 
7eb97c9d2SStefano Zampini    User interface
80f202f7eSStefano Zampini    - ** DM attached to pc?
9eb97c9d2SStefano Zampini 
10eb97c9d2SStefano Zampini    Debugging output
11b9b85e73SStefano Zampini    - * Better management of verbosity levels of debugging output
12eb97c9d2SStefano Zampini 
13eb97c9d2SStefano Zampini    Extra
14b9b85e73SStefano Zampini    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15eb97c9d2SStefano Zampini    - BDDC with MG framework?
16eb97c9d2SStefano Zampini 
17eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
18eb97c9d2SStefano Zampini    - Provide general case for subassembling
19eb97c9d2SStefano Zampini 
2053cdbc3dSStefano Zampini */
210c7d97c5SJed Brown 
22ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
23ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
243b03a366Sstefano_zampini #include <petscblaslapack.h>
25674ae819SStefano Zampini 
260c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
270c7d97c5SJed Brown #undef __FUNCT__
280c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
294416b707SBarry Smith PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc)
300c7d97c5SJed Brown {
310c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
320c7d97c5SJed Brown   PetscErrorCode ierr;
330c7d97c5SJed Brown 
340c7d97c5SJed Brown   PetscFunctionBegin;
35e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
368eeda7d8SStefano Zampini   /* Verbose debugging */
37c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_dirichlet_approximate","Inform PCBDDC that we are using approximate Dirichlet solvers","none",pcbddc->NullSpace_corr[0],&pcbddc->NullSpace_corr[0],NULL);CHKERRQ(ierr);
38c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_dirichlet_approximate_scale","Inform PCBDDC that we need to scale the Dirichlet solve","none",pcbddc->NullSpace_corr[1],&pcbddc->NullSpace_corr[1],NULL);CHKERRQ(ierr);
39c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_neumann_approximate","Inform PCBDDC that we are using approximate Neumann solvers","none",pcbddc->NullSpace_corr[2],&pcbddc->NullSpace_corr[2],NULL);CHKERRQ(ierr);
40c7017625SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_neumann_approximate_scale","Inform PCBDDC that we need to scale the Neumann solve","none",pcbddc->NullSpace_corr[3],&pcbddc->NullSpace_corr[3],NULL);CHKERRQ(ierr);
418eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
426b78500eSPatrick Sanan   /* Primal space customization */
4308a5cf49SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr);
448eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
458eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
468eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
476661aa1dSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);CHKERRQ(ierr);
48ac632422SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr);
498eeda7d8SStefano Zampini   /* Change of basis */
50b9b85e73SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
51b9b85e73SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
52674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
53674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
54674ae819SStefano Zampini   }
558eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
568eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);CHKERRQ(ierr);
5774e2c79eSStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);CHKERRQ(ierr);
580298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
592b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
60323d395dSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr);
61674ae819SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
62b96c3477SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);CHKERRQ(ierr);
63b96c3477SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);CHKERRQ(ierr);
64b96c3477SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);CHKERRQ(ierr);
65ac632422SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);CHKERRQ(ierr);
664c6709b3SStefano Zampini   ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr);
6708122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
6808122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
693301b35fSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
70b0c7d250SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);CHKERRQ(ierr);
710c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
720c7d97c5SJed Brown   PetscFunctionReturn(0);
730c7d97c5SJed Brown }
746b78500eSPatrick Sanan 
756b78500eSPatrick Sanan /* -------------------------------------------------------------------------- */
766b78500eSPatrick Sanan #undef __FUNCT__
776b78500eSPatrick Sanan #define __FUNCT__ "PCView_BDDC"
786b78500eSPatrick Sanan static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer)
796b78500eSPatrick Sanan {
806b78500eSPatrick Sanan   PC_BDDC              *pcbddc = (PC_BDDC*)pc->data;
816b78500eSPatrick Sanan   PetscErrorCode       ierr;
826b78500eSPatrick Sanan   PetscBool            isascii,isstring;
836b78500eSPatrick Sanan 
846b78500eSPatrick Sanan   PetscFunctionBegin;
856b78500eSPatrick Sanan   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
866b78500eSPatrick Sanan   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
876b78500eSPatrick Sanan 
886b78500eSPatrick Sanan   /* In a braindead way, print out anything which the user can control from the command line,
896b78500eSPatrick Sanan      cribbing from PCSetFromOptions_BDDC */
906b78500eSPatrick Sanan 
916b78500eSPatrick Sanan   /* Nothing printed for the String viewer */
926b78500eSPatrick Sanan 
936b78500eSPatrick Sanan   /* ASCII viewer */
946b78500eSPatrick Sanan   if (isascii) {
956b78500eSPatrick Sanan     /* Verbose debugging */
966b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use verbose output: %d\n",pcbddc->dbg_flag);CHKERRQ(ierr);
976b78500eSPatrick Sanan 
986b78500eSPatrick Sanan     /* Primal space customization */
996b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use local mat graph: %d\n",pcbddc->use_local_adj);CHKERRQ(ierr);
1006b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use vertices: %d\n",pcbddc->use_vertices);CHKERRQ(ierr);
1016b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use edges: %d\n",pcbddc->use_edges);CHKERRQ(ierr);
1026b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use faces: %d\n",pcbddc->use_faces);CHKERRQ(ierr);
1036b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use true near null space: %d\n",pcbddc->use_nnsp_true);CHKERRQ(ierr);
1046b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);CHKERRQ(ierr);
1056b78500eSPatrick Sanan 
1066b78500eSPatrick Sanan     /* Change of basis */
1076b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);CHKERRQ(ierr);
1086b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);CHKERRQ(ierr);
1096b78500eSPatrick Sanan 
1106b78500eSPatrick Sanan     /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
1116b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);CHKERRQ(ierr);
1126b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Coarse problem restribute procs: %d\n",pcbddc->redistribute_coarse);CHKERRQ(ierr);
1136b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Multilevel coarsening ratio: %d\n",pcbddc->coarsening_ratio);CHKERRQ(ierr);
1146b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Multilevel max levels: %d\n",pcbddc->max_levels);CHKERRQ(ierr);
1156b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);CHKERRQ(ierr);
1166b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);CHKERRQ(ierr);
1176b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);CHKERRQ(ierr);
1186b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Number of dofs' layers for the computation of principal minors: %d\n",pcbddc->sub_schurs_layers);CHKERRQ(ierr);
1196b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);CHKERRQ(ierr);
1206b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Fast deluxe scaling: %d\n",pcbddc->faster_deluxe);CHKERRQ(ierr);
1216b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Adaptive constraint selection threshold: %g\n",pcbddc->adaptive_threshold);CHKERRQ(ierr);
1226b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Min constraints / connected component: %d\n",pcbddc->adaptive_nmin);CHKERRQ(ierr);
1236b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Max constraints / connected component: %d\n",pcbddc->adaptive_nmax);CHKERRQ(ierr);
1246b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);CHKERRQ(ierr);
1256b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Num. Procs. to map coarse adjacency list: %d\n",pcbddc->coarse_adj_red);CHKERRQ(ierr);
1266b78500eSPatrick Sanan   }
1276b78500eSPatrick Sanan 
1286b78500eSPatrick Sanan   PetscFunctionReturn(0);
1296b78500eSPatrick Sanan }
1306b78500eSPatrick Sanan 
1310c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
132674ae819SStefano Zampini #undef __FUNCT__
133906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
134906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
135b9b85e73SStefano Zampini {
136b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
137b9b85e73SStefano Zampini   PetscErrorCode ierr;
138b9b85e73SStefano Zampini 
139b9b85e73SStefano Zampini   PetscFunctionBegin;
140b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
141b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
142b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
143b9b85e73SStefano Zampini   PetscFunctionReturn(0);
144b9b85e73SStefano Zampini }
145b9b85e73SStefano Zampini #undef __FUNCT__
146906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
147b9b85e73SStefano Zampini /*@
148906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
149b9b85e73SStefano Zampini 
150b9b85e73SStefano Zampini    Collective on PC
151b9b85e73SStefano Zampini 
152b9b85e73SStefano Zampini    Input Parameters:
153b9b85e73SStefano Zampini +  pc - the preconditioning context
154906d46d4SStefano Zampini -  change - the change of basis matrix
155b9b85e73SStefano Zampini 
156b9b85e73SStefano Zampini    Level: intermediate
157b9b85e73SStefano Zampini 
158b9b85e73SStefano Zampini    Notes:
159b9b85e73SStefano Zampini 
160b9b85e73SStefano Zampini .seealso: PCBDDC
161b9b85e73SStefano Zampini @*/
162906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
163b9b85e73SStefano Zampini {
164b9b85e73SStefano Zampini   PetscErrorCode ierr;
165b9b85e73SStefano Zampini 
166b9b85e73SStefano Zampini   PetscFunctionBegin;
167b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
168b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
169906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
170906d46d4SStefano Zampini   if (pc->mat) {
171906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
172906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
173906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
1746c4ed002SBarry Smith     if (rows_c != rows) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
1756c4ed002SBarry Smith     if (cols_c != cols) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
176906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
177906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
1786c4ed002SBarry Smith     if (rows_c != rows) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
1796c4ed002SBarry Smith     if (cols_c != cols) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
180906d46d4SStefano Zampini   }
181906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
182b9b85e73SStefano Zampini   PetscFunctionReturn(0);
183b9b85e73SStefano Zampini }
184b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
185b9b85e73SStefano Zampini #undef __FUNCT__
186674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
187674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
188674ae819SStefano Zampini {
189674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
190674ae819SStefano Zampini   PetscErrorCode ierr;
1911e6b0712SBarry Smith 
192674ae819SStefano Zampini   PetscFunctionBegin;
193674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
194674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
195674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
196674ae819SStefano Zampini   PetscFunctionReturn(0);
197674ae819SStefano Zampini }
198674ae819SStefano Zampini #undef __FUNCT__
199674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
200674ae819SStefano Zampini /*@
20128509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
202674ae819SStefano Zampini 
20317eb1463SStefano Zampini    Collective
204674ae819SStefano Zampini 
205674ae819SStefano Zampini    Input Parameters:
206674ae819SStefano Zampini +  pc - the preconditioning context
20717eb1463SStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
208674ae819SStefano Zampini 
209674ae819SStefano Zampini    Level: intermediate
210674ae819SStefano Zampini 
211674ae819SStefano Zampini    Notes:
212674ae819SStefano Zampini 
213674ae819SStefano Zampini .seealso: PCBDDC
214674ae819SStefano Zampini @*/
215674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
216674ae819SStefano Zampini {
217674ae819SStefano Zampini   PetscErrorCode ierr;
218674ae819SStefano Zampini 
219674ae819SStefano Zampini   PetscFunctionBegin;
220674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
221674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
22217eb1463SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
223674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
224674ae819SStefano Zampini   PetscFunctionReturn(0);
225674ae819SStefano Zampini }
226674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
2270c7d97c5SJed Brown #undef __FUNCT__
2284fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
2294fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
2304fad6a16SStefano Zampini {
2314fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2324fad6a16SStefano Zampini 
2334fad6a16SStefano Zampini   PetscFunctionBegin;
2344fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
2354fad6a16SStefano Zampini   PetscFunctionReturn(0);
2364fad6a16SStefano Zampini }
2371e6b0712SBarry Smith 
2384fad6a16SStefano Zampini #undef __FUNCT__
2394fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
2404fad6a16SStefano Zampini /*@
24128509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
2424fad6a16SStefano Zampini 
2434fad6a16SStefano Zampini    Logically collective on PC
2444fad6a16SStefano Zampini 
2454fad6a16SStefano Zampini    Input Parameters:
2464fad6a16SStefano Zampini +  pc - the preconditioning context
24728509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2484fad6a16SStefano Zampini 
2490f202f7eSStefano Zampini    Options Database Keys:
2500f202f7eSStefano Zampini .    -pc_bddc_coarsening_ratio
2514fad6a16SStefano Zampini 
2524fad6a16SStefano Zampini    Level: intermediate
2534fad6a16SStefano Zampini 
2544fad6a16SStefano Zampini    Notes:
2550f202f7eSStefano Zampini      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2564fad6a16SStefano Zampini 
2570f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetLevels()
2584fad6a16SStefano Zampini @*/
2594fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2604fad6a16SStefano Zampini {
2614fad6a16SStefano Zampini   PetscErrorCode ierr;
2624fad6a16SStefano Zampini 
2634fad6a16SStefano Zampini   PetscFunctionBegin;
2644fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2652b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2664fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2674fad6a16SStefano Zampini   PetscFunctionReturn(0);
2684fad6a16SStefano Zampini }
2692b510759SStefano Zampini 
270b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2712b510759SStefano Zampini #undef __FUNCT__
272b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
273b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
274b8ffe317SStefano Zampini {
275b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
276b8ffe317SStefano Zampini 
277b8ffe317SStefano Zampini   PetscFunctionBegin;
27885c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
279b8ffe317SStefano Zampini   PetscFunctionReturn(0);
280b8ffe317SStefano Zampini }
281b8ffe317SStefano Zampini 
282b8ffe317SStefano Zampini #undef __FUNCT__
283b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
284b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2852b510759SStefano Zampini {
2862b510759SStefano Zampini   PetscErrorCode ierr;
2872b510759SStefano Zampini 
2882b510759SStefano Zampini   PetscFunctionBegin;
2892b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
290b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
291b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2922b510759SStefano Zampini   PetscFunctionReturn(0);
2932b510759SStefano Zampini }
2941e6b0712SBarry Smith 
2954fad6a16SStefano Zampini #undef __FUNCT__
2962b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2972b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2984fad6a16SStefano Zampini {
2994fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3004fad6a16SStefano Zampini 
3014fad6a16SStefano Zampini   PetscFunctionBegin;
3022b510759SStefano Zampini   pcbddc->current_level = level;
3034fad6a16SStefano Zampini   PetscFunctionReturn(0);
3044fad6a16SStefano Zampini }
3051e6b0712SBarry Smith 
3064fad6a16SStefano Zampini #undef __FUNCT__
307b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
308b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
309b8ffe317SStefano Zampini {
310b8ffe317SStefano Zampini   PetscErrorCode ierr;
311b8ffe317SStefano Zampini 
312b8ffe317SStefano Zampini   PetscFunctionBegin;
313b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
314b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
315b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
316b8ffe317SStefano Zampini   PetscFunctionReturn(0);
317b8ffe317SStefano Zampini }
318b8ffe317SStefano Zampini 
319b8ffe317SStefano Zampini #undef __FUNCT__
3202b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
3212b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
3222b510759SStefano Zampini {
3232b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3242b510759SStefano Zampini 
3252b510759SStefano Zampini   PetscFunctionBegin;
3262b510759SStefano Zampini   pcbddc->max_levels = levels;
3272b510759SStefano Zampini   PetscFunctionReturn(0);
3282b510759SStefano Zampini }
3292b510759SStefano Zampini 
3302b510759SStefano Zampini #undef __FUNCT__
3312b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
3324fad6a16SStefano Zampini /*@
33328509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
3344fad6a16SStefano Zampini 
3354fad6a16SStefano Zampini    Logically collective on PC
3364fad6a16SStefano Zampini 
3374fad6a16SStefano Zampini    Input Parameters:
3384fad6a16SStefano Zampini +  pc - the preconditioning context
33928509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
3404fad6a16SStefano Zampini 
3410f202f7eSStefano Zampini    Options Database Keys:
3420f202f7eSStefano Zampini .    -pc_bddc_levels
3434fad6a16SStefano Zampini 
3444fad6a16SStefano Zampini    Level: intermediate
3454fad6a16SStefano Zampini 
3464fad6a16SStefano Zampini    Notes:
3470f202f7eSStefano Zampini      Default value is 0, i.e. traditional one-level BDDC
3484fad6a16SStefano Zampini 
3490f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
3504fad6a16SStefano Zampini @*/
3512b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3524fad6a16SStefano Zampini {
3534fad6a16SStefano Zampini   PetscErrorCode ierr;
3544fad6a16SStefano Zampini 
3554fad6a16SStefano Zampini   PetscFunctionBegin;
3564fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3572b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
358312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3592b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3604fad6a16SStefano Zampini   PetscFunctionReturn(0);
3614fad6a16SStefano Zampini }
3624fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3631e6b0712SBarry Smith 
3644fad6a16SStefano Zampini #undef __FUNCT__
3653b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3663b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3673b03a366Sstefano_zampini {
3683b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3693b03a366Sstefano_zampini   PetscErrorCode ierr;
3703b03a366Sstefano_zampini 
3713b03a366Sstefano_zampini   PetscFunctionBegin;
372785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
373785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3743b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
37536e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
37636e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
377fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3783b03a366Sstefano_zampini   PetscFunctionReturn(0);
3793b03a366Sstefano_zampini }
3801e6b0712SBarry Smith 
3813b03a366Sstefano_zampini #undef __FUNCT__
3823b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3833b03a366Sstefano_zampini /*@
38428509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3853b03a366Sstefano_zampini 
386785d1243SStefano Zampini    Collective
3873b03a366Sstefano_zampini 
3883b03a366Sstefano_zampini    Input Parameters:
3893b03a366Sstefano_zampini +  pc - the preconditioning context
390785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3913b03a366Sstefano_zampini 
3923b03a366Sstefano_zampini    Level: intermediate
3933b03a366Sstefano_zampini 
3940f202f7eSStefano Zampini    Notes:
3950f202f7eSStefano Zampini      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
3963b03a366Sstefano_zampini 
3970f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
3983b03a366Sstefano_zampini @*/
3993b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
4003b03a366Sstefano_zampini {
4013b03a366Sstefano_zampini   PetscErrorCode ierr;
4023b03a366Sstefano_zampini 
4033b03a366Sstefano_zampini   PetscFunctionBegin;
4043b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
405674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
406785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4073b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4083b03a366Sstefano_zampini   PetscFunctionReturn(0);
4093b03a366Sstefano_zampini }
4103b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4111e6b0712SBarry Smith 
4123b03a366Sstefano_zampini #undef __FUNCT__
41382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
41482ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4150c7d97c5SJed Brown {
4160c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4170c7d97c5SJed Brown   PetscErrorCode ierr;
4180c7d97c5SJed Brown 
4190c7d97c5SJed Brown   PetscFunctionBegin;
420785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
421785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4220c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4230c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
424785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4257d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4260c7d97c5SJed Brown   PetscFunctionReturn(0);
4270c7d97c5SJed Brown }
4280c7d97c5SJed Brown 
4290c7d97c5SJed Brown #undef __FUNCT__
43082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4319c0446d6SStefano Zampini /*@
43282ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4339c0446d6SStefano Zampini 
434785d1243SStefano Zampini    Collective
43553cdbc3dSStefano Zampini 
43653cdbc3dSStefano Zampini    Input Parameters:
43753cdbc3dSStefano Zampini +  pc - the preconditioning context
43882ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
43953cdbc3dSStefano Zampini 
44053cdbc3dSStefano Zampini    Level: intermediate
44153cdbc3dSStefano Zampini 
4429c0446d6SStefano Zampini    Notes:
44353cdbc3dSStefano Zampini 
4440f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
44553cdbc3dSStefano Zampini @*/
44682ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4470c7d97c5SJed Brown {
4480c7d97c5SJed Brown   PetscErrorCode ierr;
4490c7d97c5SJed Brown 
4500c7d97c5SJed Brown   PetscFunctionBegin;
4510c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4520c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
45382ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
45482ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4550c7d97c5SJed Brown   PetscFunctionReturn(0);
4560c7d97c5SJed Brown }
4570c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4580c7d97c5SJed Brown 
4590c7d97c5SJed Brown #undef __FUNCT__
4600c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
46153cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4620c7d97c5SJed Brown {
4630c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
46453cdbc3dSStefano Zampini   PetscErrorCode ierr;
4650c7d97c5SJed Brown 
4660c7d97c5SJed Brown   PetscFunctionBegin;
467785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
468785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
46953cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
47036e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
47136e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
472fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4730c7d97c5SJed Brown   PetscFunctionReturn(0);
4740c7d97c5SJed Brown }
4751e6b0712SBarry Smith 
4760c7d97c5SJed Brown #undef __FUNCT__
4770c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
47857527edcSJed Brown /*@
47928509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
48057527edcSJed Brown 
481785d1243SStefano Zampini    Collective
48257527edcSJed Brown 
48357527edcSJed Brown    Input Parameters:
48457527edcSJed Brown +  pc - the preconditioning context
485785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
48657527edcSJed Brown 
48757527edcSJed Brown    Level: intermediate
48857527edcSJed Brown 
4890f202f7eSStefano Zampini    Notes:
4900f202f7eSStefano Zampini      Any process can list any global node
49157527edcSJed Brown 
4920f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
49357527edcSJed Brown @*/
49453cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4950c7d97c5SJed Brown {
4960c7d97c5SJed Brown   PetscErrorCode ierr;
4970c7d97c5SJed Brown 
4980c7d97c5SJed Brown   PetscFunctionBegin;
4990c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
500674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
501785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
50253cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
50353cdbc3dSStefano Zampini   PetscFunctionReturn(0);
50453cdbc3dSStefano Zampini }
50553cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5061e6b0712SBarry Smith 
50753cdbc3dSStefano Zampini #undef __FUNCT__
50882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
50982ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5100c7d97c5SJed Brown {
5110c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5120c7d97c5SJed Brown   PetscErrorCode ierr;
5130c7d97c5SJed Brown 
5140c7d97c5SJed Brown   PetscFunctionBegin;
515785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
516785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5170c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5180c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
519785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5207d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5210c7d97c5SJed Brown   PetscFunctionReturn(0);
5220c7d97c5SJed Brown }
5230c7d97c5SJed Brown 
5240c7d97c5SJed Brown #undef __FUNCT__
52582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5260c7d97c5SJed Brown /*@
52782ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5280c7d97c5SJed Brown 
529785d1243SStefano Zampini    Collective
5300c7d97c5SJed Brown 
5310c7d97c5SJed Brown    Input Parameters:
5320c7d97c5SJed Brown +  pc - the preconditioning context
53382ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5340c7d97c5SJed Brown 
5350c7d97c5SJed Brown    Level: intermediate
5360c7d97c5SJed Brown 
5370c7d97c5SJed Brown    Notes:
5380c7d97c5SJed Brown 
5390f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
5400c7d97c5SJed Brown @*/
54182ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5420c7d97c5SJed Brown {
5430c7d97c5SJed Brown   PetscErrorCode ierr;
5440c7d97c5SJed Brown 
5450c7d97c5SJed Brown   PetscFunctionBegin;
5460c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5470c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
54882ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54982ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
55053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
55153cdbc3dSStefano Zampini }
55253cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
55353cdbc3dSStefano Zampini 
55453cdbc3dSStefano Zampini #undef __FUNCT__
555da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
556da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
557da1bb401SStefano Zampini {
558da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
559da1bb401SStefano Zampini 
560da1bb401SStefano Zampini   PetscFunctionBegin;
561da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
562da1bb401SStefano Zampini   PetscFunctionReturn(0);
563da1bb401SStefano Zampini }
5641e6b0712SBarry Smith 
565da1bb401SStefano Zampini #undef __FUNCT__
566da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
567da1bb401SStefano Zampini /*@
568785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
569da1bb401SStefano Zampini 
570785d1243SStefano Zampini    Collective
571785d1243SStefano Zampini 
572785d1243SStefano Zampini    Input Parameters:
573785d1243SStefano Zampini .  pc - the preconditioning context
574785d1243SStefano Zampini 
575785d1243SStefano Zampini    Output Parameters:
576785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
577785d1243SStefano Zampini 
578785d1243SStefano Zampini    Level: intermediate
579785d1243SStefano Zampini 
5800f202f7eSStefano Zampini    Notes:
5810f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
582785d1243SStefano Zampini 
583785d1243SStefano Zampini .seealso: PCBDDC
584785d1243SStefano Zampini @*/
585785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
586785d1243SStefano Zampini {
587785d1243SStefano Zampini   PetscErrorCode ierr;
588785d1243SStefano Zampini 
589785d1243SStefano Zampini   PetscFunctionBegin;
590785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
591785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
592785d1243SStefano Zampini   PetscFunctionReturn(0);
593785d1243SStefano Zampini }
594785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
595785d1243SStefano Zampini 
596785d1243SStefano Zampini #undef __FUNCT__
597785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
598785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
599785d1243SStefano Zampini {
600785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
601785d1243SStefano Zampini 
602785d1243SStefano Zampini   PetscFunctionBegin;
603785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
604785d1243SStefano Zampini   PetscFunctionReturn(0);
605785d1243SStefano Zampini }
606785d1243SStefano Zampini 
607785d1243SStefano Zampini #undef __FUNCT__
60882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
609da1bb401SStefano Zampini /*@
61082ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
611da1bb401SStefano Zampini 
612785d1243SStefano Zampini    Collective
613da1bb401SStefano Zampini 
614da1bb401SStefano Zampini    Input Parameters:
61528509bceSStefano Zampini .  pc - the preconditioning context
616da1bb401SStefano Zampini 
617da1bb401SStefano Zampini    Output Parameters:
61828509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
619da1bb401SStefano Zampini 
620da1bb401SStefano Zampini    Level: intermediate
621da1bb401SStefano Zampini 
622da1bb401SStefano Zampini    Notes:
6230f202f7eSStefano Zampini      The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries).
6240f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
625da1bb401SStefano Zampini 
626da1bb401SStefano Zampini .seealso: PCBDDC
627da1bb401SStefano Zampini @*/
62882ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
629da1bb401SStefano Zampini {
630da1bb401SStefano Zampini   PetscErrorCode ierr;
631da1bb401SStefano Zampini 
632da1bb401SStefano Zampini   PetscFunctionBegin;
633da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
63482ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
635da1bb401SStefano Zampini   PetscFunctionReturn(0);
636da1bb401SStefano Zampini }
637da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6381e6b0712SBarry Smith 
639da1bb401SStefano Zampini #undef __FUNCT__
64053cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
64153cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
64253cdbc3dSStefano Zampini {
64353cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
64453cdbc3dSStefano Zampini 
64553cdbc3dSStefano Zampini   PetscFunctionBegin;
64653cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
64753cdbc3dSStefano Zampini   PetscFunctionReturn(0);
64853cdbc3dSStefano Zampini }
6491e6b0712SBarry Smith 
65053cdbc3dSStefano Zampini #undef __FUNCT__
65153cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
65253cdbc3dSStefano Zampini /*@
653785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
65453cdbc3dSStefano Zampini 
655785d1243SStefano Zampini    Collective
656785d1243SStefano Zampini 
657785d1243SStefano Zampini    Input Parameters:
658785d1243SStefano Zampini .  pc - the preconditioning context
659785d1243SStefano Zampini 
660785d1243SStefano Zampini    Output Parameters:
661785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
662785d1243SStefano Zampini 
663785d1243SStefano Zampini    Level: intermediate
664785d1243SStefano Zampini 
6650f202f7eSStefano Zampini    Notes:
6660f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
667785d1243SStefano Zampini 
668785d1243SStefano Zampini .seealso: PCBDDC
669785d1243SStefano Zampini @*/
670785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
671785d1243SStefano Zampini {
672785d1243SStefano Zampini   PetscErrorCode ierr;
673785d1243SStefano Zampini 
674785d1243SStefano Zampini   PetscFunctionBegin;
675785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
676785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
677785d1243SStefano Zampini   PetscFunctionReturn(0);
678785d1243SStefano Zampini }
679785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
680785d1243SStefano Zampini 
681785d1243SStefano Zampini #undef __FUNCT__
682785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
683785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
684785d1243SStefano Zampini {
685785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
686785d1243SStefano Zampini 
687785d1243SStefano Zampini   PetscFunctionBegin;
688785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
689785d1243SStefano Zampini   PetscFunctionReturn(0);
690785d1243SStefano Zampini }
691785d1243SStefano Zampini 
692785d1243SStefano Zampini #undef __FUNCT__
69382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
69453cdbc3dSStefano Zampini /*@
69582ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
69653cdbc3dSStefano Zampini 
697785d1243SStefano Zampini    Collective
69853cdbc3dSStefano Zampini 
69953cdbc3dSStefano Zampini    Input Parameters:
70028509bceSStefano Zampini .  pc - the preconditioning context
70153cdbc3dSStefano Zampini 
70253cdbc3dSStefano Zampini    Output Parameters:
70328509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
70453cdbc3dSStefano Zampini 
70553cdbc3dSStefano Zampini    Level: intermediate
70653cdbc3dSStefano Zampini 
70753cdbc3dSStefano Zampini    Notes:
7080f202f7eSStefano Zampini      The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries).
7090f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
71053cdbc3dSStefano Zampini 
71153cdbc3dSStefano Zampini .seealso: PCBDDC
71253cdbc3dSStefano Zampini @*/
71382ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
71453cdbc3dSStefano Zampini {
71553cdbc3dSStefano Zampini   PetscErrorCode ierr;
71653cdbc3dSStefano Zampini 
71753cdbc3dSStefano Zampini   PetscFunctionBegin;
71853cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
71982ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7200c7d97c5SJed Brown   PetscFunctionReturn(0);
7210c7d97c5SJed Brown }
72236e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7231e6b0712SBarry Smith 
72436e030ebSStefano Zampini #undef __FUNCT__
725da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7261a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
72736e030ebSStefano Zampini {
72836e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
729da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
730da1bb401SStefano Zampini   PetscErrorCode ierr;
73136e030ebSStefano Zampini 
73236e030ebSStefano Zampini   PetscFunctionBegin;
733674ae819SStefano Zampini   /* free old CSR */
734674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
735fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
736674ae819SStefano Zampini   /* get CSR into graph structure */
737da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
738854ce69bSBarry Smith     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
739785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
740674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
741674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
742da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7431a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7441a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
745674ae819SStefano Zampini   }
746575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
74736e030ebSStefano Zampini   PetscFunctionReturn(0);
74836e030ebSStefano Zampini }
7491e6b0712SBarry Smith 
75036e030ebSStefano Zampini #undef __FUNCT__
751da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
75236e030ebSStefano Zampini /*@
7530f202f7eSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
75436e030ebSStefano Zampini 
75536e030ebSStefano Zampini    Not collective
75636e030ebSStefano Zampini 
75736e030ebSStefano Zampini    Input Parameters:
75836e030ebSStefano Zampini +  pc - the preconditioning context
7590f202f7eSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
76028509bceSStefano Zampini .  xadj, adjncy - the CSR graph
76128509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
76236e030ebSStefano Zampini 
76336e030ebSStefano Zampini    Level: intermediate
76436e030ebSStefano Zampini 
76536e030ebSStefano Zampini    Notes:
76636e030ebSStefano Zampini 
76728509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
76836e030ebSStefano Zampini @*/
7691a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
77036e030ebSStefano Zampini {
771575ad6abSStefano Zampini   void (*f)(void) = 0;
77236e030ebSStefano Zampini   PetscErrorCode ierr;
77336e030ebSStefano Zampini 
77436e030ebSStefano Zampini   PetscFunctionBegin;
77536e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
776674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
777d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
7786c4ed002SBarry Smith   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER)  SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d",copymode);
7791a83f524SJed Brown   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
780575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
781575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
782575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
783575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
784575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
785da1bb401SStefano Zampini   }
78636e030ebSStefano Zampini   PetscFunctionReturn(0);
78736e030ebSStefano Zampini }
7889c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7891e6b0712SBarry Smith 
7909c0446d6SStefano Zampini #undef __FUNCT__
79163602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
79263602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
79363602bcaSStefano Zampini {
79463602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
79563602bcaSStefano Zampini   PetscInt i;
79663602bcaSStefano Zampini   PetscErrorCode ierr;
79763602bcaSStefano Zampini 
79863602bcaSStefano Zampini   PetscFunctionBegin;
79963602bcaSStefano Zampini   /* Destroy ISes if they were already set */
80063602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
80163602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
80263602bcaSStefano Zampini   }
80363602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
80463602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
80563602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
80663602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
80763602bcaSStefano Zampini   }
80863602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
80963602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
81063602bcaSStefano Zampini   /* allocate space then set */
811d02579f5SStefano Zampini   if (n_is) {
812d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
813d02579f5SStefano Zampini   }
81463602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
81563602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
81663602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
81763602bcaSStefano Zampini   }
81863602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
81963602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8201a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
82163602bcaSStefano Zampini   PetscFunctionReturn(0);
82263602bcaSStefano Zampini }
82363602bcaSStefano Zampini 
82463602bcaSStefano Zampini #undef __FUNCT__
82563602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
82663602bcaSStefano Zampini /*@
82763602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
82863602bcaSStefano Zampini 
82963602bcaSStefano Zampini    Collective
83063602bcaSStefano Zampini 
83163602bcaSStefano Zampini    Input Parameters:
83263602bcaSStefano Zampini +  pc - the preconditioning context
8330f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
8340f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in local ordering
83563602bcaSStefano Zampini 
83663602bcaSStefano Zampini    Level: intermediate
83763602bcaSStefano Zampini 
8380f202f7eSStefano Zampini    Notes:
8390f202f7eSStefano Zampini      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
84063602bcaSStefano Zampini 
84163602bcaSStefano Zampini .seealso: PCBDDC
84263602bcaSStefano Zampini @*/
84363602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
84463602bcaSStefano Zampini {
84563602bcaSStefano Zampini   PetscInt       i;
84663602bcaSStefano Zampini   PetscErrorCode ierr;
84763602bcaSStefano Zampini 
84863602bcaSStefano Zampini   PetscFunctionBegin;
84963602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
85063602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
85163602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
85263602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
85363602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
85463602bcaSStefano Zampini   }
855e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
85663602bcaSStefano Zampini   PetscFunctionReturn(0);
85763602bcaSStefano Zampini }
85863602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
85963602bcaSStefano Zampini 
86063602bcaSStefano Zampini #undef __FUNCT__
8619c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8629c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8639c0446d6SStefano Zampini {
8649c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8659c0446d6SStefano Zampini   PetscInt i;
8669c0446d6SStefano Zampini   PetscErrorCode ierr;
8679c0446d6SStefano Zampini 
8689c0446d6SStefano Zampini   PetscFunctionBegin;
869da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8709c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8719c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8729c0446d6SStefano Zampini   }
873d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
87463602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
87563602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
87663602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
87763602bcaSStefano Zampini   }
87863602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
87963602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
880da1bb401SStefano Zampini   /* allocate space then set */
881d02579f5SStefano Zampini   if (n_is) {
882785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
883d02579f5SStefano Zampini   }
8849c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
885da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
886da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8879c0446d6SStefano Zampini   }
8889c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
88963602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8901a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8919c0446d6SStefano Zampini   PetscFunctionReturn(0);
8929c0446d6SStefano Zampini }
8931e6b0712SBarry Smith 
8949c0446d6SStefano Zampini #undef __FUNCT__
8959c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8969c0446d6SStefano Zampini /*@
89763602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8989c0446d6SStefano Zampini 
89963602bcaSStefano Zampini    Collective
9009c0446d6SStefano Zampini 
9019c0446d6SStefano Zampini    Input Parameters:
9029c0446d6SStefano Zampini +  pc - the preconditioning context
9030f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
9040f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in global ordering
9059c0446d6SStefano Zampini 
9069c0446d6SStefano Zampini    Level: intermediate
9079c0446d6SStefano Zampini 
9080f202f7eSStefano Zampini    Notes:
9090f202f7eSStefano Zampini      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
9109c0446d6SStefano Zampini 
9119c0446d6SStefano Zampini .seealso: PCBDDC
9129c0446d6SStefano Zampini @*/
9139c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9149c0446d6SStefano Zampini {
9152b510759SStefano Zampini   PetscInt       i;
9169c0446d6SStefano Zampini   PetscErrorCode ierr;
9179c0446d6SStefano Zampini 
9189c0446d6SStefano Zampini   PetscFunctionBegin;
9199c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
92063602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9212b510759SStefano Zampini   for (i=0;i<n_is;i++) {
92263602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
92363602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9242b510759SStefano Zampini   }
9259c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9269c0446d6SStefano Zampini   PetscFunctionReturn(0);
9279c0446d6SStefano Zampini }
928906d46d4SStefano Zampini 
929da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
930534831adSStefano Zampini #undef __FUNCT__
931534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
932534831adSStefano Zampini /* -------------------------------------------------------------------------- */
933534831adSStefano Zampini /*
934534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
935534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9369c0446d6SStefano Zampini 
937534831adSStefano Zampini    Input Parameter:
938534831adSStefano Zampini +  pc - the preconditioner contex
939534831adSStefano Zampini 
940534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
941534831adSStefano Zampini 
942534831adSStefano Zampini    Notes:
943534831adSStefano Zampini      The interface routine PCPreSolve() is not usually called directly by
944534831adSStefano Zampini    the user, but instead is called by KSPSolve().
945534831adSStefano Zampini */
946534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
947534831adSStefano Zampini {
948534831adSStefano Zampini   PetscErrorCode ierr;
949534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
950534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9513972b0daSStefano Zampini   Vec            used_vec;
9528d00608fSStefano Zampini   PetscBool      copy_rhs = PETSC_TRUE;
953534831adSStefano Zampini 
954534831adSStefano Zampini   PetscFunctionBegin;
95585c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
95685c4d303SStefano Zampini   if (ksp) {
95785c4d303SStefano Zampini     PetscBool iscg;
95885c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
95985c4d303SStefano Zampini     if (!iscg) {
96085c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
96185c4d303SStefano Zampini     }
96285c4d303SStefano Zampini   }
963c7017625SStefano Zampini 
96485c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
96562a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
96662a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
96762a6ff1dSStefano Zampini   }
96862a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
96962a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
97062a6ff1dSStefano Zampini   }
9718d00608fSStefano Zampini 
9723972b0daSStefano Zampini   if (x) {
9733972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9743972b0daSStefano Zampini     used_vec = x;
9758d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
9763972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9773972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9783972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9793972b0daSStefano Zampini   }
9808efcfb23SStefano Zampini 
9818efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
9823972b0daSStefano Zampini   if (ksp) {
983a0cb1b98SStefano Zampini     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
9848efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
9858efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
9863972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9873972b0daSStefano Zampini     }
9883972b0daSStefano Zampini   }
9893308cffdSStefano Zampini 
9908d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
9913972b0daSStefano Zampini 
9923972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
993a07ea27aSStefano Zampini   if (rhs) {
9943975b054SStefano Zampini     IS dirIS = NULL;
9953975b054SStefano Zampini 
996a07ea27aSStefano Zampini     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
9973975b054SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
9983975b054SStefano Zampini     if (dirIS) {
999906d46d4SStefano Zampini       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
1000785d1243SStefano Zampini       PetscInt          dirsize,i,*is_indices;
10012b095fd8SStefano Zampini       PetscScalar       *array_x;
10022b095fd8SStefano Zampini       const PetscScalar *array_diagonal;
1003785d1243SStefano Zampini 
10043972b0daSStefano Zampini       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
10053972b0daSStefano Zampini       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1006e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1007e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1008e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1009e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
101082ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
10113972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10122b095fd8SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10133972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10142fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
10153972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10162b095fd8SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10173972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1018e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1019e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10208d00608fSStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
10211b968477SStefano Zampini       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
10228efcfb23SStefano Zampini     }
1023a07ea27aSStefano Zampini   }
1024b76ba322SStefano Zampini 
10258efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
10268d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
10278d00608fSStefano Zampini     /* store the original rhs */
10288d00608fSStefano Zampini     if (copy_rhs) {
10298d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10308d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10318d00608fSStefano Zampini     }
10328d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10333972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10343972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10353972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10368efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10377acc28cbSStefano Zampini     if (ksp) {
10387acc28cbSStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
10397acc28cbSStefano Zampini     }
10403308cffdSStefano Zampini   }
10418efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1042b76ba322SStefano Zampini 
1043b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
1044b097fa66SStefano Zampini   if (x && pcbddc->use_exact_dirichlet_trick) {
10458efcfb23SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1046b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1047b76ba322SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1048b76ba322SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
10498efcfb23SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10508efcfb23SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1051a0cb1b98SStefano Zampini     if (ksp) {
1052b76ba322SStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1053b76ba322SStefano Zampini     }
1054b76ba322SStefano Zampini   }
1055b76ba322SStefano Zampini 
1056b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1057906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1058906d46d4SStefano Zampini 
1059906d46d4SStefano Zampini     /* get change ctx */
1060906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1061906d46d4SStefano Zampini 
1062906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1063906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1064906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1065906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1066906d46d4SStefano Zampini 
1067906d46d4SStefano Zampini     /* change iteration matrix */
1068906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1069906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1070906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1071906d46d4SStefano Zampini 
10728d00608fSStefano Zampini     /* store the original rhs */
10738d00608fSStefano Zampini     if (copy_rhs) {
10748d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10758d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10768d00608fSStefano Zampini     }
10778d00608fSStefano Zampini 
1078906d46d4SStefano Zampini     /* change rhs */
1079906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1080906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
10818d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
108292e3dcfbSStefano Zampini   }
1083534831adSStefano Zampini   PetscFunctionReturn(0);
1084534831adSStefano Zampini }
1085906d46d4SStefano Zampini 
1086534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1087534831adSStefano Zampini #undef __FUNCT__
1088534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1089534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1090534831adSStefano Zampini /*
1091534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1092534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1093534831adSStefano Zampini 
1094534831adSStefano Zampini    Input Parameter:
1095534831adSStefano Zampini +  pc - the preconditioner contex
1096534831adSStefano Zampini 
1097534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1098534831adSStefano Zampini 
1099534831adSStefano Zampini    Notes:
1100534831adSStefano Zampini      The interface routine PCPostSolve() is not usually called directly by
1101534831adSStefano Zampini      the user, but instead is called by KSPSolve().
1102534831adSStefano Zampini */
1103534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1104534831adSStefano Zampini {
1105534831adSStefano Zampini   PetscErrorCode ierr;
1106534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1107534831adSStefano Zampini 
1108534831adSStefano Zampini   PetscFunctionBegin;
1109b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1110906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1111906d46d4SStefano Zampini 
1112906d46d4SStefano Zampini     /* get change ctx */
1113906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1114906d46d4SStefano Zampini 
1115906d46d4SStefano Zampini     /* restore iteration matrix */
1116906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1117906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1118906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1119906d46d4SStefano Zampini 
1120906d46d4SStefano Zampini     /* get solution in original basis */
1121906d46d4SStefano Zampini     if (x) {
1122906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1123906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1124906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
11253425bc38SStefano Zampini     }
1126534831adSStefano Zampini   }
1127906d46d4SStefano Zampini 
11283972b0daSStefano Zampini   /* add solution removed in presolve */
11296bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
11303425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
11313425bc38SStefano Zampini   }
1132906d46d4SStefano Zampini 
1133fb223d50SStefano Zampini   /* restore rhs to its original state */
11348d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
1135fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1136fb223d50SStefano Zampini   }
11378d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
11388efcfb23SStefano Zampini 
11398efcfb23SStefano Zampini   /* restore ksp guess state */
11408efcfb23SStefano Zampini   if (ksp) {
11418efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
11428efcfb23SStefano Zampini   }
1143534831adSStefano Zampini   PetscFunctionReturn(0);
1144534831adSStefano Zampini }
1145534831adSStefano Zampini /* -------------------------------------------------------------------------- */
114653cdbc3dSStefano Zampini #undef __FUNCT__
114753cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
11480c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
11490c7d97c5SJed Brown /*
11500c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
11510c7d97c5SJed Brown                   by setting data structures and options.
11520c7d97c5SJed Brown 
11530c7d97c5SJed Brown    Input Parameter:
115453cdbc3dSStefano Zampini +  pc - the preconditioner context
11550c7d97c5SJed Brown 
11560c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
11570c7d97c5SJed Brown 
11580c7d97c5SJed Brown    Notes:
11590c7d97c5SJed Brown      The interface routine PCSetUp() is not usually called directly by
11600c7d97c5SJed Brown      the user, but instead is called by PCApply() if necessary.
11610c7d97c5SJed Brown */
116253cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
11630c7d97c5SJed Brown {
11640c7d97c5SJed Brown   PetscErrorCode ierr;
11650c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
11665e8657edSStefano Zampini   Mat_IS*        matis;
116708122e43SStefano Zampini   MatNullSpace   nearnullspace;
116891e8d312SStefano Zampini   PetscInt       nrows,ncols;
116908122e43SStefano Zampini   PetscBool      computetopography,computesolvers,computesubschurs;
11708de1fae6SStefano Zampini   PetscBool      computeconstraintsmatrix;
11715e8657edSStefano Zampini   PetscBool      new_nearnullspace_provided,ismatis;
11720c7d97c5SJed Brown 
11730c7d97c5SJed Brown   PetscFunctionBegin;
11745e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
11756c4ed002SBarry Smith   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
117691e8d312SStefano Zampini   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
11776c4ed002SBarry Smith   if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
11785e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1179f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
11803b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1181fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1182f4ddd8eeSStefano Zampini   /* split work */
1183674ae819SStefano Zampini   if (pc->setupcalled) {
1184d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1185674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1186674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1187674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1188674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1189674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1190674ae819SStefano Zampini     }
1191674ae819SStefano Zampini   } else {
1192674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1193674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1194674ae819SStefano Zampini   }
1195fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1196fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1197fb180af4SStefano Zampini   }
11988de1fae6SStefano Zampini   computeconstraintsmatrix = PETSC_FALSE;
11996c4ed002SBarry Smith   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
12005a95e1ceSStefano Zampini   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1201862806e4SStefano Zampini   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1202862806e4SStefano Zampini 
12035a95e1ceSStefano Zampini   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
12046c4ed002SBarry Smith   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute faster deluxe if adaptivity and change of basis are both requested. Rerun with -pc_bddc_deluxe_faster false");
12056c4ed002SBarry Smith 
1206f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
120770cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
120870cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
120958a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
12101575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);CHKERRQ(ierr);
1211f4ddd8eeSStefano Zampini     }
121258a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1213f4ddd8eeSStefano Zampini   }
1214f4ddd8eeSStefano Zampini 
12155e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
12165e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
12175e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
12185e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
12195e8657edSStefano Zampini   } else {
1220b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
12215e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
12225e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
1223674ae819SStefano Zampini   }
1224f4ddd8eeSStefano Zampini 
1225e496cd5dSStefano Zampini   /* workaround for reals */
1226e496cd5dSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
12273301b35fSStefano Zampini   if (matis->A->symmetric_set) {
12283301b35fSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
12293301b35fSStefano Zampini   }
1230e496cd5dSStefano Zampini #endif
1231e496cd5dSStefano Zampini 
12325e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
12335e8657edSStefano Zampini   {
12345e8657edSStefano Zampini     Mat temp_mat;
12355e8657edSStefano Zampini 
12365e8657edSStefano Zampini     temp_mat = matis->A;
12375e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
12385e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
12395e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
12405e8657edSStefano Zampini     matis->A = temp_mat;
12415e8657edSStefano Zampini   }
1242684f6988SStefano Zampini 
1243b96c3477SStefano Zampini   /* Analyze interface and setup sub_schurs data */
1244674ae819SStefano Zampini   if (computetopography) {
1245674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
12468de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1247674ae819SStefano Zampini   }
1248fb8d54d4SStefano Zampini 
1249b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1250684f6988SStefano Zampini   if (computesolvers) {
1251d5574798SStefano Zampini     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1252d5574798SStefano Zampini 
1253d5574798SStefano Zampini     if (computesubschurs && computetopography) {
125408122e43SStefano Zampini       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1255b1b3d7a2SStefano Zampini     }
1256d5574798SStefano Zampini     if (sub_schurs->use_mumps) {
12572070dbb6SStefano Zampini       if (computesubschurs) {
125808122e43SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
12592070dbb6SStefano Zampini       }
1260d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1261d5574798SStefano Zampini     } else {
1262d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
12632070dbb6SStefano Zampini       if (computesubschurs) {
1264d5574798SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1265d5574798SStefano Zampini       }
12662070dbb6SStefano Zampini     }
126708122e43SStefano Zampini     if (pcbddc->adaptive_selection) {
126808122e43SStefano Zampini       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
12698de1fae6SStefano Zampini       computeconstraintsmatrix = PETSC_TRUE;
1270b7eb3628SStefano Zampini     }
1271b1b3d7a2SStefano Zampini   }
1272684f6988SStefano Zampini 
1273f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1274fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1275f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1276f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1277f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1278f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1279f4ddd8eeSStefano Zampini     } else {
1280f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1281f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1282f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1283165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1284f4ddd8eeSStefano Zampini         PetscInt         i;
1285165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1286165b64e2SStefano Zampini         PetscObjectState state;
1287165b64e2SStefano Zampini         PetscInt         nnsp_size;
1288165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1289f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1290f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1291165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1292f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1293f4ddd8eeSStefano Zampini             break;
1294f4ddd8eeSStefano Zampini           }
1295f4ddd8eeSStefano Zampini         }
1296f4ddd8eeSStefano Zampini       }
1297f4ddd8eeSStefano Zampini     }
1298f4ddd8eeSStefano Zampini   } else {
1299f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1300f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1301f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1302f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1303f4ddd8eeSStefano Zampini     }
1304f4ddd8eeSStefano Zampini   }
1305f4ddd8eeSStefano Zampini 
1306f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1307727cdba6SStefano Zampini   /* reset primal space flags */
1308f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1309727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
13108de1fae6SStefano Zampini   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1311727cdba6SStefano Zampini     /* It also sets the primal space flags */
1312674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1313e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1314f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
13153975b054SStefano Zampini   }
13165e8657edSStefano Zampini 
13173975b054SStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
13185e8657edSStefano Zampini     if (pcbddc->use_change_of_basis) {
13195e8657edSStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
13205e8657edSStefano Zampini 
13215e8657edSStefano Zampini       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
13225e8657edSStefano Zampini       /* get submatrices */
13235e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
13245e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
13255e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
13265e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
13275e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
13285e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
13293975b054SStefano Zampini       /* set flag in pcis to not reuse submatrices during PCISCreate */
13303975b054SStefano Zampini       pcis->reusesubmatrices = PETSC_FALSE;
13315e8657edSStefano Zampini     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1332b96c3477SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
13335e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
13345e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
13355e8657edSStefano Zampini     }
1336b96c3477SStefano Zampini     /* SetUp coarse and local Neumann solvers */
133799cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1338b96c3477SStefano Zampini     /* SetUp Scaling operator */
1339674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
13400c7d97c5SJed Brown   }
134170cf5478SStefano Zampini 
134258a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
134358a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
13442b510759SStefano Zampini   }
13450c7d97c5SJed Brown   PetscFunctionReturn(0);
13460c7d97c5SJed Brown }
13470c7d97c5SJed Brown 
13480c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
13490c7d97c5SJed Brown /*
135050efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
13510c7d97c5SJed Brown 
13520c7d97c5SJed Brown    Input Parameters:
13530f202f7eSStefano Zampini +  pc - the preconditioner context
13540f202f7eSStefano Zampini -  r - input vector (global)
13550c7d97c5SJed Brown 
13560c7d97c5SJed Brown    Output Parameter:
13570c7d97c5SJed Brown .  z - output vector (global)
13580c7d97c5SJed Brown 
13590c7d97c5SJed Brown    Application Interface Routine: PCApply()
13600c7d97c5SJed Brown  */
13610c7d97c5SJed Brown #undef __FUNCT__
13620c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
136353cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
13640c7d97c5SJed Brown {
13650c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
13660c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1367b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
13680c7d97c5SJed Brown   PetscErrorCode    ierr;
13693b03a366Sstefano_zampini   const PetscScalar one = 1.0;
13703b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
13712617d88aSStefano Zampini   const PetscScalar zero = 0.0;
13720c7d97c5SJed Brown 
13730c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
13740c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
1375b097fa66SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
13760c7d97c5SJed Brown 
13770c7d97c5SJed Brown   PetscFunctionBegin;
137885c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1379b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
13800c7d97c5SJed Brown     /* First Dirichlet solve */
13810c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13820c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13830c7d97c5SJed Brown     /*
13840c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1385b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1386674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
13870c7d97c5SJed Brown     */
1388b097fa66SStefano Zampini     if (n_D) {
1389b097fa66SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
13900c7d97c5SJed Brown       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
13918eeda7d8SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1392b097fa66SStefano Zampini       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1393b097fa66SStefano Zampini     } else {
1394b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1395b097fa66SStefano Zampini     }
13960c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13970c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1398674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1399b76ba322SStefano Zampini   } else {
1400b097fa66SStefano Zampini     if (pcbddc->switch_static) {
14010bdf917eSStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1402b097fa66SStefano Zampini     }
1403674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1404b76ba322SStefano Zampini   }
1405b76ba322SStefano Zampini 
14062617d88aSStefano Zampini   /* Apply interface preconditioner
14072617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1408dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
14092617d88aSStefano Zampini 
1410674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1411674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
14120c7d97c5SJed Brown 
14133b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
14140c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14150c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1416b097fa66SStefano Zampini   if (n_B) {
14170c7d97c5SJed Brown     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
14188eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1419b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1420b097fa66SStefano Zampini     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1421b097fa66SStefano Zampini   }
1422df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1423b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1424b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1425b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1426b097fa66SStefano Zampini     } else {
1427b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1428b097fa66SStefano Zampini     }
14290c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14300c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1431b097fa66SStefano Zampini   } else {
1432b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1433b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1434b097fa66SStefano Zampini     } else {
1435b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1436b097fa66SStefano Zampini     }
1437b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1438b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1439b097fa66SStefano Zampini   }
14400c7d97c5SJed Brown   PetscFunctionReturn(0);
14410c7d97c5SJed Brown }
144250efa1b5SStefano Zampini 
144350efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
144450efa1b5SStefano Zampini /*
144550efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
144650efa1b5SStefano Zampini 
144750efa1b5SStefano Zampini    Input Parameters:
14480f202f7eSStefano Zampini +  pc - the preconditioner context
14490f202f7eSStefano Zampini -  r - input vector (global)
145050efa1b5SStefano Zampini 
145150efa1b5SStefano Zampini    Output Parameter:
145250efa1b5SStefano Zampini .  z - output vector (global)
145350efa1b5SStefano Zampini 
145450efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
145550efa1b5SStefano Zampini  */
145650efa1b5SStefano Zampini #undef __FUNCT__
145750efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
145850efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
145950efa1b5SStefano Zampini {
146050efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
146150efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1462b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
146350efa1b5SStefano Zampini   PetscErrorCode    ierr;
146450efa1b5SStefano Zampini   const PetscScalar one = 1.0;
146550efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
146650efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
146750efa1b5SStefano Zampini 
146850efa1b5SStefano Zampini   PetscFunctionBegin;
146950efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1470b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
147150efa1b5SStefano Zampini     /* First Dirichlet solve */
147250efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
147350efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
147450efa1b5SStefano Zampini     /*
147550efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
1476b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
147750efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
147850efa1b5SStefano Zampini     */
1479b097fa66SStefano Zampini     if (n_D) {
1480b097fa66SStefano Zampini       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
148150efa1b5SStefano Zampini       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
148250efa1b5SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1483b097fa66SStefano Zampini       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1484b097fa66SStefano Zampini     } else {
1485b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1486b097fa66SStefano Zampini     }
148750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
148850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
148950efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
149050efa1b5SStefano Zampini   } else {
1491b097fa66SStefano Zampini     if (pcbddc->switch_static) {
149250efa1b5SStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1493b097fa66SStefano Zampini     }
149450efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
149550efa1b5SStefano Zampini   }
149650efa1b5SStefano Zampini 
149750efa1b5SStefano Zampini   /* Apply interface preconditioner
149850efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1499dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
150050efa1b5SStefano Zampini 
150150efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
150250efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
150350efa1b5SStefano Zampini 
150450efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
150550efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
150650efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1507b097fa66SStefano Zampini   if (n_B) {
150850efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
150950efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1510b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1511b097fa66SStefano Zampini     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1512b097fa66SStefano Zampini   }
1513b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1514b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1515b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1516b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1517b097fa66SStefano Zampini     } else {
1518b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1519b097fa66SStefano Zampini     }
152050efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
152150efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1522b097fa66SStefano Zampini   } else {
1523b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1524b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1525b097fa66SStefano Zampini     } else {
1526b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1527b097fa66SStefano Zampini     }
1528b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1529b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1530b097fa66SStefano Zampini   }
153150efa1b5SStefano Zampini   PetscFunctionReturn(0);
153250efa1b5SStefano Zampini }
1533da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1534674ae819SStefano Zampini 
1535da1bb401SStefano Zampini #undef __FUNCT__
1536da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1537da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1538da1bb401SStefano Zampini {
1539da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1540da1bb401SStefano Zampini   PetscErrorCode ierr;
1541da1bb401SStefano Zampini 
1542da1bb401SStefano Zampini   PetscFunctionBegin;
1543da1bb401SStefano Zampini   /* free data created by PCIS */
1544da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1545674ae819SStefano Zampini   /* free BDDC custom data  */
1546674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1547674ae819SStefano Zampini   /* destroy objects related to topography */
1548674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1549674ae819SStefano Zampini   /* free allocated graph structure */
1550da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1551b96c3477SStefano Zampini   /* free allocated sub schurs structure */
1552b96c3477SStefano Zampini   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
155334a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1554674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
155534a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1556674ae819SStefano Zampini   /* free solvers stuff */
1557674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
155862a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
155962a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
156062a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1561906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1562906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1563906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1564906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1565906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1566906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1567906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1568906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1569906d46d4SStefano Zampini   }
1570906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
15713425bc38SStefano Zampini   /* remove functions */
1572906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1573674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1574bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
15752b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1576b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
15772b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1578bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
157982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1580bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
158182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1582bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
158382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1584bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1585785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1586bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
158763602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1588bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1589bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1590bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1591bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1592a06fd7f2SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);CHKERRQ(ierr);
1593674ae819SStefano Zampini   /* Free the private data structure */
1594674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1595da1bb401SStefano Zampini   PetscFunctionReturn(0);
1596da1bb401SStefano Zampini }
15973425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
15981e6b0712SBarry Smith 
15993425bc38SStefano Zampini #undef __FUNCT__
1600a06fd7f2SStefano Zampini #define __FUNCT__ "PCPreSolveChangeRHS_BDDC"
1601a06fd7f2SStefano Zampini static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool* change)
1602a06fd7f2SStefano Zampini {
1603a06fd7f2SStefano Zampini   PetscFunctionBegin;
1604a06fd7f2SStefano Zampini   *change = PETSC_TRUE;
1605a06fd7f2SStefano Zampini   PetscFunctionReturn(0);
1606a06fd7f2SStefano Zampini }
1607a06fd7f2SStefano Zampini 
1608a06fd7f2SStefano Zampini #undef __FUNCT__
16093425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
16103425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
16113425bc38SStefano Zampini {
1612674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
1613c08af4c6SStefano Zampini   Vec            copy_standard_rhs;
16143425bc38SStefano Zampini   PC_IS*         pcis;
16153425bc38SStefano Zampini   PC_BDDC*       pcbddc;
16163425bc38SStefano Zampini   PetscErrorCode ierr;
16170c7d97c5SJed Brown 
16183425bc38SStefano Zampini   PetscFunctionBegin;
16193425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16203425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
16213425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
16223425bc38SStefano Zampini 
1623c08af4c6SStefano Zampini   /*
1624c08af4c6SStefano Zampini      change of basis for physical rhs if needed
1625c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
1626c08af4c6SStefano Zampini      TODO: better management when FETIDP will have its own class
1627c08af4c6SStefano Zampini   */
1628c08af4c6SStefano Zampini   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1629c08af4c6SStefano Zampini   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1630c08af4c6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
16313425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
1632c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1633c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1634fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1635fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
1636c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1637c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1638674ae819SStefano Zampini   /* Apply partition of unity */
16393425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1640c08af4c6SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
16418eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
16423425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
16433425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16443425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
16453425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1646c08af4c6SStefano Zampini     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1647c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1648c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1649c08af4c6SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1650c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1651c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16523425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
16533425bc38SStefano Zampini   }
1654c08af4c6SStefano Zampini   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
16553425bc38SStefano Zampini   /* BDDC rhs */
16563425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
16578eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
16583425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16593425bc38SStefano Zampini   }
16603425bc38SStefano Zampini   /* apply BDDC */
1661dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
16623425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
16633425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
16643425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
16653425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16663425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16673425bc38SStefano Zampini   PetscFunctionReturn(0);
16683425bc38SStefano Zampini }
16691e6b0712SBarry Smith 
16703425bc38SStefano Zampini #undef __FUNCT__
16713425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
16723425bc38SStefano Zampini /*@
16730f202f7eSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
16743425bc38SStefano Zampini 
16753425bc38SStefano Zampini    Collective
16763425bc38SStefano Zampini 
16773425bc38SStefano Zampini    Input Parameters:
16780f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
16790f202f7eSStefano Zampini -  standard_rhs    - the right-hand side of the original linear system
16803425bc38SStefano Zampini 
16813425bc38SStefano Zampini    Output Parameters:
16820f202f7eSStefano Zampini .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
16833425bc38SStefano Zampini 
16843425bc38SStefano Zampini    Level: developer
16853425bc38SStefano Zampini 
16863425bc38SStefano Zampini    Notes:
16873425bc38SStefano Zampini 
16880f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
16893425bc38SStefano Zampini @*/
16903425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
16913425bc38SStefano Zampini {
1692674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
16933425bc38SStefano Zampini   PetscErrorCode ierr;
16943425bc38SStefano Zampini 
16953425bc38SStefano Zampini   PetscFunctionBegin;
16963425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1697163d334eSBarry Smith   ierr = PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
16983425bc38SStefano Zampini   PetscFunctionReturn(0);
16993425bc38SStefano Zampini }
17003425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
17011e6b0712SBarry Smith 
17023425bc38SStefano Zampini #undef __FUNCT__
17033425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
17043425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
17053425bc38SStefano Zampini {
1706674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17073425bc38SStefano Zampini   PC_IS*         pcis;
17083425bc38SStefano Zampini   PC_BDDC*       pcbddc;
17093425bc38SStefano Zampini   PetscErrorCode ierr;
17103425bc38SStefano Zampini 
17113425bc38SStefano Zampini   PetscFunctionBegin;
17123425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
17133425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
17143425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
17153425bc38SStefano Zampini 
17163425bc38SStefano Zampini   /* apply B_delta^T */
17173425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17183425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17193425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
17203425bc38SStefano Zampini   /* compute rhs for BDDC application */
17213425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
17228eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
17233425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17243425bc38SStefano Zampini   }
17253425bc38SStefano Zampini   /* apply BDDC */
1726dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
17273425bc38SStefano Zampini   /* put values into standard global vector */
17283425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17293425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17308eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
17313425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
17323425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
17333425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
17343425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17353425bc38SStefano Zampini   }
17363425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17373425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17383425bc38SStefano Zampini   /* final change of basis if needed
17393425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
17403308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
17413425bc38SStefano Zampini   PetscFunctionReturn(0);
17423425bc38SStefano Zampini }
17431e6b0712SBarry Smith 
17443425bc38SStefano Zampini #undef __FUNCT__
17453425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
17463425bc38SStefano Zampini /*@
17470f202f7eSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
17483425bc38SStefano Zampini 
17493425bc38SStefano Zampini    Collective
17503425bc38SStefano Zampini 
17513425bc38SStefano Zampini    Input Parameters:
17520f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
17530f202f7eSStefano Zampini -  fetidp_flux_sol - the solution of the FETI-DP linear system
17543425bc38SStefano Zampini 
17553425bc38SStefano Zampini    Output Parameters:
17560f202f7eSStefano Zampini .  standard_sol    - the solution defined on the physical domain
17573425bc38SStefano Zampini 
17583425bc38SStefano Zampini    Level: developer
17593425bc38SStefano Zampini 
17603425bc38SStefano Zampini    Notes:
17613425bc38SStefano Zampini 
17620f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
17633425bc38SStefano Zampini @*/
17643425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
17653425bc38SStefano Zampini {
1766674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17673425bc38SStefano Zampini   PetscErrorCode ierr;
17683425bc38SStefano Zampini 
17693425bc38SStefano Zampini   PetscFunctionBegin;
17703425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1771163d334eSBarry Smith   ierr = PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
17723425bc38SStefano Zampini   PetscFunctionReturn(0);
17733425bc38SStefano Zampini }
17743425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
17751e6b0712SBarry Smith 
1776f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1777edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1778f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1779f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1780edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1781f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1782674ae819SStefano Zampini 
17833425bc38SStefano Zampini #undef __FUNCT__
17843425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
1785*1720468bSStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, Mat *fetidp_mat, PC *fetidp_pc)
17863425bc38SStefano Zampini {
1787674ae819SStefano Zampini 
1788674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
17893425bc38SStefano Zampini   Mat            newmat;
1790674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
17913425bc38SStefano Zampini   PC             newpc;
1792ce94432eSBarry Smith   MPI_Comm       comm;
17933425bc38SStefano Zampini   PetscErrorCode ierr;
17943425bc38SStefano Zampini 
17953425bc38SStefano Zampini   PetscFunctionBegin;
1796ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
17973425bc38SStefano Zampini   /* FETIDP linear matrix */
17983425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
1799*1720468bSStefano Zampini   fetidpmat_ctx->fully_redundant = fully_redundant;
18003425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
18013425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
18023425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1803edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
18043425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
18053425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
18063425bc38SStefano Zampini   /* FETIDP preconditioner */
18073425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
18083425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
18093425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
18103425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
18113425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
18123425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1813edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
18143425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
181523ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
18163425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
18173425bc38SStefano Zampini   /* return pointers for objects created */
18183425bc38SStefano Zampini   *fetidp_mat=newmat;
18193425bc38SStefano Zampini   *fetidp_pc=newpc;
18203425bc38SStefano Zampini   PetscFunctionReturn(0);
18213425bc38SStefano Zampini }
18221e6b0712SBarry Smith 
18233425bc38SStefano Zampini #undef __FUNCT__
18243425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
18253425bc38SStefano Zampini /*@
18260f202f7eSStefano Zampini  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
18273425bc38SStefano Zampini 
18283425bc38SStefano Zampini    Collective
18293425bc38SStefano Zampini 
18303425bc38SStefano Zampini    Input Parameters:
1831*1720468bSStefano Zampini +  pc - the BDDC preconditioning context (setup should have been called before)
1832*1720468bSStefano Zampini -  fully_redundant - true for a fully redundant set of Lagrange multipliers
183328509bceSStefano Zampini 
183428509bceSStefano Zampini    Output Parameters:
18350f202f7eSStefano Zampini +  fetidp_mat - shell FETI-DP matrix object
18360f202f7eSStefano Zampini -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
183728509bceSStefano Zampini 
18383425bc38SStefano Zampini    Level: developer
18393425bc38SStefano Zampini 
18403425bc38SStefano Zampini    Notes:
18410f202f7eSStefano Zampini      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
18423425bc38SStefano Zampini 
18430f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
18443425bc38SStefano Zampini @*/
1845*1720468bSStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, Mat *fetidp_mat, PC *fetidp_pc)
18463425bc38SStefano Zampini {
18473425bc38SStefano Zampini   PetscErrorCode ierr;
18483425bc38SStefano Zampini 
18493425bc38SStefano Zampini   PetscFunctionBegin;
18503425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18513425bc38SStefano Zampini   if (pc->setupcalled) {
1852*1720468bSStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,PetscBool,Mat*,PC*),(pc,fully_redundant,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1853f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
18543425bc38SStefano Zampini   PetscFunctionReturn(0);
18553425bc38SStefano Zampini }
18560c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1857da1bb401SStefano Zampini /*MC
1858da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
18590c7d97c5SJed Brown 
186028509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
186128509bceSStefano Zampini 
186228509bceSStefano Zampini .vb
186328509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
186428509bceSStefano Zampini    [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf
186528509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
18660f202f7eSStefano Zampini    [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf
186728509bceSStefano Zampini .ve
186828509bceSStefano Zampini 
186928509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
187028509bceSStefano Zampini 
18710f202f7eSStefano Zampini    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
187228509bceSStefano Zampini 
187328509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
187428509bceSStefano Zampini 
1875b6fdb6dfSStefano Zampini    Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains.
1876b6fdb6dfSStefano Zampini 
1877c7017625SStefano Zampini    Approximate local solvers are automatically adapted (see [1]) if the user has attached a nullspace object to the subdomain matrices, and informed BDDC of using approximate solvers (via the command line).
187828509bceSStefano Zampini 
18790f202f7eSStefano Zampini    Boundary nodes are split in vertices, edges and faces classes using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph()
18800f202f7eSStefano Zampini    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
188128509bceSStefano Zampini 
18820f202f7eSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
188328509bceSStefano Zampini 
18840f202f7eSStefano Zampini    Change of basis is performed similarly to [2] when requested. When more than one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used.
18850f202f7eSStefano Zampini    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
188628509bceSStefano Zampini 
18870f202f7eSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
188828509bceSStefano Zampini 
18890f202f7eSStefano Zampini    Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS is present. Future versions of the code will also consider using MKL_PARDISO or PASTIX.
189028509bceSStefano Zampini 
18910f202f7eSStefano Zampini    An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases.
18920f202f7eSStefano Zampini    Deluxe scaling is not supported yet for FETI-DP.
18930f202f7eSStefano Zampini 
18940f202f7eSStefano Zampini    Options Database Keys (some of them, run with -h for a complete list):
18950f202f7eSStefano Zampini 
18960f202f7eSStefano Zampini .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
18970f202f7eSStefano Zampini .    -pc_bddc_use_edges <true> - use or not edges in primal space
18980f202f7eSStefano Zampini .    -pc_bddc_use_faces <false> - use or not faces in primal space
18990f202f7eSStefano Zampini .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
19000f202f7eSStefano Zampini .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
19010f202f7eSStefano Zampini .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
19020f202f7eSStefano Zampini .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
190328509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
19040f202f7eSStefano Zampini .    -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case)
19050f202f7eSStefano Zampini .    -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level)
19060f202f7eSStefano Zampini .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
19070f202f7eSStefano Zampini .    -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling)
19080f202f7eSStefano Zampini .    -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS installed)
190928509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
191028509bceSStefano Zampini 
191128509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
191228509bceSStefano Zampini .vb
191328509bceSStefano Zampini       -pc_bddc_dirichlet_
191428509bceSStefano Zampini       -pc_bddc_neumann_
191528509bceSStefano Zampini       -pc_bddc_coarse_
191628509bceSStefano Zampini .ve
19170f202f7eSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
191828509bceSStefano Zampini 
19190f202f7eSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
192028509bceSStefano Zampini .vb
1921312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1922312be037SStefano Zampini       -pc_bddc_neumann_lN_
1923312be037SStefano Zampini       -pc_bddc_coarse_lN_
192428509bceSStefano Zampini .ve
19250f202f7eSStefano Zampini    Note that level number ranges from the finest (0) to the coarsest (N).
19260f202f7eSStefano Zampini    In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g.
19270f202f7eSStefano Zampini .vb
19280f202f7eSStefano Zampini      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
19290f202f7eSStefano Zampini .ve
19300f202f7eSStefano Zampini    will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors
1931da1bb401SStefano Zampini 
1932da1bb401SStefano Zampini    Level: intermediate
1933da1bb401SStefano Zampini 
1934b6fdb6dfSStefano Zampini    Developer notes:
1935da1bb401SStefano Zampini 
1936da1bb401SStefano Zampini    Contributed by Stefano Zampini
1937da1bb401SStefano Zampini 
1938da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1939da1bb401SStefano Zampini M*/
1940b2573a8aSBarry Smith 
1941da1bb401SStefano Zampini #undef __FUNCT__
1942da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
19438cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1944da1bb401SStefano Zampini {
1945da1bb401SStefano Zampini   PetscErrorCode      ierr;
1946da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1947da1bb401SStefano Zampini 
1948da1bb401SStefano Zampini   PetscFunctionBegin;
1949da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1950b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1951da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1952da1bb401SStefano Zampini 
1953da1bb401SStefano Zampini   /* create PCIS data structure */
1954da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1955da1bb401SStefano Zampini 
195647d04d0dSStefano Zampini   /* BDDC customization */
195708a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
195847d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
195947d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
196047d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
196147d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
196247d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
196347d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
1964fa434743SStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE;
1965fa434743SStefano Zampini   pcbddc->use_qr_single       = PETSC_FALSE;
19663301b35fSStefano Zampini   pcbddc->symmetric_primal    = PETSC_TRUE;
196747d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1968b9d89cd5SStefano Zampini   /* private */
1969727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
19700e6343abSStefano Zampini   pcbddc->local_primal_size_cc       = 0;
19710e6343abSStefano Zampini   pcbddc->local_primal_ref_node      = 0;
19720e6343abSStefano Zampini   pcbddc->local_primal_ref_mult      = 0;
1973e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1974727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1975fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
197668457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
1977f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1978727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1979f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1980f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1981f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1982674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
19833972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1984534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1985534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1986534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1987b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
1988906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
1989da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1990da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1991da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1992da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
199315aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
199415aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1995da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1996da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1997da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1998da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1999da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
2000da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
2001da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
2002da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
2003da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
2004da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
2005785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
2006785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
2007785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
200860d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
200960d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
201063602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
2011da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
201263602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
2013da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
201485c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
201547d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
201647d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
2017b0c7d250SStefano Zampini   pcbddc->coarse_adj_red             = 0;
20184fad6a16SStefano Zampini   pcbddc->current_level              = 0;
20192b510759SStefano Zampini   pcbddc->max_levels                 = 0;
2020323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
202174e2c79eSStefano Zampini   pcbddc->redistribute_coarse        = 0;
2022f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
2023323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
2024da1bb401SStefano Zampini 
2025674ae819SStefano Zampini   /* create local graph structure */
2026674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2027674ae819SStefano Zampini 
2028674ae819SStefano Zampini   /* scaling */
2029674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
203034a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2031ac632422SStefano Zampini   pcbddc->faster_deluxe         = PETSC_FALSE;
2032b96c3477SStefano Zampini 
2033b96c3477SStefano Zampini   /* create sub schurs structure */
2034b96c3477SStefano Zampini   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2035b96c3477SStefano Zampini   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2036b96c3477SStefano Zampini   pcbddc->sub_schurs_layers      = -1;
2037b96c3477SStefano Zampini   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2038b96c3477SStefano Zampini 
2039b96c3477SStefano Zampini   pcbddc->computed_rowadj = PETSC_FALSE;
2040da1bb401SStefano Zampini 
2041b7eb3628SStefano Zampini   /* adaptivity */
2042f6f667cfSStefano Zampini   pcbddc->adaptive_threshold      = 0.0;
204308122e43SStefano Zampini   pcbddc->adaptive_nmax           = 0;
2044f6f667cfSStefano Zampini   pcbddc->adaptive_nmin           = 0;
2045b7eb3628SStefano Zampini 
2046da1bb401SStefano Zampini   /* function pointers */
2047da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
204893bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2049da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
2050da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
2051da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
20526b78500eSPatrick Sanan   pc->ops->view                = PCView_BDDC;
2053da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
2054da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
2055da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
2056534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
2057534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
2058da1bb401SStefano Zampini 
2059da1bb401SStefano Zampini   /* composing function */
2060906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2061674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2062bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
20632b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2064b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
20652b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2066bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
206782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2068bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
206982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2070bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
207182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2072bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
207382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2074bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
207563602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2076bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2077bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2078bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2079bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2080a06fd7f2SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_BDDC);CHKERRQ(ierr);
2081da1bb401SStefano Zampini   PetscFunctionReturn(0);
2082da1bb401SStefano Zampini }
20833425bc38SStefano Zampini 
2084