xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 163d334ea1df412efbe58c2e8fce360a2afeaf87)
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    FETIDP
18eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
19eb97c9d2SStefano Zampini 
20eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
21eb97c9d2SStefano Zampini    - Provide general case for subassembling
22eb97c9d2SStefano Zampini 
2353cdbc3dSStefano Zampini */
240c7d97c5SJed Brown 
25ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
26ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
273b03a366Sstefano_zampini #include <petscblaslapack.h>
28674ae819SStefano Zampini 
290c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
300c7d97c5SJed Brown #undef __FUNCT__
310c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
324416b707SBarry Smith PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc)
330c7d97c5SJed Brown {
340c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
350c7d97c5SJed Brown   PetscErrorCode ierr;
360c7d97c5SJed Brown 
370c7d97c5SJed Brown   PetscFunctionBegin;
38e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
398eeda7d8SStefano Zampini   /* Verbose debugging */
408eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
416b78500eSPatrick Sanan   /* Primal space customization */
4208a5cf49SStefano 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);
438eeda7d8SStefano 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);
448eeda7d8SStefano 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);
458eeda7d8SStefano 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);
466661aa1dSStefano 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);
47ac632422SStefano 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);
488eeda7d8SStefano Zampini   /* Change of basis */
49b9b85e73SStefano 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);
50b9b85e73SStefano 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);
51674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
52674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
53674ae819SStefano Zampini   }
548eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
558eeda7d8SStefano 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);
5674e2c79eSStefano 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);
570298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
582b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
59323d395dSStefano 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);
60674ae819SStefano 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);
61b96c3477SStefano 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);
62b96c3477SStefano 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);
63b96c3477SStefano 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);
64ac632422SStefano 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);
654c6709b3SStefano 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);
6608122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
6708122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
683301b35fSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
69b0c7d250SStefano 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);
700c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
710c7d97c5SJed Brown   PetscFunctionReturn(0);
720c7d97c5SJed Brown }
736b78500eSPatrick Sanan 
746b78500eSPatrick Sanan /* -------------------------------------------------------------------------- */
756b78500eSPatrick Sanan #undef __FUNCT__
766b78500eSPatrick Sanan #define __FUNCT__ "PCView_BDDC"
776b78500eSPatrick Sanan static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer)
786b78500eSPatrick Sanan {
796b78500eSPatrick Sanan   PC_BDDC              *pcbddc = (PC_BDDC*)pc->data;
806b78500eSPatrick Sanan   PetscErrorCode       ierr;
816b78500eSPatrick Sanan   PetscBool            isascii,isstring;
826b78500eSPatrick Sanan 
836b78500eSPatrick Sanan   PetscFunctionBegin;
846b78500eSPatrick Sanan   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
856b78500eSPatrick Sanan   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
866b78500eSPatrick Sanan 
876b78500eSPatrick Sanan   /* In a braindead way, print out anything which the user can control from the command line,
886b78500eSPatrick Sanan      cribbing from PCSetFromOptions_BDDC */
896b78500eSPatrick Sanan 
906b78500eSPatrick Sanan   /* Nothing printed for the String viewer */
916b78500eSPatrick Sanan 
926b78500eSPatrick Sanan   /* ASCII viewer */
936b78500eSPatrick Sanan   if (isascii) {
946b78500eSPatrick Sanan     /* Verbose debugging */
956b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use verbose output: %d\n",pcbddc->dbg_flag);CHKERRQ(ierr);
966b78500eSPatrick Sanan 
976b78500eSPatrick Sanan     /* Primal space customization */
986b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use local mat graph: %d\n",pcbddc->use_local_adj);CHKERRQ(ierr);
996b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use vertices: %d\n",pcbddc->use_vertices);CHKERRQ(ierr);
1006b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use edges: %d\n",pcbddc->use_edges);CHKERRQ(ierr);
1016b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use faces: %d\n",pcbddc->use_faces);CHKERRQ(ierr);
1026b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use true near null space: %d\n",pcbddc->use_nnsp_true);CHKERRQ(ierr);
1036b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);CHKERRQ(ierr);
1046b78500eSPatrick Sanan 
1056b78500eSPatrick Sanan     /* Change of basis */
1066b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);CHKERRQ(ierr);
1076b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);CHKERRQ(ierr);
1086b78500eSPatrick Sanan 
1096b78500eSPatrick Sanan     /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
1106b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);CHKERRQ(ierr);
1116b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Coarse problem restribute procs: %d\n",pcbddc->redistribute_coarse);CHKERRQ(ierr);
1126b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Multilevel coarsening ratio: %d\n",pcbddc->coarsening_ratio);CHKERRQ(ierr);
1136b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Multilevel max levels: %d\n",pcbddc->max_levels);CHKERRQ(ierr);
1146b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);CHKERRQ(ierr);
1156b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);CHKERRQ(ierr);
1166b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);CHKERRQ(ierr);
1176b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Number of dofs' layers for the computation of principal minors: %d\n",pcbddc->sub_schurs_layers);CHKERRQ(ierr);
1186b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);CHKERRQ(ierr);
1196b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Fast deluxe scaling: %d\n",pcbddc->faster_deluxe);CHKERRQ(ierr);
1206b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Adaptive constraint selection threshold: %g\n",pcbddc->adaptive_threshold);CHKERRQ(ierr);
1216b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Min constraints / connected component: %d\n",pcbddc->adaptive_nmin);CHKERRQ(ierr);
1226b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Max constraints / connected component: %d\n",pcbddc->adaptive_nmax);CHKERRQ(ierr);
1236b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);CHKERRQ(ierr);
1246b78500eSPatrick Sanan     ierr = PetscViewerASCIIPrintf(viewer,    "  BDDC: Num. Procs. to map coarse adjacency list: %d\n",pcbddc->coarse_adj_red);CHKERRQ(ierr);
1256b78500eSPatrick Sanan   }
1266b78500eSPatrick Sanan 
1276b78500eSPatrick Sanan   PetscFunctionReturn(0);
1286b78500eSPatrick Sanan }
1296b78500eSPatrick Sanan 
1300c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
131674ae819SStefano Zampini #undef __FUNCT__
132906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
133906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
134b9b85e73SStefano Zampini {
135b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
136b9b85e73SStefano Zampini   PetscErrorCode ierr;
137b9b85e73SStefano Zampini 
138b9b85e73SStefano Zampini   PetscFunctionBegin;
139b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
140b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
141b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
142b9b85e73SStefano Zampini   PetscFunctionReturn(0);
143b9b85e73SStefano Zampini }
144b9b85e73SStefano Zampini #undef __FUNCT__
145906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
146b9b85e73SStefano Zampini /*@
147906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
148b9b85e73SStefano Zampini 
149b9b85e73SStefano Zampini    Collective on PC
150b9b85e73SStefano Zampini 
151b9b85e73SStefano Zampini    Input Parameters:
152b9b85e73SStefano Zampini +  pc - the preconditioning context
153906d46d4SStefano Zampini -  change - the change of basis matrix
154b9b85e73SStefano Zampini 
155b9b85e73SStefano Zampini    Level: intermediate
156b9b85e73SStefano Zampini 
157b9b85e73SStefano Zampini    Notes:
158b9b85e73SStefano Zampini 
159b9b85e73SStefano Zampini .seealso: PCBDDC
160b9b85e73SStefano Zampini @*/
161906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
162b9b85e73SStefano Zampini {
163b9b85e73SStefano Zampini   PetscErrorCode ierr;
164b9b85e73SStefano Zampini 
165b9b85e73SStefano Zampini   PetscFunctionBegin;
166b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
167b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
168906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
169906d46d4SStefano Zampini   if (pc->mat) {
170906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
171906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
172906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
1736c4ed002SBarry 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);
1746c4ed002SBarry 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);
175906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
176906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
1776c4ed002SBarry 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);
1786c4ed002SBarry 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);
179906d46d4SStefano Zampini   }
180906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
181b9b85e73SStefano Zampini   PetscFunctionReturn(0);
182b9b85e73SStefano Zampini }
183b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
184b9b85e73SStefano Zampini #undef __FUNCT__
185674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
186674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
187674ae819SStefano Zampini {
188674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
189674ae819SStefano Zampini   PetscErrorCode ierr;
1901e6b0712SBarry Smith 
191674ae819SStefano Zampini   PetscFunctionBegin;
192674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
193674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
194674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
195674ae819SStefano Zampini   PetscFunctionReturn(0);
196674ae819SStefano Zampini }
197674ae819SStefano Zampini #undef __FUNCT__
198674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
199674ae819SStefano Zampini /*@
20028509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
201674ae819SStefano Zampini 
20217eb1463SStefano Zampini    Collective
203674ae819SStefano Zampini 
204674ae819SStefano Zampini    Input Parameters:
205674ae819SStefano Zampini +  pc - the preconditioning context
20617eb1463SStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
207674ae819SStefano Zampini 
208674ae819SStefano Zampini    Level: intermediate
209674ae819SStefano Zampini 
210674ae819SStefano Zampini    Notes:
211674ae819SStefano Zampini 
212674ae819SStefano Zampini .seealso: PCBDDC
213674ae819SStefano Zampini @*/
214674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
215674ae819SStefano Zampini {
216674ae819SStefano Zampini   PetscErrorCode ierr;
217674ae819SStefano Zampini 
218674ae819SStefano Zampini   PetscFunctionBegin;
219674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
220674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
22117eb1463SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
222674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
223674ae819SStefano Zampini   PetscFunctionReturn(0);
224674ae819SStefano Zampini }
225674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
2260c7d97c5SJed Brown #undef __FUNCT__
2274fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
2284fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
2294fad6a16SStefano Zampini {
2304fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2314fad6a16SStefano Zampini 
2324fad6a16SStefano Zampini   PetscFunctionBegin;
2334fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
2344fad6a16SStefano Zampini   PetscFunctionReturn(0);
2354fad6a16SStefano Zampini }
2361e6b0712SBarry Smith 
2374fad6a16SStefano Zampini #undef __FUNCT__
2384fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
2394fad6a16SStefano Zampini /*@
24028509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
2414fad6a16SStefano Zampini 
2424fad6a16SStefano Zampini    Logically collective on PC
2434fad6a16SStefano Zampini 
2444fad6a16SStefano Zampini    Input Parameters:
2454fad6a16SStefano Zampini +  pc - the preconditioning context
24628509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2474fad6a16SStefano Zampini 
2480f202f7eSStefano Zampini    Options Database Keys:
2490f202f7eSStefano Zampini .    -pc_bddc_coarsening_ratio
2504fad6a16SStefano Zampini 
2514fad6a16SStefano Zampini    Level: intermediate
2524fad6a16SStefano Zampini 
2534fad6a16SStefano Zampini    Notes:
2540f202f7eSStefano Zampini      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2554fad6a16SStefano Zampini 
2560f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetLevels()
2574fad6a16SStefano Zampini @*/
2584fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2594fad6a16SStefano Zampini {
2604fad6a16SStefano Zampini   PetscErrorCode ierr;
2614fad6a16SStefano Zampini 
2624fad6a16SStefano Zampini   PetscFunctionBegin;
2634fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2642b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2654fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2664fad6a16SStefano Zampini   PetscFunctionReturn(0);
2674fad6a16SStefano Zampini }
2682b510759SStefano Zampini 
269b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2702b510759SStefano Zampini #undef __FUNCT__
271b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
272b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
273b8ffe317SStefano Zampini {
274b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
275b8ffe317SStefano Zampini 
276b8ffe317SStefano Zampini   PetscFunctionBegin;
27785c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
278b8ffe317SStefano Zampini   PetscFunctionReturn(0);
279b8ffe317SStefano Zampini }
280b8ffe317SStefano Zampini 
281b8ffe317SStefano Zampini #undef __FUNCT__
282b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
283b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2842b510759SStefano Zampini {
2852b510759SStefano Zampini   PetscErrorCode ierr;
2862b510759SStefano Zampini 
2872b510759SStefano Zampini   PetscFunctionBegin;
2882b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
289b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
290b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2912b510759SStefano Zampini   PetscFunctionReturn(0);
2922b510759SStefano Zampini }
2931e6b0712SBarry Smith 
2944fad6a16SStefano Zampini #undef __FUNCT__
2952b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2962b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2974fad6a16SStefano Zampini {
2984fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2994fad6a16SStefano Zampini 
3004fad6a16SStefano Zampini   PetscFunctionBegin;
3012b510759SStefano Zampini   pcbddc->current_level = level;
3024fad6a16SStefano Zampini   PetscFunctionReturn(0);
3034fad6a16SStefano Zampini }
3041e6b0712SBarry Smith 
3054fad6a16SStefano Zampini #undef __FUNCT__
306b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
307b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
308b8ffe317SStefano Zampini {
309b8ffe317SStefano Zampini   PetscErrorCode ierr;
310b8ffe317SStefano Zampini 
311b8ffe317SStefano Zampini   PetscFunctionBegin;
312b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
313b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
314b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
315b8ffe317SStefano Zampini   PetscFunctionReturn(0);
316b8ffe317SStefano Zampini }
317b8ffe317SStefano Zampini 
318b8ffe317SStefano Zampini #undef __FUNCT__
3192b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
3202b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
3212b510759SStefano Zampini {
3222b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3232b510759SStefano Zampini 
3242b510759SStefano Zampini   PetscFunctionBegin;
3252b510759SStefano Zampini   pcbddc->max_levels = levels;
3262b510759SStefano Zampini   PetscFunctionReturn(0);
3272b510759SStefano Zampini }
3282b510759SStefano Zampini 
3292b510759SStefano Zampini #undef __FUNCT__
3302b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
3314fad6a16SStefano Zampini /*@
33228509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
3334fad6a16SStefano Zampini 
3344fad6a16SStefano Zampini    Logically collective on PC
3354fad6a16SStefano Zampini 
3364fad6a16SStefano Zampini    Input Parameters:
3374fad6a16SStefano Zampini +  pc - the preconditioning context
33828509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
3394fad6a16SStefano Zampini 
3400f202f7eSStefano Zampini    Options Database Keys:
3410f202f7eSStefano Zampini .    -pc_bddc_levels
3424fad6a16SStefano Zampini 
3434fad6a16SStefano Zampini    Level: intermediate
3444fad6a16SStefano Zampini 
3454fad6a16SStefano Zampini    Notes:
3460f202f7eSStefano Zampini      Default value is 0, i.e. traditional one-level BDDC
3474fad6a16SStefano Zampini 
3480f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
3494fad6a16SStefano Zampini @*/
3502b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3514fad6a16SStefano Zampini {
3524fad6a16SStefano Zampini   PetscErrorCode ierr;
3534fad6a16SStefano Zampini 
3544fad6a16SStefano Zampini   PetscFunctionBegin;
3554fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3562b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
357312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3582b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3594fad6a16SStefano Zampini   PetscFunctionReturn(0);
3604fad6a16SStefano Zampini }
3614fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3621e6b0712SBarry Smith 
3634fad6a16SStefano Zampini #undef __FUNCT__
3640bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3650bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3660bdf917eSStefano Zampini {
3670bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3680bdf917eSStefano Zampini   PetscErrorCode ierr;
3690bdf917eSStefano Zampini 
3700bdf917eSStefano Zampini   PetscFunctionBegin;
3710bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3720bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3730bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3740bdf917eSStefano Zampini   PetscFunctionReturn(0);
3750bdf917eSStefano Zampini }
3761e6b0712SBarry Smith 
3770bdf917eSStefano Zampini #undef __FUNCT__
3780bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3790bdf917eSStefano Zampini /*@
38028509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3810bdf917eSStefano Zampini 
3820bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3830bdf917eSStefano Zampini 
3840bdf917eSStefano Zampini    Input Parameters:
3850bdf917eSStefano Zampini +  pc - the preconditioning context
38628509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3870bdf917eSStefano Zampini 
3880bdf917eSStefano Zampini    Level: intermediate
3890bdf917eSStefano Zampini 
3900bdf917eSStefano Zampini    Notes:
3910bdf917eSStefano Zampini 
3920bdf917eSStefano Zampini .seealso: PCBDDC
3930bdf917eSStefano Zampini @*/
3940bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3950bdf917eSStefano Zampini {
3960bdf917eSStefano Zampini   PetscErrorCode ierr;
3970bdf917eSStefano Zampini 
3980bdf917eSStefano Zampini   PetscFunctionBegin;
3990bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
400674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
4012b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
4020bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
4030bdf917eSStefano Zampini   PetscFunctionReturn(0);
4040bdf917eSStefano Zampini }
4050bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
4061e6b0712SBarry Smith 
4070bdf917eSStefano Zampini #undef __FUNCT__
4083b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
4093b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
4103b03a366Sstefano_zampini {
4113b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4123b03a366Sstefano_zampini   PetscErrorCode ierr;
4133b03a366Sstefano_zampini 
4143b03a366Sstefano_zampini   PetscFunctionBegin;
415785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
416785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4173b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
41836e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
41936e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
420fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4213b03a366Sstefano_zampini   PetscFunctionReturn(0);
4223b03a366Sstefano_zampini }
4231e6b0712SBarry Smith 
4243b03a366Sstefano_zampini #undef __FUNCT__
4253b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
4263b03a366Sstefano_zampini /*@
42728509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
4283b03a366Sstefano_zampini 
429785d1243SStefano Zampini    Collective
4303b03a366Sstefano_zampini 
4313b03a366Sstefano_zampini    Input Parameters:
4323b03a366Sstefano_zampini +  pc - the preconditioning context
433785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
4343b03a366Sstefano_zampini 
4353b03a366Sstefano_zampini    Level: intermediate
4363b03a366Sstefano_zampini 
4370f202f7eSStefano Zampini    Notes:
4380f202f7eSStefano Zampini      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
4393b03a366Sstefano_zampini 
4400f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
4413b03a366Sstefano_zampini @*/
4423b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
4433b03a366Sstefano_zampini {
4443b03a366Sstefano_zampini   PetscErrorCode ierr;
4453b03a366Sstefano_zampini 
4463b03a366Sstefano_zampini   PetscFunctionBegin;
4473b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
448674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
449785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4503b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4513b03a366Sstefano_zampini   PetscFunctionReturn(0);
4523b03a366Sstefano_zampini }
4533b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4541e6b0712SBarry Smith 
4553b03a366Sstefano_zampini #undef __FUNCT__
45682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
45782ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4580c7d97c5SJed Brown {
4590c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4600c7d97c5SJed Brown   PetscErrorCode ierr;
4610c7d97c5SJed Brown 
4620c7d97c5SJed Brown   PetscFunctionBegin;
463785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
464785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4650c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4660c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
467785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4687d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4690c7d97c5SJed Brown   PetscFunctionReturn(0);
4700c7d97c5SJed Brown }
4710c7d97c5SJed Brown 
4720c7d97c5SJed Brown #undef __FUNCT__
47382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4749c0446d6SStefano Zampini /*@
47582ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4769c0446d6SStefano Zampini 
477785d1243SStefano Zampini    Collective
47853cdbc3dSStefano Zampini 
47953cdbc3dSStefano Zampini    Input Parameters:
48053cdbc3dSStefano Zampini +  pc - the preconditioning context
48182ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
48253cdbc3dSStefano Zampini 
48353cdbc3dSStefano Zampini    Level: intermediate
48453cdbc3dSStefano Zampini 
4859c0446d6SStefano Zampini    Notes:
48653cdbc3dSStefano Zampini 
4870f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
48853cdbc3dSStefano Zampini @*/
48982ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4900c7d97c5SJed Brown {
4910c7d97c5SJed Brown   PetscErrorCode ierr;
4920c7d97c5SJed Brown 
4930c7d97c5SJed Brown   PetscFunctionBegin;
4940c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4950c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
49682ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
49782ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4980c7d97c5SJed Brown   PetscFunctionReturn(0);
4990c7d97c5SJed Brown }
5000c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
5010c7d97c5SJed Brown 
5020c7d97c5SJed Brown #undef __FUNCT__
5030c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
50453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
5050c7d97c5SJed Brown {
5060c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
50753cdbc3dSStefano Zampini   PetscErrorCode ierr;
5080c7d97c5SJed Brown 
5090c7d97c5SJed Brown   PetscFunctionBegin;
510785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
511785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
51253cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
51336e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
51436e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
515fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5160c7d97c5SJed Brown   PetscFunctionReturn(0);
5170c7d97c5SJed Brown }
5181e6b0712SBarry Smith 
5190c7d97c5SJed Brown #undef __FUNCT__
5200c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
52157527edcSJed Brown /*@
52228509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
52357527edcSJed Brown 
524785d1243SStefano Zampini    Collective
52557527edcSJed Brown 
52657527edcSJed Brown    Input Parameters:
52757527edcSJed Brown +  pc - the preconditioning context
528785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
52957527edcSJed Brown 
53057527edcSJed Brown    Level: intermediate
53157527edcSJed Brown 
5320f202f7eSStefano Zampini    Notes:
5330f202f7eSStefano Zampini      Any process can list any global node
53457527edcSJed Brown 
5350f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
53657527edcSJed Brown @*/
53753cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
5380c7d97c5SJed Brown {
5390c7d97c5SJed Brown   PetscErrorCode ierr;
5400c7d97c5SJed Brown 
5410c7d97c5SJed Brown   PetscFunctionBegin;
5420c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
543674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
544785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54553cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
54653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
54753cdbc3dSStefano Zampini }
54853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5491e6b0712SBarry Smith 
55053cdbc3dSStefano Zampini #undef __FUNCT__
55182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
55282ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5530c7d97c5SJed Brown {
5540c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5550c7d97c5SJed Brown   PetscErrorCode ierr;
5560c7d97c5SJed Brown 
5570c7d97c5SJed Brown   PetscFunctionBegin;
558785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
559785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5600c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5610c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
562785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5637d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5640c7d97c5SJed Brown   PetscFunctionReturn(0);
5650c7d97c5SJed Brown }
5660c7d97c5SJed Brown 
5670c7d97c5SJed Brown #undef __FUNCT__
56882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5690c7d97c5SJed Brown /*@
57082ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5710c7d97c5SJed Brown 
572785d1243SStefano Zampini    Collective
5730c7d97c5SJed Brown 
5740c7d97c5SJed Brown    Input Parameters:
5750c7d97c5SJed Brown +  pc - the preconditioning context
57682ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5770c7d97c5SJed Brown 
5780c7d97c5SJed Brown    Level: intermediate
5790c7d97c5SJed Brown 
5800c7d97c5SJed Brown    Notes:
5810c7d97c5SJed Brown 
5820f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
5830c7d97c5SJed Brown @*/
58482ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5850c7d97c5SJed Brown {
5860c7d97c5SJed Brown   PetscErrorCode ierr;
5870c7d97c5SJed Brown 
5880c7d97c5SJed Brown   PetscFunctionBegin;
5890c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5900c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
59182ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
59282ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
59353cdbc3dSStefano Zampini   PetscFunctionReturn(0);
59453cdbc3dSStefano Zampini }
59553cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
59653cdbc3dSStefano Zampini 
59753cdbc3dSStefano Zampini #undef __FUNCT__
598da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
599da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
600da1bb401SStefano Zampini {
601da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
602da1bb401SStefano Zampini 
603da1bb401SStefano Zampini   PetscFunctionBegin;
604da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
605da1bb401SStefano Zampini   PetscFunctionReturn(0);
606da1bb401SStefano Zampini }
6071e6b0712SBarry Smith 
608da1bb401SStefano Zampini #undef __FUNCT__
609da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
610da1bb401SStefano Zampini /*@
611785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
612da1bb401SStefano Zampini 
613785d1243SStefano Zampini    Collective
614785d1243SStefano Zampini 
615785d1243SStefano Zampini    Input Parameters:
616785d1243SStefano Zampini .  pc - the preconditioning context
617785d1243SStefano Zampini 
618785d1243SStefano Zampini    Output Parameters:
619785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
620785d1243SStefano Zampini 
621785d1243SStefano Zampini    Level: intermediate
622785d1243SStefano Zampini 
6230f202f7eSStefano Zampini    Notes:
6240f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
625785d1243SStefano Zampini 
626785d1243SStefano Zampini .seealso: PCBDDC
627785d1243SStefano Zampini @*/
628785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
629785d1243SStefano Zampini {
630785d1243SStefano Zampini   PetscErrorCode ierr;
631785d1243SStefano Zampini 
632785d1243SStefano Zampini   PetscFunctionBegin;
633785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
634785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
635785d1243SStefano Zampini   PetscFunctionReturn(0);
636785d1243SStefano Zampini }
637785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
638785d1243SStefano Zampini 
639785d1243SStefano Zampini #undef __FUNCT__
640785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
641785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
642785d1243SStefano Zampini {
643785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
644785d1243SStefano Zampini 
645785d1243SStefano Zampini   PetscFunctionBegin;
646785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
647785d1243SStefano Zampini   PetscFunctionReturn(0);
648785d1243SStefano Zampini }
649785d1243SStefano Zampini 
650785d1243SStefano Zampini #undef __FUNCT__
65182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
652da1bb401SStefano Zampini /*@
65382ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
654da1bb401SStefano Zampini 
655785d1243SStefano Zampini    Collective
656da1bb401SStefano Zampini 
657da1bb401SStefano Zampini    Input Parameters:
65828509bceSStefano Zampini .  pc - the preconditioning context
659da1bb401SStefano Zampini 
660da1bb401SStefano Zampini    Output Parameters:
66128509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
662da1bb401SStefano Zampini 
663da1bb401SStefano Zampini    Level: intermediate
664da1bb401SStefano Zampini 
665da1bb401SStefano Zampini    Notes:
6660f202f7eSStefano 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).
6670f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
668da1bb401SStefano Zampini 
669da1bb401SStefano Zampini .seealso: PCBDDC
670da1bb401SStefano Zampini @*/
67182ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
672da1bb401SStefano Zampini {
673da1bb401SStefano Zampini   PetscErrorCode ierr;
674da1bb401SStefano Zampini 
675da1bb401SStefano Zampini   PetscFunctionBegin;
676da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
67782ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
678da1bb401SStefano Zampini   PetscFunctionReturn(0);
679da1bb401SStefano Zampini }
680da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6811e6b0712SBarry Smith 
682da1bb401SStefano Zampini #undef __FUNCT__
68353cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
68453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
68553cdbc3dSStefano Zampini {
68653cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
68753cdbc3dSStefano Zampini 
68853cdbc3dSStefano Zampini   PetscFunctionBegin;
68953cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
69053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
69153cdbc3dSStefano Zampini }
6921e6b0712SBarry Smith 
69353cdbc3dSStefano Zampini #undef __FUNCT__
69453cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
69553cdbc3dSStefano Zampini /*@
696785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
69753cdbc3dSStefano Zampini 
698785d1243SStefano Zampini    Collective
699785d1243SStefano Zampini 
700785d1243SStefano Zampini    Input Parameters:
701785d1243SStefano Zampini .  pc - the preconditioning context
702785d1243SStefano Zampini 
703785d1243SStefano Zampini    Output Parameters:
704785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
705785d1243SStefano Zampini 
706785d1243SStefano Zampini    Level: intermediate
707785d1243SStefano Zampini 
7080f202f7eSStefano Zampini    Notes:
7090f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
710785d1243SStefano Zampini 
711785d1243SStefano Zampini .seealso: PCBDDC
712785d1243SStefano Zampini @*/
713785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
714785d1243SStefano Zampini {
715785d1243SStefano Zampini   PetscErrorCode ierr;
716785d1243SStefano Zampini 
717785d1243SStefano Zampini   PetscFunctionBegin;
718785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
719785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
720785d1243SStefano Zampini   PetscFunctionReturn(0);
721785d1243SStefano Zampini }
722785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
723785d1243SStefano Zampini 
724785d1243SStefano Zampini #undef __FUNCT__
725785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
726785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
727785d1243SStefano Zampini {
728785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
729785d1243SStefano Zampini 
730785d1243SStefano Zampini   PetscFunctionBegin;
731785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
732785d1243SStefano Zampini   PetscFunctionReturn(0);
733785d1243SStefano Zampini }
734785d1243SStefano Zampini 
735785d1243SStefano Zampini #undef __FUNCT__
73682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
73753cdbc3dSStefano Zampini /*@
73882ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
73953cdbc3dSStefano Zampini 
740785d1243SStefano Zampini    Collective
74153cdbc3dSStefano Zampini 
74253cdbc3dSStefano Zampini    Input Parameters:
74328509bceSStefano Zampini .  pc - the preconditioning context
74453cdbc3dSStefano Zampini 
74553cdbc3dSStefano Zampini    Output Parameters:
74628509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
74753cdbc3dSStefano Zampini 
74853cdbc3dSStefano Zampini    Level: intermediate
74953cdbc3dSStefano Zampini 
75053cdbc3dSStefano Zampini    Notes:
7510f202f7eSStefano 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).
7520f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
75353cdbc3dSStefano Zampini 
75453cdbc3dSStefano Zampini .seealso: PCBDDC
75553cdbc3dSStefano Zampini @*/
75682ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
75753cdbc3dSStefano Zampini {
75853cdbc3dSStefano Zampini   PetscErrorCode ierr;
75953cdbc3dSStefano Zampini 
76053cdbc3dSStefano Zampini   PetscFunctionBegin;
76153cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
76282ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7630c7d97c5SJed Brown   PetscFunctionReturn(0);
7640c7d97c5SJed Brown }
76536e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7661e6b0712SBarry Smith 
76736e030ebSStefano Zampini #undef __FUNCT__
768da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7691a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
77036e030ebSStefano Zampini {
77136e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
772da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
773da1bb401SStefano Zampini   PetscErrorCode ierr;
77436e030ebSStefano Zampini 
77536e030ebSStefano Zampini   PetscFunctionBegin;
776674ae819SStefano Zampini   /* free old CSR */
777674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
778fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
779674ae819SStefano Zampini   /* get CSR into graph structure */
780da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
781854ce69bSBarry Smith     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
782785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
783674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
784674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
785da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7861a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7871a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
788674ae819SStefano Zampini   }
789575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
79036e030ebSStefano Zampini   PetscFunctionReturn(0);
79136e030ebSStefano Zampini }
7921e6b0712SBarry Smith 
79336e030ebSStefano Zampini #undef __FUNCT__
794da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
79536e030ebSStefano Zampini /*@
7960f202f7eSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
79736e030ebSStefano Zampini 
79836e030ebSStefano Zampini    Not collective
79936e030ebSStefano Zampini 
80036e030ebSStefano Zampini    Input Parameters:
80136e030ebSStefano Zampini +  pc - the preconditioning context
8020f202f7eSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
80328509bceSStefano Zampini .  xadj, adjncy - the CSR graph
80428509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
80536e030ebSStefano Zampini 
80636e030ebSStefano Zampini    Level: intermediate
80736e030ebSStefano Zampini 
80836e030ebSStefano Zampini    Notes:
80936e030ebSStefano Zampini 
81028509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
81136e030ebSStefano Zampini @*/
8121a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
81336e030ebSStefano Zampini {
814575ad6abSStefano Zampini   void (*f)(void) = 0;
81536e030ebSStefano Zampini   PetscErrorCode ierr;
81636e030ebSStefano Zampini 
81736e030ebSStefano Zampini   PetscFunctionBegin;
81836e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
819674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
820d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
8216c4ed002SBarry Smith   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER)  SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d",copymode);
8221a83f524SJed Brown   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
823575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
824575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
825575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
826575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
827575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
828da1bb401SStefano Zampini   }
82936e030ebSStefano Zampini   PetscFunctionReturn(0);
83036e030ebSStefano Zampini }
8319c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
8321e6b0712SBarry Smith 
8339c0446d6SStefano Zampini #undef __FUNCT__
83463602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
83563602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
83663602bcaSStefano Zampini {
83763602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
83863602bcaSStefano Zampini   PetscInt i;
83963602bcaSStefano Zampini   PetscErrorCode ierr;
84063602bcaSStefano Zampini 
84163602bcaSStefano Zampini   PetscFunctionBegin;
84263602bcaSStefano Zampini   /* Destroy ISes if they were already set */
84363602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
84463602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
84563602bcaSStefano Zampini   }
84663602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
84763602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
84863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
84963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
85063602bcaSStefano Zampini   }
85163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
85263602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
85363602bcaSStefano Zampini   /* allocate space then set */
854d02579f5SStefano Zampini   if (n_is) {
855d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
856d02579f5SStefano Zampini   }
85763602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
85863602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
85963602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
86063602bcaSStefano Zampini   }
86163602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
86263602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8631a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
86463602bcaSStefano Zampini   PetscFunctionReturn(0);
86563602bcaSStefano Zampini }
86663602bcaSStefano Zampini 
86763602bcaSStefano Zampini #undef __FUNCT__
86863602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
86963602bcaSStefano Zampini /*@
87063602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
87163602bcaSStefano Zampini 
87263602bcaSStefano Zampini    Collective
87363602bcaSStefano Zampini 
87463602bcaSStefano Zampini    Input Parameters:
87563602bcaSStefano Zampini +  pc - the preconditioning context
8760f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
8770f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in local ordering
87863602bcaSStefano Zampini 
87963602bcaSStefano Zampini    Level: intermediate
88063602bcaSStefano Zampini 
8810f202f7eSStefano Zampini    Notes:
8820f202f7eSStefano Zampini      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
88363602bcaSStefano Zampini 
88463602bcaSStefano Zampini .seealso: PCBDDC
88563602bcaSStefano Zampini @*/
88663602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
88763602bcaSStefano Zampini {
88863602bcaSStefano Zampini   PetscInt       i;
88963602bcaSStefano Zampini   PetscErrorCode ierr;
89063602bcaSStefano Zampini 
89163602bcaSStefano Zampini   PetscFunctionBegin;
89263602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
89363602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
89463602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
89563602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
89663602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
89763602bcaSStefano Zampini   }
898e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
89963602bcaSStefano Zampini   PetscFunctionReturn(0);
90063602bcaSStefano Zampini }
90163602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
90263602bcaSStefano Zampini 
90363602bcaSStefano Zampini #undef __FUNCT__
9049c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
9059c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
9069c0446d6SStefano Zampini {
9079c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
9089c0446d6SStefano Zampini   PetscInt i;
9099c0446d6SStefano Zampini   PetscErrorCode ierr;
9109c0446d6SStefano Zampini 
9119c0446d6SStefano Zampini   PetscFunctionBegin;
912da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
9139c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
9149c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
9159c0446d6SStefano Zampini   }
916d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
91763602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
91863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
91963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
92063602bcaSStefano Zampini   }
92163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
92263602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
923da1bb401SStefano Zampini   /* allocate space then set */
924d02579f5SStefano Zampini   if (n_is) {
925785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
926d02579f5SStefano Zampini   }
9279c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
928da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
929da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
9309c0446d6SStefano Zampini   }
9319c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
93263602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
9331a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
9349c0446d6SStefano Zampini   PetscFunctionReturn(0);
9359c0446d6SStefano Zampini }
9361e6b0712SBarry Smith 
9379c0446d6SStefano Zampini #undef __FUNCT__
9389c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
9399c0446d6SStefano Zampini /*@
94063602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
9419c0446d6SStefano Zampini 
94263602bcaSStefano Zampini    Collective
9439c0446d6SStefano Zampini 
9449c0446d6SStefano Zampini    Input Parameters:
9459c0446d6SStefano Zampini +  pc - the preconditioning context
9460f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
9470f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in global ordering
9489c0446d6SStefano Zampini 
9499c0446d6SStefano Zampini    Level: intermediate
9509c0446d6SStefano Zampini 
9510f202f7eSStefano Zampini    Notes:
9520f202f7eSStefano Zampini      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
9539c0446d6SStefano Zampini 
9549c0446d6SStefano Zampini .seealso: PCBDDC
9559c0446d6SStefano Zampini @*/
9569c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9579c0446d6SStefano Zampini {
9582b510759SStefano Zampini   PetscInt       i;
9599c0446d6SStefano Zampini   PetscErrorCode ierr;
9609c0446d6SStefano Zampini 
9619c0446d6SStefano Zampini   PetscFunctionBegin;
9629c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
96363602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9642b510759SStefano Zampini   for (i=0;i<n_is;i++) {
96563602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
96663602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9672b510759SStefano Zampini   }
9689c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9699c0446d6SStefano Zampini   PetscFunctionReturn(0);
9709c0446d6SStefano Zampini }
971906d46d4SStefano Zampini 
972da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
973534831adSStefano Zampini #undef __FUNCT__
974534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
975534831adSStefano Zampini /* -------------------------------------------------------------------------- */
976534831adSStefano Zampini /*
977534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
978534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9799c0446d6SStefano Zampini 
980534831adSStefano Zampini    Input Parameter:
981534831adSStefano Zampini +  pc - the preconditioner contex
982534831adSStefano Zampini 
983534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
984534831adSStefano Zampini 
985534831adSStefano Zampini    Notes:
986534831adSStefano Zampini      The interface routine PCPreSolve() is not usually called directly by
987534831adSStefano Zampini    the user, but instead is called by KSPSolve().
988534831adSStefano Zampini */
989534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
990534831adSStefano Zampini {
991534831adSStefano Zampini   PetscErrorCode ierr;
992534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
993534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9943972b0daSStefano Zampini   Vec            used_vec;
9958d00608fSStefano Zampini   PetscBool      copy_rhs = PETSC_TRUE;
996534831adSStefano Zampini 
997534831adSStefano Zampini   PetscFunctionBegin;
99885c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
99985c4d303SStefano Zampini   if (ksp) {
100085c4d303SStefano Zampini     PetscBool iscg;
100185c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
100285c4d303SStefano Zampini     if (!iscg) {
100385c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
100485c4d303SStefano Zampini     }
100585c4d303SStefano Zampini   }
100685c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
100762a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
100862a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
100962a6ff1dSStefano Zampini   }
101062a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
101162a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
101262a6ff1dSStefano Zampini   }
10138d00608fSStefano Zampini 
10143972b0daSStefano Zampini   if (x) {
10153972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
10163972b0daSStefano Zampini     used_vec = x;
10178d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
10183972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
10193972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
10203972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
10213972b0daSStefano Zampini   }
10228efcfb23SStefano Zampini 
10238efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
10243972b0daSStefano Zampini   if (ksp) {
1025a0cb1b98SStefano Zampini     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
10268efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
10278efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
10283972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
10293972b0daSStefano Zampini     }
10303972b0daSStefano Zampini   }
10313308cffdSStefano Zampini 
10328d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
10333972b0daSStefano Zampini 
10343972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
1035a07ea27aSStefano Zampini   if (rhs) {
10363975b054SStefano Zampini     IS dirIS = NULL;
10373975b054SStefano Zampini 
1038a07ea27aSStefano Zampini     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
10393975b054SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
10403975b054SStefano Zampini     if (dirIS) {
1041906d46d4SStefano Zampini       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
1042785d1243SStefano Zampini       PetscInt          dirsize,i,*is_indices;
10432b095fd8SStefano Zampini       PetscScalar       *array_x;
10442b095fd8SStefano Zampini       const PetscScalar *array_diagonal;
1045785d1243SStefano Zampini 
10463972b0daSStefano Zampini       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
10473972b0daSStefano Zampini       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1048e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1049e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1050e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1051e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
105282ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
10533972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10542b095fd8SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10553972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10562fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
10573972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10582b095fd8SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10593972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1060e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1061e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10628d00608fSStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
10631b968477SStefano Zampini       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
10648efcfb23SStefano Zampini     }
1065a07ea27aSStefano Zampini   }
1066b76ba322SStefano Zampini 
10678efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
10688d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
10698d00608fSStefano Zampini     /* store the original rhs */
10708d00608fSStefano Zampini     if (copy_rhs) {
10718d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10728d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10738d00608fSStefano Zampini     }
10748d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10753972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10763972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10773972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10788efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10797acc28cbSStefano Zampini     if (ksp) {
10807acc28cbSStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
10817acc28cbSStefano Zampini     }
10823308cffdSStefano Zampini   }
10838efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1084b76ba322SStefano Zampini 
1085b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
1086b097fa66SStefano Zampini   if (x && pcbddc->use_exact_dirichlet_trick) {
10878efcfb23SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1088b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1089b76ba322SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1090b76ba322SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
10918efcfb23SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10928efcfb23SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1093a0cb1b98SStefano Zampini     if (ksp) {
1094b76ba322SStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1095b76ba322SStefano Zampini     }
1096b76ba322SStefano Zampini   }
1097b76ba322SStefano Zampini 
1098b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1099906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1100906d46d4SStefano Zampini 
1101906d46d4SStefano Zampini     /* get change ctx */
1102906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1103906d46d4SStefano Zampini 
1104906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1105906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1106906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1107906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1108906d46d4SStefano Zampini 
1109906d46d4SStefano Zampini     /* change iteration matrix */
1110906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1111906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1112906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1113906d46d4SStefano Zampini 
11148d00608fSStefano Zampini     /* store the original rhs */
11158d00608fSStefano Zampini     if (copy_rhs) {
11168d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
11178d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
11188d00608fSStefano Zampini     }
11198d00608fSStefano Zampini 
1120906d46d4SStefano Zampini     /* change rhs */
1121906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1122906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
11238d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
112492e3dcfbSStefano Zampini   }
112592e3dcfbSStefano Zampini 
112692e3dcfbSStefano Zampini   /* remove nullspace if present */
11278efcfb23SStefano Zampini   if (ksp && x && pcbddc->NullSpace) {
11288efcfb23SStefano Zampini     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
11298d00608fSStefano Zampini     /* store the original rhs */
11308d00608fSStefano Zampini     if (copy_rhs) {
11318d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
11328d00608fSStefano Zampini     }
11338d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
1134d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1135b76ba322SStefano Zampini   }
1136534831adSStefano Zampini   PetscFunctionReturn(0);
1137534831adSStefano Zampini }
1138906d46d4SStefano Zampini 
1139534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1140534831adSStefano Zampini #undef __FUNCT__
1141534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1142534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1143534831adSStefano Zampini /*
1144534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1145534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1146534831adSStefano Zampini 
1147534831adSStefano Zampini    Input Parameter:
1148534831adSStefano Zampini +  pc - the preconditioner contex
1149534831adSStefano Zampini 
1150534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1151534831adSStefano Zampini 
1152534831adSStefano Zampini    Notes:
1153534831adSStefano Zampini      The interface routine PCPostSolve() is not usually called directly by
1154534831adSStefano Zampini      the user, but instead is called by KSPSolve().
1155534831adSStefano Zampini */
1156534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1157534831adSStefano Zampini {
1158534831adSStefano Zampini   PetscErrorCode ierr;
1159534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1160534831adSStefano Zampini 
1161534831adSStefano Zampini   PetscFunctionBegin;
1162b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1163906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1164906d46d4SStefano Zampini 
1165906d46d4SStefano Zampini     /* get change ctx */
1166906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1167906d46d4SStefano Zampini 
1168906d46d4SStefano Zampini     /* restore iteration matrix */
1169906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1170906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1171906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1172906d46d4SStefano Zampini 
1173906d46d4SStefano Zampini     /* get solution in original basis */
1174906d46d4SStefano Zampini     if (x) {
1175906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1176906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1177906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
11783425bc38SStefano Zampini     }
1179534831adSStefano Zampini   }
1180906d46d4SStefano Zampini 
11813972b0daSStefano Zampini   /* add solution removed in presolve */
11826bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
11833425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
11843425bc38SStefano Zampini   }
1185906d46d4SStefano Zampini 
1186fb223d50SStefano Zampini   /* restore rhs to its original state */
11878d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
1188fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1189fb223d50SStefano Zampini   }
11908d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
11918efcfb23SStefano Zampini 
11928efcfb23SStefano Zampini   /* restore ksp guess state */
11938efcfb23SStefano Zampini   if (ksp) {
11948efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
11958efcfb23SStefano Zampini   }
1196534831adSStefano Zampini   PetscFunctionReturn(0);
1197534831adSStefano Zampini }
1198534831adSStefano Zampini /* -------------------------------------------------------------------------- */
119953cdbc3dSStefano Zampini #undef __FUNCT__
120053cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
12010c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
12020c7d97c5SJed Brown /*
12030c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
12040c7d97c5SJed Brown                   by setting data structures and options.
12050c7d97c5SJed Brown 
12060c7d97c5SJed Brown    Input Parameter:
120753cdbc3dSStefano Zampini +  pc - the preconditioner context
12080c7d97c5SJed Brown 
12090c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
12100c7d97c5SJed Brown 
12110c7d97c5SJed Brown    Notes:
12120c7d97c5SJed Brown      The interface routine PCSetUp() is not usually called directly by
12130c7d97c5SJed Brown      the user, but instead is called by PCApply() if necessary.
12140c7d97c5SJed Brown */
121553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
12160c7d97c5SJed Brown {
12170c7d97c5SJed Brown   PetscErrorCode ierr;
12180c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
12195e8657edSStefano Zampini   Mat_IS*        matis;
122008122e43SStefano Zampini   MatNullSpace   nearnullspace;
122191e8d312SStefano Zampini   PetscInt       nrows,ncols;
122208122e43SStefano Zampini   PetscBool      computetopography,computesolvers,computesubschurs;
12238de1fae6SStefano Zampini   PetscBool      computeconstraintsmatrix;
12245e8657edSStefano Zampini   PetscBool      new_nearnullspace_provided,ismatis;
12250c7d97c5SJed Brown 
12260c7d97c5SJed Brown   PetscFunctionBegin;
12275e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
12286c4ed002SBarry Smith   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
122991e8d312SStefano Zampini   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
12306c4ed002SBarry Smith   if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
12315e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1232f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
12333b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1234fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1235f4ddd8eeSStefano Zampini   /* split work */
1236674ae819SStefano Zampini   if (pc->setupcalled) {
1237d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1238674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1239674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1240674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1241674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1242674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1243674ae819SStefano Zampini     }
1244674ae819SStefano Zampini   } else {
1245674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1246674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1247674ae819SStefano Zampini   }
1248fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1249fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1250fb180af4SStefano Zampini   }
12518de1fae6SStefano Zampini   computeconstraintsmatrix = PETSC_FALSE;
12526c4ed002SBarry 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");
12535a95e1ceSStefano Zampini   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1254862806e4SStefano Zampini   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1255862806e4SStefano Zampini 
12565a95e1ceSStefano Zampini   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
12576c4ed002SBarry 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");
12586c4ed002SBarry Smith 
1259f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
126070cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
126170cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
126258a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
12631575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);CHKERRQ(ierr);
1264f4ddd8eeSStefano Zampini     }
126558a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1266f4ddd8eeSStefano Zampini   }
1267f4ddd8eeSStefano Zampini 
12685e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
12695e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
12705e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
12715e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
12725e8657edSStefano Zampini   } else {
1273b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
12745e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
12755e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
1276674ae819SStefano Zampini   }
1277f4ddd8eeSStefano Zampini 
1278e496cd5dSStefano Zampini   /* workaround for reals */
1279e496cd5dSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
12803301b35fSStefano Zampini   if (matis->A->symmetric_set) {
12813301b35fSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
12823301b35fSStefano Zampini   }
1283e496cd5dSStefano Zampini #endif
1284e496cd5dSStefano Zampini 
12855e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
12865e8657edSStefano Zampini   {
12875e8657edSStefano Zampini     Mat temp_mat;
12885e8657edSStefano Zampini 
12895e8657edSStefano Zampini     temp_mat = matis->A;
12905e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
12915e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
12925e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
12935e8657edSStefano Zampini     matis->A = temp_mat;
12945e8657edSStefano Zampini   }
1295684f6988SStefano Zampini 
1296b96c3477SStefano Zampini   /* Analyze interface and setup sub_schurs data */
1297674ae819SStefano Zampini   if (computetopography) {
1298674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
12998de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1300674ae819SStefano Zampini   }
1301fb8d54d4SStefano Zampini 
1302b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1303684f6988SStefano Zampini   if (computesolvers) {
1304d5574798SStefano Zampini     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1305d5574798SStefano Zampini 
1306d5574798SStefano Zampini     if (computesubschurs && computetopography) {
130708122e43SStefano Zampini       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1308b1b3d7a2SStefano Zampini     }
1309d5574798SStefano Zampini     if (sub_schurs->use_mumps) {
13102070dbb6SStefano Zampini       if (computesubschurs) {
131108122e43SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
13122070dbb6SStefano Zampini       }
1313d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1314d5574798SStefano Zampini     } else {
1315d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
13162070dbb6SStefano Zampini       if (computesubschurs) {
1317d5574798SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1318d5574798SStefano Zampini       }
13192070dbb6SStefano Zampini     }
132008122e43SStefano Zampini     if (pcbddc->adaptive_selection) {
132108122e43SStefano Zampini       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
13228de1fae6SStefano Zampini       computeconstraintsmatrix = PETSC_TRUE;
1323b7eb3628SStefano Zampini     }
1324b1b3d7a2SStefano Zampini   }
1325684f6988SStefano Zampini 
1326f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1327fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1328f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1329f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1330f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1331f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1332f4ddd8eeSStefano Zampini     } else {
1333f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1334f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1335f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1336165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1337f4ddd8eeSStefano Zampini         PetscInt         i;
1338165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1339165b64e2SStefano Zampini         PetscObjectState state;
1340165b64e2SStefano Zampini         PetscInt         nnsp_size;
1341165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1342f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1343f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1344165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1345f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1346f4ddd8eeSStefano Zampini             break;
1347f4ddd8eeSStefano Zampini           }
1348f4ddd8eeSStefano Zampini         }
1349f4ddd8eeSStefano Zampini       }
1350f4ddd8eeSStefano Zampini     }
1351f4ddd8eeSStefano Zampini   } else {
1352f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1353f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1354f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1355f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1356f4ddd8eeSStefano Zampini     }
1357f4ddd8eeSStefano Zampini   }
1358f4ddd8eeSStefano Zampini 
1359f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1360727cdba6SStefano Zampini   /* reset primal space flags */
1361f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1362727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
13638de1fae6SStefano Zampini   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1364727cdba6SStefano Zampini     /* It also sets the primal space flags */
1365674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1366e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1367f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
13683975b054SStefano Zampini   }
13695e8657edSStefano Zampini 
13703975b054SStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
13715e8657edSStefano Zampini     if (pcbddc->use_change_of_basis) {
13725e8657edSStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
13735e8657edSStefano Zampini 
13745e8657edSStefano Zampini       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
13755e8657edSStefano Zampini       /* get submatrices */
13765e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
13775e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
13785e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
13795e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
13805e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
13815e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
13823975b054SStefano Zampini       /* set flag in pcis to not reuse submatrices during PCISCreate */
13833975b054SStefano Zampini       pcis->reusesubmatrices = PETSC_FALSE;
13845e8657edSStefano Zampini     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1385b96c3477SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
13865e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
13875e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
13885e8657edSStefano Zampini     }
1389b96c3477SStefano Zampini     /* SetUp coarse and local Neumann solvers */
139099cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1391b96c3477SStefano Zampini     /* SetUp Scaling operator */
1392674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
13930c7d97c5SJed Brown   }
139470cf5478SStefano Zampini 
139558a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
139658a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
13972b510759SStefano Zampini   }
13980c7d97c5SJed Brown   PetscFunctionReturn(0);
13990c7d97c5SJed Brown }
14000c7d97c5SJed Brown 
14010c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
14020c7d97c5SJed Brown /*
140350efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
14040c7d97c5SJed Brown 
14050c7d97c5SJed Brown    Input Parameters:
14060f202f7eSStefano Zampini +  pc - the preconditioner context
14070f202f7eSStefano Zampini -  r - input vector (global)
14080c7d97c5SJed Brown 
14090c7d97c5SJed Brown    Output Parameter:
14100c7d97c5SJed Brown .  z - output vector (global)
14110c7d97c5SJed Brown 
14120c7d97c5SJed Brown    Application Interface Routine: PCApply()
14130c7d97c5SJed Brown  */
14140c7d97c5SJed Brown #undef __FUNCT__
14150c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
141653cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
14170c7d97c5SJed Brown {
14180c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
14190c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1420b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
14210c7d97c5SJed Brown   PetscErrorCode    ierr;
14223b03a366Sstefano_zampini   const PetscScalar one = 1.0;
14233b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
14242617d88aSStefano Zampini   const PetscScalar zero = 0.0;
14250c7d97c5SJed Brown 
14260c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
14270c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
1428b097fa66SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
14290c7d97c5SJed Brown 
14300c7d97c5SJed Brown   PetscFunctionBegin;
143185c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1432b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
14330c7d97c5SJed Brown     /* First Dirichlet solve */
14340c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14350c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14360c7d97c5SJed Brown     /*
14370c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1438b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1439674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
14400c7d97c5SJed Brown     */
1441b097fa66SStefano Zampini     if (n_D) {
1442b097fa66SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
14430c7d97c5SJed Brown       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
14448eeda7d8SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1445b097fa66SStefano Zampini       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1446b097fa66SStefano Zampini     } else {
1447b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1448b097fa66SStefano Zampini     }
14490c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14500c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1451674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1452b76ba322SStefano Zampini   } else {
1453b097fa66SStefano Zampini     if (pcbddc->switch_static) {
14540bdf917eSStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1455b097fa66SStefano Zampini     }
1456674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1457b76ba322SStefano Zampini   }
1458b76ba322SStefano Zampini 
14592617d88aSStefano Zampini   /* Apply interface preconditioner
14602617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1461dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
14622617d88aSStefano Zampini 
1463674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1464674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
14650c7d97c5SJed Brown 
14663b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
14670c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14680c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1469b097fa66SStefano Zampini   if (n_B) {
14700c7d97c5SJed Brown     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
14718eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1472b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1473b097fa66SStefano Zampini     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1474b097fa66SStefano Zampini   }
1475df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1476b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1477b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1478b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1479b097fa66SStefano Zampini     } else {
1480b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1481b097fa66SStefano Zampini     }
14820c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14830c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1484b097fa66SStefano Zampini   } else {
1485b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1486b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1487b097fa66SStefano Zampini     } else {
1488b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1489b097fa66SStefano Zampini     }
1490b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1491b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1492b097fa66SStefano Zampini   }
14930c7d97c5SJed Brown   PetscFunctionReturn(0);
14940c7d97c5SJed Brown }
149550efa1b5SStefano Zampini 
149650efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
149750efa1b5SStefano Zampini /*
149850efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
149950efa1b5SStefano Zampini 
150050efa1b5SStefano Zampini    Input Parameters:
15010f202f7eSStefano Zampini +  pc - the preconditioner context
15020f202f7eSStefano Zampini -  r - input vector (global)
150350efa1b5SStefano Zampini 
150450efa1b5SStefano Zampini    Output Parameter:
150550efa1b5SStefano Zampini .  z - output vector (global)
150650efa1b5SStefano Zampini 
150750efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
150850efa1b5SStefano Zampini  */
150950efa1b5SStefano Zampini #undef __FUNCT__
151050efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
151150efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
151250efa1b5SStefano Zampini {
151350efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
151450efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1515b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
151650efa1b5SStefano Zampini   PetscErrorCode    ierr;
151750efa1b5SStefano Zampini   const PetscScalar one = 1.0;
151850efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
151950efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
152050efa1b5SStefano Zampini 
152150efa1b5SStefano Zampini   PetscFunctionBegin;
152250efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1523b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
152450efa1b5SStefano Zampini     /* First Dirichlet solve */
152550efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
152650efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
152750efa1b5SStefano Zampini     /*
152850efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
1529b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
153050efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
153150efa1b5SStefano Zampini     */
1532b097fa66SStefano Zampini     if (n_D) {
1533b097fa66SStefano Zampini       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
153450efa1b5SStefano Zampini       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
153550efa1b5SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1536b097fa66SStefano Zampini       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1537b097fa66SStefano Zampini     } else {
1538b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1539b097fa66SStefano Zampini     }
154050efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
154150efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
154250efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
154350efa1b5SStefano Zampini   } else {
1544b097fa66SStefano Zampini     if (pcbddc->switch_static) {
154550efa1b5SStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1546b097fa66SStefano Zampini     }
154750efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
154850efa1b5SStefano Zampini   }
154950efa1b5SStefano Zampini 
155050efa1b5SStefano Zampini   /* Apply interface preconditioner
155150efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1552dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
155350efa1b5SStefano Zampini 
155450efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
155550efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
155650efa1b5SStefano Zampini 
155750efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
155850efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
155950efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1560b097fa66SStefano Zampini   if (n_B) {
156150efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
156250efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1563b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1564b097fa66SStefano Zampini     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1565b097fa66SStefano Zampini   }
1566b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1567b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1568b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1569b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1570b097fa66SStefano Zampini     } else {
1571b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1572b097fa66SStefano Zampini     }
157350efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
157450efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1575b097fa66SStefano Zampini   } else {
1576b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1577b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1578b097fa66SStefano Zampini     } else {
1579b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1580b097fa66SStefano Zampini     }
1581b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1582b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1583b097fa66SStefano Zampini   }
158450efa1b5SStefano Zampini   PetscFunctionReturn(0);
158550efa1b5SStefano Zampini }
1586da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1587674ae819SStefano Zampini 
1588da1bb401SStefano Zampini #undef __FUNCT__
1589da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1590da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1591da1bb401SStefano Zampini {
1592da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1593da1bb401SStefano Zampini   PetscErrorCode ierr;
1594da1bb401SStefano Zampini 
1595da1bb401SStefano Zampini   PetscFunctionBegin;
1596da1bb401SStefano Zampini   /* free data created by PCIS */
1597da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1598674ae819SStefano Zampini   /* free BDDC custom data  */
1599674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1600674ae819SStefano Zampini   /* destroy objects related to topography */
1601674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1602674ae819SStefano Zampini   /* free allocated graph structure */
1603da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1604b96c3477SStefano Zampini   /* free allocated sub schurs structure */
1605b96c3477SStefano Zampini   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
160634a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1607674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
160834a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1609674ae819SStefano Zampini   /* free solvers stuff */
1610674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
161162a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
161262a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
161362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1614906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1615906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1616906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1617906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1618906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1619906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1620906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1621906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1622906d46d4SStefano Zampini   }
1623906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
16243425bc38SStefano Zampini   /* remove functions */
1625906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1626674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1627bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
16282b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1629b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
16302b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1631bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1632bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
163382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1634bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
163582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1636bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
163782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1638bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1639785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1640bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
164163602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1642bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1643bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1644bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1645bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1646674ae819SStefano Zampini   /* Free the private data structure */
1647674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1648da1bb401SStefano Zampini   PetscFunctionReturn(0);
1649da1bb401SStefano Zampini }
16503425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
16511e6b0712SBarry Smith 
16523425bc38SStefano Zampini #undef __FUNCT__
16533425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
16543425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
16553425bc38SStefano Zampini {
1656674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
1657c08af4c6SStefano Zampini   Vec            copy_standard_rhs;
16583425bc38SStefano Zampini   PC_IS*         pcis;
16593425bc38SStefano Zampini   PC_BDDC*       pcbddc;
16603425bc38SStefano Zampini   PetscErrorCode ierr;
16610c7d97c5SJed Brown 
16623425bc38SStefano Zampini   PetscFunctionBegin;
16633425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16643425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
16653425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
16663425bc38SStefano Zampini 
1667c08af4c6SStefano Zampini   /*
1668c08af4c6SStefano Zampini      change of basis for physical rhs if needed
1669c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
1670c08af4c6SStefano Zampini      TODO: better management when FETIDP will have its own class
1671c08af4c6SStefano Zampini   */
1672c08af4c6SStefano Zampini   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1673c08af4c6SStefano Zampini   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1674c08af4c6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
16753425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
1676c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1677c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1678fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1679fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
1680c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1681c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1682674ae819SStefano Zampini   /* Apply partition of unity */
16833425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1684c08af4c6SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
16858eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
16863425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
16873425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16883425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
16893425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1690c08af4c6SStefano Zampini     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1691c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1692c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1693c08af4c6SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1694c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1695c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16963425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
16973425bc38SStefano Zampini   }
1698c08af4c6SStefano Zampini   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
16993425bc38SStefano Zampini   /* BDDC rhs */
17003425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
17018eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
17023425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17033425bc38SStefano Zampini   }
17043425bc38SStefano Zampini   /* apply BDDC */
1705dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
17063425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
17073425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
17083425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
17093425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17103425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17113425bc38SStefano Zampini   PetscFunctionReturn(0);
17123425bc38SStefano Zampini }
17131e6b0712SBarry Smith 
17143425bc38SStefano Zampini #undef __FUNCT__
17153425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
17163425bc38SStefano Zampini /*@
17170f202f7eSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
17183425bc38SStefano Zampini 
17193425bc38SStefano Zampini    Collective
17203425bc38SStefano Zampini 
17213425bc38SStefano Zampini    Input Parameters:
17220f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
17230f202f7eSStefano Zampini -  standard_rhs    - the right-hand side of the original linear system
17243425bc38SStefano Zampini 
17253425bc38SStefano Zampini    Output Parameters:
17260f202f7eSStefano Zampini .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
17273425bc38SStefano Zampini 
17283425bc38SStefano Zampini    Level: developer
17293425bc38SStefano Zampini 
17303425bc38SStefano Zampini    Notes:
17313425bc38SStefano Zampini 
17320f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
17333425bc38SStefano Zampini @*/
17343425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
17353425bc38SStefano Zampini {
1736674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17373425bc38SStefano Zampini   PetscErrorCode ierr;
17383425bc38SStefano Zampini 
17393425bc38SStefano Zampini   PetscFunctionBegin;
17403425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1741*163d334eSBarry Smith   ierr = PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
17423425bc38SStefano Zampini   PetscFunctionReturn(0);
17433425bc38SStefano Zampini }
17443425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
17451e6b0712SBarry Smith 
17463425bc38SStefano Zampini #undef __FUNCT__
17473425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
17483425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
17493425bc38SStefano Zampini {
1750674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17513425bc38SStefano Zampini   PC_IS*         pcis;
17523425bc38SStefano Zampini   PC_BDDC*       pcbddc;
17533425bc38SStefano Zampini   PetscErrorCode ierr;
17543425bc38SStefano Zampini 
17553425bc38SStefano Zampini   PetscFunctionBegin;
17563425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
17573425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
17583425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
17593425bc38SStefano Zampini 
17603425bc38SStefano Zampini   /* apply B_delta^T */
17613425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17623425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17633425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
17643425bc38SStefano Zampini   /* compute rhs for BDDC application */
17653425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
17668eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
17673425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17683425bc38SStefano Zampini   }
17693425bc38SStefano Zampini   /* apply BDDC */
1770dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
17713425bc38SStefano Zampini   /* put values into standard global vector */
17723425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17733425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17748eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
17753425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
17763425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
17773425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
17783425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17793425bc38SStefano Zampini   }
17803425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17813425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17823425bc38SStefano Zampini   /* final change of basis if needed
17833425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
17843308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
17853425bc38SStefano Zampini   PetscFunctionReturn(0);
17863425bc38SStefano Zampini }
17871e6b0712SBarry Smith 
17883425bc38SStefano Zampini #undef __FUNCT__
17893425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
17903425bc38SStefano Zampini /*@
17910f202f7eSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
17923425bc38SStefano Zampini 
17933425bc38SStefano Zampini    Collective
17943425bc38SStefano Zampini 
17953425bc38SStefano Zampini    Input Parameters:
17960f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
17970f202f7eSStefano Zampini -  fetidp_flux_sol - the solution of the FETI-DP linear system
17983425bc38SStefano Zampini 
17993425bc38SStefano Zampini    Output Parameters:
18000f202f7eSStefano Zampini .  standard_sol    - the solution defined on the physical domain
18013425bc38SStefano Zampini 
18023425bc38SStefano Zampini    Level: developer
18033425bc38SStefano Zampini 
18043425bc38SStefano Zampini    Notes:
18053425bc38SStefano Zampini 
18060f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
18073425bc38SStefano Zampini @*/
18083425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
18093425bc38SStefano Zampini {
1810674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
18113425bc38SStefano Zampini   PetscErrorCode ierr;
18123425bc38SStefano Zampini 
18133425bc38SStefano Zampini   PetscFunctionBegin;
18143425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
1815*163d334eSBarry Smith   ierr = PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
18163425bc38SStefano Zampini   PetscFunctionReturn(0);
18173425bc38SStefano Zampini }
18183425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
18191e6b0712SBarry Smith 
1820f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1821edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1822f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1823f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1824edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1825f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1826674ae819SStefano Zampini 
18273425bc38SStefano Zampini #undef __FUNCT__
18283425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
18293425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
18303425bc38SStefano Zampini {
1831674ae819SStefano Zampini 
1832674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
18333425bc38SStefano Zampini   Mat            newmat;
1834674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
18353425bc38SStefano Zampini   PC             newpc;
1836ce94432eSBarry Smith   MPI_Comm       comm;
18373425bc38SStefano Zampini   PetscErrorCode ierr;
18383425bc38SStefano Zampini 
18393425bc38SStefano Zampini   PetscFunctionBegin;
1840ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
18413425bc38SStefano Zampini   /* FETIDP linear matrix */
18423425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
18433425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
18443425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
18453425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1846edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
18473425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
18483425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
18493425bc38SStefano Zampini   /* FETIDP preconditioner */
18503425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
18513425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
18523425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
18533425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
18543425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
18553425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1856edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
18573425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
185823ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
18593425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
18603425bc38SStefano Zampini   /* return pointers for objects created */
18613425bc38SStefano Zampini   *fetidp_mat=newmat;
18623425bc38SStefano Zampini   *fetidp_pc=newpc;
18633425bc38SStefano Zampini   PetscFunctionReturn(0);
18643425bc38SStefano Zampini }
18651e6b0712SBarry Smith 
18663425bc38SStefano Zampini #undef __FUNCT__
18673425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
18683425bc38SStefano Zampini /*@
18690f202f7eSStefano Zampini  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
18703425bc38SStefano Zampini 
18713425bc38SStefano Zampini    Collective
18723425bc38SStefano Zampini 
18733425bc38SStefano Zampini    Input Parameters:
18740f202f7eSStefano Zampini .  pc - the BDDC preconditioning context (setup should have been called before)
187528509bceSStefano Zampini 
187628509bceSStefano Zampini    Output Parameters:
18770f202f7eSStefano Zampini +  fetidp_mat - shell FETI-DP matrix object
18780f202f7eSStefano Zampini -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
187928509bceSStefano Zampini 
188028509bceSStefano Zampini    Options Database Keys:
18810f202f7eSStefano Zampini .    -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
18823425bc38SStefano Zampini 
18833425bc38SStefano Zampini    Level: developer
18843425bc38SStefano Zampini 
18853425bc38SStefano Zampini    Notes:
18860f202f7eSStefano Zampini      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
18873425bc38SStefano Zampini 
18880f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
18893425bc38SStefano Zampini @*/
18903425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
18913425bc38SStefano Zampini {
18923425bc38SStefano Zampini   PetscErrorCode ierr;
18933425bc38SStefano Zampini 
18943425bc38SStefano Zampini   PetscFunctionBegin;
18953425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18963425bc38SStefano Zampini   if (pc->setupcalled) {
1897516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1898f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
18993425bc38SStefano Zampini   PetscFunctionReturn(0);
19003425bc38SStefano Zampini }
19010c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1902da1bb401SStefano Zampini /*MC
1903da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
19040c7d97c5SJed Brown 
190528509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
190628509bceSStefano Zampini 
190728509bceSStefano Zampini .vb
190828509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
190928509bceSStefano 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
191028509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
19110f202f7eSStefano 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
191228509bceSStefano Zampini .ve
191328509bceSStefano Zampini 
191428509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
191528509bceSStefano Zampini 
19160f202f7eSStefano Zampini    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
191728509bceSStefano Zampini 
191828509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
191928509bceSStefano Zampini 
1920b6fdb6dfSStefano 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.
1921b6fdb6dfSStefano Zampini 
19220f202f7eSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
192328509bceSStefano Zampini 
19240f202f7eSStefano 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()
19250f202f7eSStefano Zampini    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
192628509bceSStefano Zampini 
19270f202f7eSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
192828509bceSStefano Zampini 
19290f202f7eSStefano 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.
19300f202f7eSStefano Zampini    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
193128509bceSStefano Zampini 
19320f202f7eSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
193328509bceSStefano Zampini 
19340f202f7eSStefano 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.
193528509bceSStefano Zampini 
19360f202f7eSStefano 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.
19370f202f7eSStefano Zampini    Deluxe scaling is not supported yet for FETI-DP.
19380f202f7eSStefano Zampini 
19390f202f7eSStefano Zampini    Options Database Keys (some of them, run with -h for a complete list):
19400f202f7eSStefano Zampini 
19410f202f7eSStefano Zampini .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
19420f202f7eSStefano Zampini .    -pc_bddc_use_edges <true> - use or not edges in primal space
19430f202f7eSStefano Zampini .    -pc_bddc_use_faces <false> - use or not faces in primal space
19440f202f7eSStefano Zampini .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
19450f202f7eSStefano Zampini .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
19460f202f7eSStefano Zampini .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
19470f202f7eSStefano Zampini .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
194828509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
19490f202f7eSStefano 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)
19500f202f7eSStefano 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)
19510f202f7eSStefano Zampini .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
19520f202f7eSStefano 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)
19530f202f7eSStefano 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)
195428509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
195528509bceSStefano Zampini 
195628509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
195728509bceSStefano Zampini .vb
195828509bceSStefano Zampini       -pc_bddc_dirichlet_
195928509bceSStefano Zampini       -pc_bddc_neumann_
196028509bceSStefano Zampini       -pc_bddc_coarse_
196128509bceSStefano Zampini .ve
19620f202f7eSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
196328509bceSStefano Zampini 
19640f202f7eSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
196528509bceSStefano Zampini .vb
1966312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1967312be037SStefano Zampini       -pc_bddc_neumann_lN_
1968312be037SStefano Zampini       -pc_bddc_coarse_lN_
196928509bceSStefano Zampini .ve
19700f202f7eSStefano Zampini    Note that level number ranges from the finest (0) to the coarsest (N).
19710f202f7eSStefano 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.
19720f202f7eSStefano Zampini .vb
19730f202f7eSStefano Zampini      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
19740f202f7eSStefano Zampini .ve
19750f202f7eSStefano 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
1976da1bb401SStefano Zampini 
1977da1bb401SStefano Zampini    Level: intermediate
1978da1bb401SStefano Zampini 
1979b6fdb6dfSStefano Zampini    Developer notes:
1980da1bb401SStefano Zampini 
1981da1bb401SStefano Zampini    Contributed by Stefano Zampini
1982da1bb401SStefano Zampini 
1983da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1984da1bb401SStefano Zampini M*/
1985b2573a8aSBarry Smith 
1986da1bb401SStefano Zampini #undef __FUNCT__
1987da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
19888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1989da1bb401SStefano Zampini {
1990da1bb401SStefano Zampini   PetscErrorCode      ierr;
1991da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1992da1bb401SStefano Zampini 
1993da1bb401SStefano Zampini   PetscFunctionBegin;
1994da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1995b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1996da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1997da1bb401SStefano Zampini 
1998da1bb401SStefano Zampini   /* create PCIS data structure */
1999da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
2000da1bb401SStefano Zampini 
200147d04d0dSStefano Zampini   /* BDDC customization */
200208a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
200347d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
200447d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
200547d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
200647d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
200747d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
200847d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
2009fa434743SStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE;
2010fa434743SStefano Zampini   pcbddc->use_qr_single       = PETSC_FALSE;
20113301b35fSStefano Zampini   pcbddc->symmetric_primal    = PETSC_TRUE;
201247d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
2013b9d89cd5SStefano Zampini   /* private */
2014727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
20150e6343abSStefano Zampini   pcbddc->local_primal_size_cc       = 0;
20160e6343abSStefano Zampini   pcbddc->local_primal_ref_node      = 0;
20170e6343abSStefano Zampini   pcbddc->local_primal_ref_mult      = 0;
2018e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
2019727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
2020fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
202168457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
2022f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
2023727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
2024f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
2025f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
2026f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
2027674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
20280bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
20293972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
2030534831adSStefano Zampini   pcbddc->original_rhs               = 0;
2031534831adSStefano Zampini   pcbddc->local_mat                  = 0;
2032534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
2033b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
2034906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
2035da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
2036da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
2037da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
2038da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
203915aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
204015aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
2041da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
2042da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
2043da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
2044da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
2045da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
2046da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
2047da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
2048da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
2049da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
2050da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
2051785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
2052785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
2053785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
205460d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
205560d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
205663602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
2057da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
205863602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
2059da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
206085c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
206147d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
206247d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
2063b0c7d250SStefano Zampini   pcbddc->coarse_adj_red             = 0;
20644fad6a16SStefano Zampini   pcbddc->current_level              = 0;
20652b510759SStefano Zampini   pcbddc->max_levels                 = 0;
2066323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
206774e2c79eSStefano Zampini   pcbddc->redistribute_coarse        = 0;
2068f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
2069323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
2070da1bb401SStefano Zampini 
2071674ae819SStefano Zampini   /* create local graph structure */
2072674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2073674ae819SStefano Zampini 
2074674ae819SStefano Zampini   /* scaling */
2075674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
207634a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2077ac632422SStefano Zampini   pcbddc->faster_deluxe         = PETSC_FALSE;
2078b96c3477SStefano Zampini 
2079b96c3477SStefano Zampini   /* create sub schurs structure */
2080b96c3477SStefano Zampini   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2081b96c3477SStefano Zampini   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2082b96c3477SStefano Zampini   pcbddc->sub_schurs_layers      = -1;
2083b96c3477SStefano Zampini   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2084b96c3477SStefano Zampini 
2085b96c3477SStefano Zampini   pcbddc->computed_rowadj = PETSC_FALSE;
2086da1bb401SStefano Zampini 
2087b7eb3628SStefano Zampini   /* adaptivity */
2088f6f667cfSStefano Zampini   pcbddc->adaptive_threshold      = 0.0;
208908122e43SStefano Zampini   pcbddc->adaptive_nmax           = 0;
2090f6f667cfSStefano Zampini   pcbddc->adaptive_nmin           = 0;
2091b7eb3628SStefano Zampini 
2092da1bb401SStefano Zampini   /* function pointers */
2093da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
209493bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2095da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
2096da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
2097da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
20986b78500eSPatrick Sanan   pc->ops->view                = PCView_BDDC;
2099da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
2100da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
2101da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
2102534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
2103534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
2104da1bb401SStefano Zampini 
2105da1bb401SStefano Zampini   /* composing function */
2106906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2107674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2108bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
21092b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2110b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
21112b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2112bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2113bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
211482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2115bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
211682ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2117bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
211882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2119bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
212082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2121bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
212263602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2123bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2124bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2125bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2126bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2127da1bb401SStefano Zampini   PetscFunctionReturn(0);
2128da1bb401SStefano Zampini }
21293425bc38SStefano Zampini 
2130