xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 30368db7a196d34abcffce59efe277e055a7d266)
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 
290369aaf7SStefano Zampini /* temporarily declare it */
300369aaf7SStefano Zampini PetscErrorCode PCApply_BDDC(PC,Vec,Vec);
310369aaf7SStefano Zampini 
320c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
330c7d97c5SJed Brown #undef __FUNCT__
340c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
358c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc)
360c7d97c5SJed Brown {
370c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
380c7d97c5SJed Brown   PetscErrorCode ierr;
390c7d97c5SJed Brown 
400c7d97c5SJed Brown   PetscFunctionBegin;
41e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
428eeda7d8SStefano Zampini   /* Verbose debugging */
438eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
448eeda7d8SStefano Zampini   /* Primal space cumstomization */
4508a5cf49SStefano 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);
468eeda7d8SStefano 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);
478eeda7d8SStefano 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);
488eeda7d8SStefano 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);
496661aa1dSStefano 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);
50ac632422SStefano 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);
518eeda7d8SStefano Zampini   /* Change of basis */
52b9b85e73SStefano 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);
53b9b85e73SStefano 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);
54674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
55674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
56674ae819SStefano Zampini   }
578eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
588eeda7d8SStefano 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);
5974e2c79eSStefano 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);
600298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
612b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
62323d395dSStefano 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);
63674ae819SStefano 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);
64b96c3477SStefano 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);
65b96c3477SStefano 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);
66b96c3477SStefano 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);
67ac632422SStefano 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);
684c6709b3SStefano 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);
6908122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
7008122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
713301b35fSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
72b0c7d250SStefano 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);
7306a4e24aSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_benign_trick","Apply the benign subspace trick to a class of saddle point problems","none",pcbddc->benign_saddle_point,&pcbddc->benign_saddle_point,NULL);CHKERRQ(ierr);
744f1b2e48SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_detect_disconnected","Detects disconnected subdomains","none",pcbddc->detect_disconnected,&pcbddc->detect_disconnected,NULL);CHKERRQ(ierr);
750c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
760c7d97c5SJed Brown   PetscFunctionReturn(0);
770c7d97c5SJed Brown }
780c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
79674ae819SStefano Zampini #undef __FUNCT__
80906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
81906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
82b9b85e73SStefano Zampini {
83b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
84b9b85e73SStefano Zampini   PetscErrorCode ierr;
85b9b85e73SStefano Zampini 
86b9b85e73SStefano Zampini   PetscFunctionBegin;
87b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
88b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
89b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
90b9b85e73SStefano Zampini   PetscFunctionReturn(0);
91b9b85e73SStefano Zampini }
92b9b85e73SStefano Zampini #undef __FUNCT__
93906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
94b9b85e73SStefano Zampini /*@
95906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
96b9b85e73SStefano Zampini 
97b9b85e73SStefano Zampini    Collective on PC
98b9b85e73SStefano Zampini 
99b9b85e73SStefano Zampini    Input Parameters:
100b9b85e73SStefano Zampini +  pc - the preconditioning context
101906d46d4SStefano Zampini -  change - the change of basis matrix
102b9b85e73SStefano Zampini 
103b9b85e73SStefano Zampini    Level: intermediate
104b9b85e73SStefano Zampini 
105b9b85e73SStefano Zampini    Notes:
106b9b85e73SStefano Zampini 
107b9b85e73SStefano Zampini .seealso: PCBDDC
108b9b85e73SStefano Zampini @*/
109906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
110b9b85e73SStefano Zampini {
111b9b85e73SStefano Zampini   PetscErrorCode ierr;
112b9b85e73SStefano Zampini 
113b9b85e73SStefano Zampini   PetscFunctionBegin;
114b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
116906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
117906d46d4SStefano Zampini   if (pc->mat) {
118906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
119906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
120906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
121906d46d4SStefano Zampini     if (rows_c != rows) {
122906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
123906d46d4SStefano Zampini     }
124906d46d4SStefano Zampini     if (cols_c != cols) {
125906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
126906d46d4SStefano Zampini     }
127906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
128906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
129906d46d4SStefano Zampini     if (rows_c != rows) {
130906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
131906d46d4SStefano Zampini     }
132906d46d4SStefano Zampini     if (cols_c != cols) {
133906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
134906d46d4SStefano Zampini     }
135906d46d4SStefano Zampini   }
136906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
137b9b85e73SStefano Zampini   PetscFunctionReturn(0);
138b9b85e73SStefano Zampini }
139b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
140b9b85e73SStefano Zampini #undef __FUNCT__
141*30368db7SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesIS_BDDC"
142*30368db7SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices)
143*30368db7SStefano Zampini {
144*30368db7SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
145*30368db7SStefano Zampini   PetscErrorCode ierr;
146*30368db7SStefano Zampini 
147*30368db7SStefano Zampini   PetscFunctionBegin;
148*30368db7SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
149*30368db7SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices_local);CHKERRQ(ierr);
150*30368db7SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
151*30368db7SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
152*30368db7SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
153*30368db7SStefano Zampini   PetscFunctionReturn(0);
154*30368db7SStefano Zampini }
155*30368db7SStefano Zampini #undef __FUNCT__
156*30368db7SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesIS"
157*30368db7SStefano Zampini /*@
158*30368db7SStefano Zampini  PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in PCBDDC
159*30368db7SStefano Zampini 
160*30368db7SStefano Zampini    Collective
161*30368db7SStefano Zampini 
162*30368db7SStefano Zampini    Input Parameters:
163*30368db7SStefano Zampini +  pc - the preconditioning context
164*30368db7SStefano Zampini -  PrimalVertices - index set of primal vertices in global numbering (can be empty)
165*30368db7SStefano Zampini 
166*30368db7SStefano Zampini    Level: intermediate
167*30368db7SStefano Zampini 
168*30368db7SStefano Zampini    Notes:
169*30368db7SStefano Zampini      Any process can list any global node
170*30368db7SStefano Zampini 
171*30368db7SStefano Zampini .seealso: PCBDDC
172*30368db7SStefano Zampini @*/
173*30368db7SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices)
174*30368db7SStefano Zampini {
175*30368db7SStefano Zampini   PetscErrorCode ierr;
176*30368db7SStefano Zampini 
177*30368db7SStefano Zampini   PetscFunctionBegin;
178*30368db7SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
179*30368db7SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
180*30368db7SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
181*30368db7SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
182*30368db7SStefano Zampini   PetscFunctionReturn(0);
183*30368db7SStefano Zampini }
184*30368db7SStefano Zampini /* -------------------------------------------------------------------------- */
185*30368db7SStefano Zampini #undef __FUNCT__
186674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
187674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
188674ae819SStefano Zampini {
189674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
190674ae819SStefano Zampini   PetscErrorCode ierr;
1911e6b0712SBarry Smith 
192674ae819SStefano Zampini   PetscFunctionBegin;
193674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
194*30368db7SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices_local);CHKERRQ(ierr);
195674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
196*30368db7SStefano Zampini   pcbddc->user_primal_vertices_local = PrimalVertices;
197*30368db7SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
198674ae819SStefano Zampini   PetscFunctionReturn(0);
199674ae819SStefano Zampini }
200674ae819SStefano Zampini #undef __FUNCT__
201674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
202674ae819SStefano Zampini /*@
20328509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
204674ae819SStefano Zampini 
20517eb1463SStefano Zampini    Collective
206674ae819SStefano Zampini 
207674ae819SStefano Zampini    Input Parameters:
208674ae819SStefano Zampini +  pc - the preconditioning context
20917eb1463SStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
210674ae819SStefano Zampini 
211674ae819SStefano Zampini    Level: intermediate
212674ae819SStefano Zampini 
213674ae819SStefano Zampini    Notes:
214674ae819SStefano Zampini 
215674ae819SStefano Zampini .seealso: PCBDDC
216674ae819SStefano Zampini @*/
217674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
218674ae819SStefano Zampini {
219674ae819SStefano Zampini   PetscErrorCode ierr;
220674ae819SStefano Zampini 
221674ae819SStefano Zampini   PetscFunctionBegin;
222674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
223674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
22417eb1463SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
225674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
226674ae819SStefano Zampini   PetscFunctionReturn(0);
227674ae819SStefano Zampini }
228674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
2290c7d97c5SJed Brown #undef __FUNCT__
2304fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
2314fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
2324fad6a16SStefano Zampini {
2334fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2344fad6a16SStefano Zampini 
2354fad6a16SStefano Zampini   PetscFunctionBegin;
2364fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
2374fad6a16SStefano Zampini   PetscFunctionReturn(0);
2384fad6a16SStefano Zampini }
2391e6b0712SBarry Smith 
2404fad6a16SStefano Zampini #undef __FUNCT__
2414fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
2424fad6a16SStefano Zampini /*@
24328509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
2444fad6a16SStefano Zampini 
2454fad6a16SStefano Zampini    Logically collective on PC
2464fad6a16SStefano Zampini 
2474fad6a16SStefano Zampini    Input Parameters:
2484fad6a16SStefano Zampini +  pc - the preconditioning context
24928509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2504fad6a16SStefano Zampini 
2510f202f7eSStefano Zampini    Options Database Keys:
2520f202f7eSStefano Zampini .    -pc_bddc_coarsening_ratio
2534fad6a16SStefano Zampini 
2544fad6a16SStefano Zampini    Level: intermediate
2554fad6a16SStefano Zampini 
2564fad6a16SStefano Zampini    Notes:
2570f202f7eSStefano Zampini      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2584fad6a16SStefano Zampini 
2590f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetLevels()
2604fad6a16SStefano Zampini @*/
2614fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2624fad6a16SStefano Zampini {
2634fad6a16SStefano Zampini   PetscErrorCode ierr;
2644fad6a16SStefano Zampini 
2654fad6a16SStefano Zampini   PetscFunctionBegin;
2664fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2672b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2684fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2694fad6a16SStefano Zampini   PetscFunctionReturn(0);
2704fad6a16SStefano Zampini }
2712b510759SStefano Zampini 
272b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2732b510759SStefano Zampini #undef __FUNCT__
274b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
275b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
276b8ffe317SStefano Zampini {
277b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
278b8ffe317SStefano Zampini 
279b8ffe317SStefano Zampini   PetscFunctionBegin;
28085c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
281b8ffe317SStefano Zampini   PetscFunctionReturn(0);
282b8ffe317SStefano Zampini }
283b8ffe317SStefano Zampini 
284b8ffe317SStefano Zampini #undef __FUNCT__
285b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
286b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2872b510759SStefano Zampini {
2882b510759SStefano Zampini   PetscErrorCode ierr;
2892b510759SStefano Zampini 
2902b510759SStefano Zampini   PetscFunctionBegin;
2912b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
292b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
293b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2942b510759SStefano Zampini   PetscFunctionReturn(0);
2952b510759SStefano Zampini }
2961e6b0712SBarry Smith 
2974fad6a16SStefano Zampini #undef __FUNCT__
2982b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2992b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
3004fad6a16SStefano Zampini {
3014fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3024fad6a16SStefano Zampini 
3034fad6a16SStefano Zampini   PetscFunctionBegin;
3042b510759SStefano Zampini   pcbddc->current_level = level;
3054fad6a16SStefano Zampini   PetscFunctionReturn(0);
3064fad6a16SStefano Zampini }
3071e6b0712SBarry Smith 
3084fad6a16SStefano Zampini #undef __FUNCT__
309b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
310b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
311b8ffe317SStefano Zampini {
312b8ffe317SStefano Zampini   PetscErrorCode ierr;
313b8ffe317SStefano Zampini 
314b8ffe317SStefano Zampini   PetscFunctionBegin;
315b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
316b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
317b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
318b8ffe317SStefano Zampini   PetscFunctionReturn(0);
319b8ffe317SStefano Zampini }
320b8ffe317SStefano Zampini 
321b8ffe317SStefano Zampini #undef __FUNCT__
3222b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
3232b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
3242b510759SStefano Zampini {
3252b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3262b510759SStefano Zampini 
3272b510759SStefano Zampini   PetscFunctionBegin;
3282b510759SStefano Zampini   pcbddc->max_levels = levels;
3292b510759SStefano Zampini   PetscFunctionReturn(0);
3302b510759SStefano Zampini }
3312b510759SStefano Zampini 
3322b510759SStefano Zampini #undef __FUNCT__
3332b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
3344fad6a16SStefano Zampini /*@
33528509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
3364fad6a16SStefano Zampini 
3374fad6a16SStefano Zampini    Logically collective on PC
3384fad6a16SStefano Zampini 
3394fad6a16SStefano Zampini    Input Parameters:
3404fad6a16SStefano Zampini +  pc - the preconditioning context
34128509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
3424fad6a16SStefano Zampini 
3430f202f7eSStefano Zampini    Options Database Keys:
3440f202f7eSStefano Zampini .    -pc_bddc_levels
3454fad6a16SStefano Zampini 
3464fad6a16SStefano Zampini    Level: intermediate
3474fad6a16SStefano Zampini 
3484fad6a16SStefano Zampini    Notes:
3490f202f7eSStefano Zampini      Default value is 0, i.e. traditional one-level BDDC
3504fad6a16SStefano Zampini 
3510f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
3524fad6a16SStefano Zampini @*/
3532b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3544fad6a16SStefano Zampini {
3554fad6a16SStefano Zampini   PetscErrorCode ierr;
3564fad6a16SStefano Zampini 
3574fad6a16SStefano Zampini   PetscFunctionBegin;
3584fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3592b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
360312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3612b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3624fad6a16SStefano Zampini   PetscFunctionReturn(0);
3634fad6a16SStefano Zampini }
3644fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3651e6b0712SBarry Smith 
3664fad6a16SStefano Zampini #undef __FUNCT__
3670bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3680bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3690bdf917eSStefano Zampini {
3700bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3710bdf917eSStefano Zampini   PetscErrorCode ierr;
3720bdf917eSStefano Zampini 
3730bdf917eSStefano Zampini   PetscFunctionBegin;
3740bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3750bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3760bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3770bdf917eSStefano Zampini   PetscFunctionReturn(0);
3780bdf917eSStefano Zampini }
3791e6b0712SBarry Smith 
3800bdf917eSStefano Zampini #undef __FUNCT__
3810bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3820bdf917eSStefano Zampini /*@
38328509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3840bdf917eSStefano Zampini 
3850bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3860bdf917eSStefano Zampini 
3870bdf917eSStefano Zampini    Input Parameters:
3880bdf917eSStefano Zampini +  pc - the preconditioning context
38928509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3900bdf917eSStefano Zampini 
3910bdf917eSStefano Zampini    Level: intermediate
3920bdf917eSStefano Zampini 
3930bdf917eSStefano Zampini    Notes:
3940bdf917eSStefano Zampini 
3950bdf917eSStefano Zampini .seealso: PCBDDC
3960bdf917eSStefano Zampini @*/
3970bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3980bdf917eSStefano Zampini {
3990bdf917eSStefano Zampini   PetscErrorCode ierr;
4000bdf917eSStefano Zampini 
4010bdf917eSStefano Zampini   PetscFunctionBegin;
4020bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
403674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
4042b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
4050bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
4060bdf917eSStefano Zampini   PetscFunctionReturn(0);
4070bdf917eSStefano Zampini }
4080bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
4091e6b0712SBarry Smith 
4100bdf917eSStefano Zampini #undef __FUNCT__
4113b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
4123b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
4133b03a366Sstefano_zampini {
4143b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4153b03a366Sstefano_zampini   PetscErrorCode ierr;
4163b03a366Sstefano_zampini 
4173b03a366Sstefano_zampini   PetscFunctionBegin;
418785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
419785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4203b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
42136e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
42236e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
423fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4243b03a366Sstefano_zampini   PetscFunctionReturn(0);
4253b03a366Sstefano_zampini }
4261e6b0712SBarry Smith 
4273b03a366Sstefano_zampini #undef __FUNCT__
4283b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
4293b03a366Sstefano_zampini /*@
43028509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
4313b03a366Sstefano_zampini 
432785d1243SStefano Zampini    Collective
4333b03a366Sstefano_zampini 
4343b03a366Sstefano_zampini    Input Parameters:
4353b03a366Sstefano_zampini +  pc - the preconditioning context
436785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
4373b03a366Sstefano_zampini 
4383b03a366Sstefano_zampini    Level: intermediate
4393b03a366Sstefano_zampini 
4400f202f7eSStefano Zampini    Notes:
4410f202f7eSStefano Zampini      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
4423b03a366Sstefano_zampini 
4430f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
4443b03a366Sstefano_zampini @*/
4453b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
4463b03a366Sstefano_zampini {
4473b03a366Sstefano_zampini   PetscErrorCode ierr;
4483b03a366Sstefano_zampini 
4493b03a366Sstefano_zampini   PetscFunctionBegin;
4503b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
451674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
452785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4533b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4543b03a366Sstefano_zampini   PetscFunctionReturn(0);
4553b03a366Sstefano_zampini }
4563b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4571e6b0712SBarry Smith 
4583b03a366Sstefano_zampini #undef __FUNCT__
45982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
46082ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4610c7d97c5SJed Brown {
4620c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4630c7d97c5SJed Brown   PetscErrorCode ierr;
4640c7d97c5SJed Brown 
4650c7d97c5SJed Brown   PetscFunctionBegin;
466785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
467785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4680c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4690c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
470785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4717d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4720c7d97c5SJed Brown   PetscFunctionReturn(0);
4730c7d97c5SJed Brown }
4740c7d97c5SJed Brown 
4750c7d97c5SJed Brown #undef __FUNCT__
47682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4779c0446d6SStefano Zampini /*@
47882ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4799c0446d6SStefano Zampini 
480785d1243SStefano Zampini    Collective
48153cdbc3dSStefano Zampini 
48253cdbc3dSStefano Zampini    Input Parameters:
48353cdbc3dSStefano Zampini +  pc - the preconditioning context
48482ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
48553cdbc3dSStefano Zampini 
48653cdbc3dSStefano Zampini    Level: intermediate
48753cdbc3dSStefano Zampini 
4889c0446d6SStefano Zampini    Notes:
48953cdbc3dSStefano Zampini 
4900f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
49153cdbc3dSStefano Zampini @*/
49282ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4930c7d97c5SJed Brown {
4940c7d97c5SJed Brown   PetscErrorCode ierr;
4950c7d97c5SJed Brown 
4960c7d97c5SJed Brown   PetscFunctionBegin;
4970c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4980c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
49982ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
50082ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
5010c7d97c5SJed Brown   PetscFunctionReturn(0);
5020c7d97c5SJed Brown }
5030c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
5040c7d97c5SJed Brown 
5050c7d97c5SJed Brown #undef __FUNCT__
5060c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
50753cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
5080c7d97c5SJed Brown {
5090c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
51053cdbc3dSStefano Zampini   PetscErrorCode ierr;
5110c7d97c5SJed Brown 
5120c7d97c5SJed Brown   PetscFunctionBegin;
513785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
514785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
51553cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
51636e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
51736e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
518fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5190c7d97c5SJed Brown   PetscFunctionReturn(0);
5200c7d97c5SJed Brown }
5211e6b0712SBarry Smith 
5220c7d97c5SJed Brown #undef __FUNCT__
5230c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
52457527edcSJed Brown /*@
52528509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
52657527edcSJed Brown 
527785d1243SStefano Zampini    Collective
52857527edcSJed Brown 
52957527edcSJed Brown    Input Parameters:
53057527edcSJed Brown +  pc - the preconditioning context
531785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
53257527edcSJed Brown 
53357527edcSJed Brown    Level: intermediate
53457527edcSJed Brown 
5350f202f7eSStefano Zampini    Notes:
5360f202f7eSStefano Zampini      Any process can list any global node
53757527edcSJed Brown 
5380f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
53957527edcSJed Brown @*/
54053cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
5410c7d97c5SJed Brown {
5420c7d97c5SJed Brown   PetscErrorCode ierr;
5430c7d97c5SJed Brown 
5440c7d97c5SJed Brown   PetscFunctionBegin;
5450c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
546674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
547785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54853cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
54953cdbc3dSStefano Zampini   PetscFunctionReturn(0);
55053cdbc3dSStefano Zampini }
55153cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5521e6b0712SBarry Smith 
55353cdbc3dSStefano Zampini #undef __FUNCT__
55482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
55582ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5560c7d97c5SJed Brown {
5570c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5580c7d97c5SJed Brown   PetscErrorCode ierr;
5590c7d97c5SJed Brown 
5600c7d97c5SJed Brown   PetscFunctionBegin;
561785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
562785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5630c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5640c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
565785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5667d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5670c7d97c5SJed Brown   PetscFunctionReturn(0);
5680c7d97c5SJed Brown }
5690c7d97c5SJed Brown 
5700c7d97c5SJed Brown #undef __FUNCT__
57182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5720c7d97c5SJed Brown /*@
57382ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5740c7d97c5SJed Brown 
575785d1243SStefano Zampini    Collective
5760c7d97c5SJed Brown 
5770c7d97c5SJed Brown    Input Parameters:
5780c7d97c5SJed Brown +  pc - the preconditioning context
57982ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5800c7d97c5SJed Brown 
5810c7d97c5SJed Brown    Level: intermediate
5820c7d97c5SJed Brown 
5830c7d97c5SJed Brown    Notes:
5840c7d97c5SJed Brown 
5850f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
5860c7d97c5SJed Brown @*/
58782ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5880c7d97c5SJed Brown {
5890c7d97c5SJed Brown   PetscErrorCode ierr;
5900c7d97c5SJed Brown 
5910c7d97c5SJed Brown   PetscFunctionBegin;
5920c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5930c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
59482ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
59582ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
59653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
59753cdbc3dSStefano Zampini }
59853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
59953cdbc3dSStefano Zampini 
60053cdbc3dSStefano Zampini #undef __FUNCT__
601da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
602da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
603da1bb401SStefano Zampini {
604da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
605da1bb401SStefano Zampini 
606da1bb401SStefano Zampini   PetscFunctionBegin;
607da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
608da1bb401SStefano Zampini   PetscFunctionReturn(0);
609da1bb401SStefano Zampini }
6101e6b0712SBarry Smith 
611da1bb401SStefano Zampini #undef __FUNCT__
612da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
613da1bb401SStefano Zampini /*@
614785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
615da1bb401SStefano Zampini 
616785d1243SStefano Zampini    Collective
617785d1243SStefano Zampini 
618785d1243SStefano Zampini    Input Parameters:
619785d1243SStefano Zampini .  pc - the preconditioning context
620785d1243SStefano Zampini 
621785d1243SStefano Zampini    Output Parameters:
622785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
623785d1243SStefano Zampini 
624785d1243SStefano Zampini    Level: intermediate
625785d1243SStefano Zampini 
6260f202f7eSStefano Zampini    Notes:
6270f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
628785d1243SStefano Zampini 
629785d1243SStefano Zampini .seealso: PCBDDC
630785d1243SStefano Zampini @*/
631785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
632785d1243SStefano Zampini {
633785d1243SStefano Zampini   PetscErrorCode ierr;
634785d1243SStefano Zampini 
635785d1243SStefano Zampini   PetscFunctionBegin;
636785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
637785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
638785d1243SStefano Zampini   PetscFunctionReturn(0);
639785d1243SStefano Zampini }
640785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
641785d1243SStefano Zampini 
642785d1243SStefano Zampini #undef __FUNCT__
643785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
644785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
645785d1243SStefano Zampini {
646785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
647785d1243SStefano Zampini 
648785d1243SStefano Zampini   PetscFunctionBegin;
649785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
650785d1243SStefano Zampini   PetscFunctionReturn(0);
651785d1243SStefano Zampini }
652785d1243SStefano Zampini 
653785d1243SStefano Zampini #undef __FUNCT__
65482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
655da1bb401SStefano Zampini /*@
65682ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
657da1bb401SStefano Zampini 
658785d1243SStefano Zampini    Collective
659da1bb401SStefano Zampini 
660da1bb401SStefano Zampini    Input Parameters:
66128509bceSStefano Zampini .  pc - the preconditioning context
662da1bb401SStefano Zampini 
663da1bb401SStefano Zampini    Output Parameters:
66428509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
665da1bb401SStefano Zampini 
666da1bb401SStefano Zampini    Level: intermediate
667da1bb401SStefano Zampini 
668da1bb401SStefano Zampini    Notes:
6690f202f7eSStefano 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).
6700f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
671da1bb401SStefano Zampini 
672da1bb401SStefano Zampini .seealso: PCBDDC
673da1bb401SStefano Zampini @*/
67482ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
675da1bb401SStefano Zampini {
676da1bb401SStefano Zampini   PetscErrorCode ierr;
677da1bb401SStefano Zampini 
678da1bb401SStefano Zampini   PetscFunctionBegin;
679da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
68082ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
681da1bb401SStefano Zampini   PetscFunctionReturn(0);
682da1bb401SStefano Zampini }
683da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6841e6b0712SBarry Smith 
685da1bb401SStefano Zampini #undef __FUNCT__
68653cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
68753cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
68853cdbc3dSStefano Zampini {
68953cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
69053cdbc3dSStefano Zampini 
69153cdbc3dSStefano Zampini   PetscFunctionBegin;
69253cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
69353cdbc3dSStefano Zampini   PetscFunctionReturn(0);
69453cdbc3dSStefano Zampini }
6951e6b0712SBarry Smith 
69653cdbc3dSStefano Zampini #undef __FUNCT__
69753cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
69853cdbc3dSStefano Zampini /*@
699785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
70053cdbc3dSStefano Zampini 
701785d1243SStefano Zampini    Collective
702785d1243SStefano Zampini 
703785d1243SStefano Zampini    Input Parameters:
704785d1243SStefano Zampini .  pc - the preconditioning context
705785d1243SStefano Zampini 
706785d1243SStefano Zampini    Output Parameters:
707785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
708785d1243SStefano Zampini 
709785d1243SStefano Zampini    Level: intermediate
710785d1243SStefano Zampini 
7110f202f7eSStefano Zampini    Notes:
7120f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
713785d1243SStefano Zampini 
714785d1243SStefano Zampini .seealso: PCBDDC
715785d1243SStefano Zampini @*/
716785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
717785d1243SStefano Zampini {
718785d1243SStefano Zampini   PetscErrorCode ierr;
719785d1243SStefano Zampini 
720785d1243SStefano Zampini   PetscFunctionBegin;
721785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
722785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
723785d1243SStefano Zampini   PetscFunctionReturn(0);
724785d1243SStefano Zampini }
725785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
726785d1243SStefano Zampini 
727785d1243SStefano Zampini #undef __FUNCT__
728785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
729785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
730785d1243SStefano Zampini {
731785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
732785d1243SStefano Zampini 
733785d1243SStefano Zampini   PetscFunctionBegin;
734785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
735785d1243SStefano Zampini   PetscFunctionReturn(0);
736785d1243SStefano Zampini }
737785d1243SStefano Zampini 
738785d1243SStefano Zampini #undef __FUNCT__
73982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
74053cdbc3dSStefano Zampini /*@
74182ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
74253cdbc3dSStefano Zampini 
743785d1243SStefano Zampini    Collective
74453cdbc3dSStefano Zampini 
74553cdbc3dSStefano Zampini    Input Parameters:
74628509bceSStefano Zampini .  pc - the preconditioning context
74753cdbc3dSStefano Zampini 
74853cdbc3dSStefano Zampini    Output Parameters:
74928509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
75053cdbc3dSStefano Zampini 
75153cdbc3dSStefano Zampini    Level: intermediate
75253cdbc3dSStefano Zampini 
75353cdbc3dSStefano Zampini    Notes:
7540f202f7eSStefano 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).
7550f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
75653cdbc3dSStefano Zampini 
75753cdbc3dSStefano Zampini .seealso: PCBDDC
75853cdbc3dSStefano Zampini @*/
75982ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
76053cdbc3dSStefano Zampini {
76153cdbc3dSStefano Zampini   PetscErrorCode ierr;
76253cdbc3dSStefano Zampini 
76353cdbc3dSStefano Zampini   PetscFunctionBegin;
76453cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
76582ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7660c7d97c5SJed Brown   PetscFunctionReturn(0);
7670c7d97c5SJed Brown }
76836e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7691e6b0712SBarry Smith 
77036e030ebSStefano Zampini #undef __FUNCT__
771da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7721a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
77336e030ebSStefano Zampini {
77436e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
775da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
776da1bb401SStefano Zampini   PetscErrorCode ierr;
77736e030ebSStefano Zampini 
77836e030ebSStefano Zampini   PetscFunctionBegin;
779674ae819SStefano Zampini   /* free old CSR */
780674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
781fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
782674ae819SStefano Zampini   /* get CSR into graph structure */
783da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
784854ce69bSBarry Smith     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
785785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
786674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
787674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
788da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7891a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7901a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
791674ae819SStefano Zampini   }
792575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
79336e030ebSStefano Zampini   PetscFunctionReturn(0);
79436e030ebSStefano Zampini }
7951e6b0712SBarry Smith 
79636e030ebSStefano Zampini #undef __FUNCT__
797da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
79836e030ebSStefano Zampini /*@
7990f202f7eSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
80036e030ebSStefano Zampini 
80136e030ebSStefano Zampini    Not collective
80236e030ebSStefano Zampini 
80336e030ebSStefano Zampini    Input Parameters:
80436e030ebSStefano Zampini +  pc - the preconditioning context
8050f202f7eSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
80628509bceSStefano Zampini .  xadj, adjncy - the CSR graph
80728509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
80836e030ebSStefano Zampini 
80936e030ebSStefano Zampini    Level: intermediate
81036e030ebSStefano Zampini 
81136e030ebSStefano Zampini    Notes:
81236e030ebSStefano Zampini 
81328509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
81436e030ebSStefano Zampini @*/
8151a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
81636e030ebSStefano Zampini {
817575ad6abSStefano Zampini   void (*f)(void) = 0;
81836e030ebSStefano Zampini   PetscErrorCode ierr;
81936e030ebSStefano Zampini 
82036e030ebSStefano Zampini   PetscFunctionBegin;
82136e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
822674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
823d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
824674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
825674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
826da1bb401SStefano Zampini   }
82736e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
828575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
829575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
830575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
831575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
832575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
83336e030ebSStefano Zampini   }
83436e030ebSStefano Zampini   PetscFunctionReturn(0);
83536e030ebSStefano Zampini }
8369c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
8371e6b0712SBarry Smith 
8389c0446d6SStefano Zampini #undef __FUNCT__
83963602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
84063602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
84163602bcaSStefano Zampini {
84263602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
84363602bcaSStefano Zampini   PetscInt i;
84463602bcaSStefano Zampini   PetscErrorCode ierr;
84563602bcaSStefano Zampini 
84663602bcaSStefano Zampini   PetscFunctionBegin;
84763602bcaSStefano Zampini   /* Destroy ISes if they were already set */
84863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
84963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
85063602bcaSStefano Zampini   }
85163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
85263602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
85363602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
85463602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
85563602bcaSStefano Zampini   }
85663602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
85763602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
85863602bcaSStefano Zampini   /* allocate space then set */
859d02579f5SStefano Zampini   if (n_is) {
860d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
861d02579f5SStefano Zampini   }
86263602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
86363602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
86463602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
86563602bcaSStefano Zampini   }
86663602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
86763602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8681a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
86963602bcaSStefano Zampini   PetscFunctionReturn(0);
87063602bcaSStefano Zampini }
87163602bcaSStefano Zampini 
87263602bcaSStefano Zampini #undef __FUNCT__
87363602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
87463602bcaSStefano Zampini /*@
87563602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
87663602bcaSStefano Zampini 
87763602bcaSStefano Zampini    Collective
87863602bcaSStefano Zampini 
87963602bcaSStefano Zampini    Input Parameters:
88063602bcaSStefano Zampini +  pc - the preconditioning context
8810f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
8820f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in local ordering
88363602bcaSStefano Zampini 
88463602bcaSStefano Zampini    Level: intermediate
88563602bcaSStefano Zampini 
8860f202f7eSStefano Zampini    Notes:
8870f202f7eSStefano Zampini      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
88863602bcaSStefano Zampini 
88963602bcaSStefano Zampini .seealso: PCBDDC
89063602bcaSStefano Zampini @*/
89163602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
89263602bcaSStefano Zampini {
89363602bcaSStefano Zampini   PetscInt       i;
89463602bcaSStefano Zampini   PetscErrorCode ierr;
89563602bcaSStefano Zampini 
89663602bcaSStefano Zampini   PetscFunctionBegin;
89763602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
89863602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
89963602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
90063602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
90163602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
90263602bcaSStefano Zampini   }
903e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
90463602bcaSStefano Zampini   PetscFunctionReturn(0);
90563602bcaSStefano Zampini }
90663602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
90763602bcaSStefano Zampini 
90863602bcaSStefano Zampini #undef __FUNCT__
9099c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
9109c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
9119c0446d6SStefano Zampini {
9129c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
9139c0446d6SStefano Zampini   PetscInt i;
9149c0446d6SStefano Zampini   PetscErrorCode ierr;
9159c0446d6SStefano Zampini 
9169c0446d6SStefano Zampini   PetscFunctionBegin;
917da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
9189c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
9199c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
9209c0446d6SStefano Zampini   }
921d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
92263602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
92363602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
92463602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
92563602bcaSStefano Zampini   }
92663602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
92763602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
928da1bb401SStefano Zampini   /* allocate space then set */
929d02579f5SStefano Zampini   if (n_is) {
930785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
931d02579f5SStefano Zampini   }
9329c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
933da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
934da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
9359c0446d6SStefano Zampini   }
9369c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
93763602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
9381a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
9399c0446d6SStefano Zampini   PetscFunctionReturn(0);
9409c0446d6SStefano Zampini }
9411e6b0712SBarry Smith 
9429c0446d6SStefano Zampini #undef __FUNCT__
9439c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
9449c0446d6SStefano Zampini /*@
94563602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
9469c0446d6SStefano Zampini 
94763602bcaSStefano Zampini    Collective
9489c0446d6SStefano Zampini 
9499c0446d6SStefano Zampini    Input Parameters:
9509c0446d6SStefano Zampini +  pc - the preconditioning context
9510f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
9520f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in global ordering
9539c0446d6SStefano Zampini 
9549c0446d6SStefano Zampini    Level: intermediate
9559c0446d6SStefano Zampini 
9560f202f7eSStefano Zampini    Notes:
9570f202f7eSStefano Zampini      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
9589c0446d6SStefano Zampini 
9599c0446d6SStefano Zampini .seealso: PCBDDC
9609c0446d6SStefano Zampini @*/
9619c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9629c0446d6SStefano Zampini {
9632b510759SStefano Zampini   PetscInt       i;
9649c0446d6SStefano Zampini   PetscErrorCode ierr;
9659c0446d6SStefano Zampini 
9669c0446d6SStefano Zampini   PetscFunctionBegin;
9679c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
96863602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9692b510759SStefano Zampini   for (i=0;i<n_is;i++) {
97063602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
97163602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9722b510759SStefano Zampini   }
9739c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9749c0446d6SStefano Zampini   PetscFunctionReturn(0);
9759c0446d6SStefano Zampini }
976906d46d4SStefano Zampini 
977da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
978534831adSStefano Zampini #undef __FUNCT__
979534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
980534831adSStefano Zampini /* -------------------------------------------------------------------------- */
981534831adSStefano Zampini /*
982534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
983534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9849c0446d6SStefano Zampini 
985534831adSStefano Zampini    Input Parameter:
986534831adSStefano Zampini +  pc - the preconditioner contex
987534831adSStefano Zampini 
988534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
989534831adSStefano Zampini 
990534831adSStefano Zampini    Notes:
991534831adSStefano Zampini      The interface routine PCPreSolve() is not usually called directly by
992534831adSStefano Zampini    the user, but instead is called by KSPSolve().
993534831adSStefano Zampini */
994534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
995534831adSStefano Zampini {
996534831adSStefano Zampini   PetscErrorCode ierr;
997534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
998534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9993972b0daSStefano Zampini   Vec            used_vec;
10008d00608fSStefano Zampini   PetscBool      copy_rhs = PETSC_TRUE;
100151694757SStefano Zampini   PetscBool      benign_correction_is_zero = PETSC_FALSE;
100251694757SStefano Zampini   PetscBool      iscg;
1003534831adSStefano Zampini 
1004534831adSStefano Zampini   PetscFunctionBegin;
100585c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
100685c4d303SStefano Zampini   if (ksp) {
100785c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
100885c4d303SStefano Zampini     if (!iscg) {
100985c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
101085c4d303SStefano Zampini     }
101185c4d303SStefano Zampini   }
101285c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
101362a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
101462a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
101562a6ff1dSStefano Zampini   }
101662a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
101762a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
101862a6ff1dSStefano Zampini   }
10198d00608fSStefano Zampini 
10203972b0daSStefano Zampini   if (x) {
10213972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
10223972b0daSStefano Zampini     used_vec = x;
10238d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
10243972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
10253972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
10263972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
10273972b0daSStefano Zampini   }
10288efcfb23SStefano Zampini 
10298efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
10303972b0daSStefano Zampini   if (ksp) {
1031a0cb1b98SStefano Zampini     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
10328efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
10338efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
10343972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
10353972b0daSStefano Zampini     }
10363972b0daSStefano Zampini   }
10373308cffdSStefano Zampini 
10388d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
10393972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
1040a07ea27aSStefano Zampini   if (rhs) {
10413975b054SStefano Zampini     IS dirIS = NULL;
10423975b054SStefano Zampini 
1043a07ea27aSStefano Zampini     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
10443975b054SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
10453975b054SStefano Zampini     if (dirIS) {
1046906d46d4SStefano Zampini       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
1047785d1243SStefano Zampini       PetscInt          dirsize,i,*is_indices;
10482b095fd8SStefano Zampini       PetscScalar       *array_x;
10492b095fd8SStefano Zampini       const PetscScalar *array_diagonal;
1050785d1243SStefano Zampini 
10513972b0daSStefano Zampini       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
10523972b0daSStefano Zampini       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1053e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1054e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1055e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1056e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
105782ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
10583972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10592b095fd8SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10603972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10612fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
10623972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10632b095fd8SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10643972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1065e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1066e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10678d00608fSStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
10681b968477SStefano Zampini       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
10698efcfb23SStefano Zampini     }
1070a07ea27aSStefano Zampini   }
1071b76ba322SStefano Zampini 
10728efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
10738d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
10748d00608fSStefano Zampini     /* store the original rhs */
10758d00608fSStefano Zampini     if (copy_rhs) {
10768d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10778d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10788d00608fSStefano Zampini     }
10798d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10803972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10810369aaf7SStefano Zampini     ierr = MatMultAdd(pc->mat,used_vec,rhs,rhs);CHKERRQ(ierr);
10823972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10838efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10847acc28cbSStefano Zampini     if (ksp) {
10857acc28cbSStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
10867acc28cbSStefano Zampini     }
10873308cffdSStefano Zampini   }
10888efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1089b76ba322SStefano Zampini 
10900369aaf7SStefano Zampini   /* When using the benign trick: (TODO: what about FETI-DP?)
10910369aaf7SStefano Zampini      - change rhs on pressures (iteration matrix is surely of type MATIS)
10920369aaf7SStefano Zampini      - compute initial vector in benign space
10930369aaf7SStefano Zampini   */
10940369aaf7SStefano Zampini   if (rhs && pcbddc->benign_saddle_point) {
109506f24817SStefano Zampini     Mat_IS   *matis = (Mat_IS*)(pc->mat->data);
10964f1b2e48SStefano Zampini     PetscInt i;
109706f24817SStefano Zampini 
109806f24817SStefano Zampini     /* store the original rhs */
109906f24817SStefano Zampini     if (copy_rhs) {
110006f24817SStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
110106f24817SStefano Zampini       copy_rhs = PETSC_FALSE;
110206f24817SStefano Zampini     }
110306f24817SStefano Zampini 
110406f24817SStefano Zampini     /* now change (locally) the basis */
110506f24817SStefano Zampini     ierr = VecScatterBegin(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
110606f24817SStefano Zampini     ierr = VecScatterEnd(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
110706f24817SStefano Zampini     if (pcbddc->benign_change) {
11081cf9b237SStefano Zampini       Mat M;
11091cf9b237SStefano Zampini 
111006f24817SStefano Zampini       ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
111106f24817SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
111206f24817SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11131cf9b237SStefano Zampini       /* change local iteration matrix TODO (case when mat == pmat) */
11141cf9b237SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
11151cf9b237SStefano Zampini       ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_INITIAL_MATRIX,2.0,&M);CHKERRQ(ierr);
11161cf9b237SStefano Zampini       ierr = MatSeqAIJCompress(M,&pcbddc->local_mat);CHKERRQ(ierr);
11171cf9b237SStefano Zampini       ierr = MatDestroy(&M);CHKERRQ(ierr);
111806f24817SStefano Zampini       ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
111906f24817SStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
112006f24817SStefano Zampini       pcbddc->benign_original_mat = matis->A;
112106f24817SStefano Zampini       ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
112206f24817SStefano Zampini       ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr);
112306f24817SStefano Zampini       matis->A = pcbddc->local_mat;
112406f24817SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
112506f24817SStefano Zampini     } else {
112606f24817SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
112706f24817SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
112806f24817SStefano Zampini     }
112906f24817SStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
113006f24817SStefano Zampini 
11310369aaf7SStefano Zampini     /* compute u^*_h as in Xuemin Tu's thesis (see Section 4.8.1) */
11320369aaf7SStefano Zampini     /* TODO: what about Stokes? */
11330369aaf7SStefano Zampini     if (!pcbddc->benign_vec) {
11340369aaf7SStefano Zampini       ierr = VecDuplicate(rhs,&pcbddc->benign_vec);CHKERRQ(ierr);
11350369aaf7SStefano Zampini     }
11364f1b2e48SStefano Zampini     ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
11374f1b2e48SStefano Zampini     if (pcbddc->benign_n) {
11380369aaf7SStefano Zampini       const PetscScalar *array;
11390369aaf7SStefano Zampini 
11400369aaf7SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array);CHKERRQ(ierr);
11414f1b2e48SStefano Zampini       for (i=0;i<pcbddc->benign_n;i++) pcbddc->benign_p0[i] = array[pcbddc->benign_p0_lidx[i]];
11420369aaf7SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array);CHKERRQ(ierr);
11430369aaf7SStefano Zampini     }
114451694757SStefano Zampini     if (pcbddc->benign_null && iscg) { /* this is a workaround, need to understand more */
11454f1b2e48SStefano Zampini       PetscBool iszero_l = PETSC_TRUE;
11464f1b2e48SStefano Zampini       for (i=0;i<pcbddc->benign_n;i++) {
11474f1b2e48SStefano Zampini         iszero_l = (iszero_l && (PetscAbsScalar(pcbddc->benign_p0[i]) < PETSC_SMALL ? PETSC_TRUE : PETSC_FALSE));
11484f1b2e48SStefano Zampini       }
114951694757SStefano Zampini       ierr = MPI_Allreduce(&iszero_l,&benign_correction_is_zero,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
115051694757SStefano Zampini     }
115151694757SStefano Zampini     if (!benign_correction_is_zero) {
115251694757SStefano Zampini       ierr = VecSet(pcis->vec1_global,0.);CHKERRQ(ierr);
1153015636ebSStefano Zampini       ierr = PCBDDCBenignGetOrSetP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr);
11540369aaf7SStefano Zampini       ierr = PCApply_BDDC(pc,pcis->vec1_global,pcbddc->benign_vec);CHKERRQ(ierr);
11554f1b2e48SStefano Zampini       ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
1156015636ebSStefano Zampini       ierr = PCBDDCBenignGetOrSetP0(pc,pcbddc->benign_vec,PETSC_FALSE);CHKERRQ(ierr);
1157b76ba322SStefano Zampini     }
115851694757SStefano Zampini   }
1159b76ba322SStefano Zampini 
11600369aaf7SStefano Zampini   /* change rhs and iteration matrix if using the change of basis */
1161b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1162906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1163906d46d4SStefano Zampini 
1164906d46d4SStefano Zampini     /* get change ctx */
1165906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1166906d46d4SStefano Zampini 
1167906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1168906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1169906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1170906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1171906d46d4SStefano Zampini 
1172906d46d4SStefano Zampini     /* change iteration matrix */
1173906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1174906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1175906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1176906d46d4SStefano Zampini 
11778d00608fSStefano Zampini     /* store the original rhs */
11788d00608fSStefano Zampini     if (copy_rhs) {
11798d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
11808d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
11818d00608fSStefano Zampini     }
11828d00608fSStefano Zampini 
1183906d46d4SStefano Zampini     /* change rhs */
1184906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1185906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
11868d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
118792e3dcfbSStefano Zampini   }
118892e3dcfbSStefano Zampini 
11890369aaf7SStefano Zampini   /* remove non-benign solution from the rhs */
11900369aaf7SStefano Zampini   if (pcbddc->benign_saddle_point) {
11910369aaf7SStefano Zampini     /* store the original rhs */
11920369aaf7SStefano Zampini     if (copy_rhs) {
11930369aaf7SStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
11940369aaf7SStefano Zampini       copy_rhs = PETSC_FALSE;
11950369aaf7SStefano Zampini     }
119651694757SStefano Zampini     if (benign_correction_is_zero) { /* still need to understand why it works great */
119751694757SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
119851694757SStefano Zampini       ierr = PCApply_BDDC(pc,rhs,pcbddc->benign_vec);CHKERRQ(ierr);
119951694757SStefano Zampini     }
12000369aaf7SStefano Zampini     ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr);
12010369aaf7SStefano Zampini     ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);CHKERRQ(ierr);
12020369aaf7SStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
12030369aaf7SStefano Zampini   }
12040369aaf7SStefano Zampini 
12050369aaf7SStefano Zampini   /* set initial guess if using PCG */
12060369aaf7SStefano Zampini   if (x && pcbddc->use_exact_dirichlet_trick) {
12070369aaf7SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
12080369aaf7SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12090369aaf7SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12100369aaf7SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
12110369aaf7SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12120369aaf7SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12130369aaf7SStefano Zampini     if (ksp) {
12140369aaf7SStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
12150369aaf7SStefano Zampini     }
12160369aaf7SStefano Zampini   }
1217b0f5fe93SStefano Zampini   if (pcbddc->benign_null) {
1218b0f5fe93SStefano Zampini     MatNullSpace null_space;
1219b0f5fe93SStefano Zampini     Vec          nullv;
1220b0f5fe93SStefano Zampini     PetscBool    isnull;
12214f1b2e48SStefano Zampini     PetscInt     i;
1222b0f5fe93SStefano Zampini 
12234f1b2e48SStefano Zampini     for (i=0;i<pcbddc->benign_n;i++) pcbddc->benign_p0[i] = 1.;
1224b0f5fe93SStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&nullv);CHKERRQ(ierr);
1225b0f5fe93SStefano Zampini     ierr = VecSet(nullv,0.);CHKERRQ(ierr);
1226b0f5fe93SStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,nullv,PETSC_FALSE);CHKERRQ(ierr);
1227b0f5fe93SStefano Zampini     ierr = VecNormalize(nullv,NULL);CHKERRQ(ierr);
1228b0f5fe93SStefano Zampini     ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)pc),PETSC_FALSE,1,&nullv,&null_space);CHKERRQ(ierr);
1229b0f5fe93SStefano Zampini     ierr = MatNullSpaceTest(null_space,pc->mat,&isnull);CHKERRQ(ierr);
12304f1b2e48SStefano Zampini     if (isnull) {
1231b0f5fe93SStefano Zampini       ierr = MatSetNullSpace(pc->mat,null_space);CHKERRQ(ierr);
12324f1b2e48SStefano Zampini     }
1233b0f5fe93SStefano Zampini     ierr = MatNullSpaceDestroy(&null_space);CHKERRQ(ierr);
1234b0f5fe93SStefano Zampini     ierr = VecDestroy(&nullv);CHKERRQ(ierr);
1235b0f5fe93SStefano Zampini   }
1236b0f5fe93SStefano Zampini 
12370369aaf7SStefano Zampini 
123892e3dcfbSStefano Zampini   /* remove nullspace if present */
12398efcfb23SStefano Zampini   if (ksp && x && pcbddc->NullSpace) {
12408efcfb23SStefano Zampini     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
12418d00608fSStefano Zampini     /* store the original rhs */
12428d00608fSStefano Zampini     if (copy_rhs) {
12438d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
12448d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
12458d00608fSStefano Zampini     }
12468d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
1247d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1248b76ba322SStefano Zampini   }
1249534831adSStefano Zampini   PetscFunctionReturn(0);
1250534831adSStefano Zampini }
1251906d46d4SStefano Zampini 
1252534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1253534831adSStefano Zampini #undef __FUNCT__
1254534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1255534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1256534831adSStefano Zampini /*
1257534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1258534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1259534831adSStefano Zampini 
1260534831adSStefano Zampini    Input Parameter:
1261534831adSStefano Zampini +  pc - the preconditioner contex
1262534831adSStefano Zampini 
1263534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1264534831adSStefano Zampini 
1265534831adSStefano Zampini    Notes:
1266534831adSStefano Zampini      The interface routine PCPostSolve() is not usually called directly by
1267534831adSStefano Zampini      the user, but instead is called by KSPSolve().
1268534831adSStefano Zampini */
1269534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1270534831adSStefano Zampini {
1271534831adSStefano Zampini   PetscErrorCode ierr;
1272534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1273534831adSStefano Zampini 
1274534831adSStefano Zampini   PetscFunctionBegin;
1275b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1276906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1277906d46d4SStefano Zampini 
1278906d46d4SStefano Zampini     /* get change ctx */
1279906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1280906d46d4SStefano Zampini 
1281906d46d4SStefano Zampini     /* restore iteration matrix */
1282906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1283906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1284906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
12859c6a02ceSStefano Zampini   }
1286906d46d4SStefano Zampini 
12878666afb4SStefano Zampini   /* need to restore the local matrices */
12889c6a02ceSStefano Zampini   if (pcbddc->benign_change) {
12898666afb4SStefano Zampini     Mat_IS *matis = (Mat_IS*)pc->mat->data;
12908666afb4SStefano Zampini 
12918666afb4SStefano Zampini     ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
12928666afb4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->benign_original_mat);CHKERRQ(ierr);
12938666afb4SStefano Zampini     matis->A = pcbddc->benign_original_mat;
12948666afb4SStefano Zampini     ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
12958666afb4SStefano Zampini   }
12968666afb4SStefano Zampini 
1297906d46d4SStefano Zampini   /* get solution in original basis */
1298906d46d4SStefano Zampini   if (x) {
1299906d46d4SStefano Zampini     PC_IS *pcis = (PC_IS*)(pc->data);
130006f24817SStefano Zampini 
13012d3346b4SStefano Zampini     /* restore solution on pressures */
130206a4e24aSStefano Zampini     if (pcbddc->benign_saddle_point) {
13032d3346b4SStefano Zampini       Mat_IS *matis = (Mat_IS*)pc->mat->data;
13042d3346b4SStefano Zampini 
13050369aaf7SStefano Zampini       /* add non-benign solution */
13068666afb4SStefano Zampini       ierr = VecAXPY(x,-1.0,pcbddc->benign_vec);CHKERRQ(ierr);
13070369aaf7SStefano Zampini 
13080369aaf7SStefano Zampini       /* change basis on pressures for x */
130906f24817SStefano Zampini       ierr = VecScatterBegin(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
131006f24817SStefano Zampini       ierr = VecScatterEnd(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13112d3346b4SStefano Zampini       if (pcbddc->benign_change) {
13122d3346b4SStefano Zampini         ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
131306f24817SStefano Zampini         ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
131406f24817SStefano Zampini         ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13152d3346b4SStefano Zampini       } else {
131606f24817SStefano Zampini         ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
131706f24817SStefano Zampini         ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13182d3346b4SStefano Zampini       }
13192d3346b4SStefano Zampini     }
13209c6a02ceSStefano Zampini 
13210369aaf7SStefano Zampini     /* change basis on x */
13229c6a02ceSStefano Zampini     if (pcbddc->ChangeOfBasisMatrix) {
13239c6a02ceSStefano Zampini       PCBDDCChange_ctx change_ctx;
13249c6a02ceSStefano Zampini 
13259c6a02ceSStefano Zampini       ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
13260369aaf7SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1327906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
13283425bc38SStefano Zampini     }
1329534831adSStefano Zampini   }
1330906d46d4SStefano Zampini 
13313972b0daSStefano Zampini   /* add solution removed in presolve */
13326bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
13333425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
13343425bc38SStefano Zampini   }
1335906d46d4SStefano Zampini 
1336fb223d50SStefano Zampini   /* restore rhs to its original state */
13378d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
1338fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1339fb223d50SStefano Zampini   }
13408d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
13418efcfb23SStefano Zampini 
13428efcfb23SStefano Zampini   /* restore ksp guess state */
13438efcfb23SStefano Zampini   if (ksp) {
13448efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
13458efcfb23SStefano Zampini   }
1346534831adSStefano Zampini   PetscFunctionReturn(0);
1347534831adSStefano Zampini }
1348534831adSStefano Zampini /* -------------------------------------------------------------------------- */
134953cdbc3dSStefano Zampini #undef __FUNCT__
135053cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
13510c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
13520c7d97c5SJed Brown /*
13530c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
13540c7d97c5SJed Brown                   by setting data structures and options.
13550c7d97c5SJed Brown 
13560c7d97c5SJed Brown    Input Parameter:
135753cdbc3dSStefano Zampini +  pc - the preconditioner context
13580c7d97c5SJed Brown 
13590c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
13600c7d97c5SJed Brown 
13610c7d97c5SJed Brown    Notes:
13620c7d97c5SJed Brown      The interface routine PCSetUp() is not usually called directly by
13630c7d97c5SJed Brown      the user, but instead is called by PCApply() if necessary.
13640c7d97c5SJed Brown */
136553cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
13660c7d97c5SJed Brown {
13670c7d97c5SJed Brown   PetscErrorCode ierr;
13680c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
13695e8657edSStefano Zampini   Mat_IS*        matis;
137008122e43SStefano Zampini   MatNullSpace   nearnullspace;
1371c9ed8603SStefano Zampini   IS             zerodiag = NULL;
137291e8d312SStefano Zampini   PetscInt       nrows,ncols;
137308122e43SStefano Zampini   PetscBool      computetopography,computesolvers,computesubschurs;
13748de1fae6SStefano Zampini   PetscBool      computeconstraintsmatrix;
13755e8657edSStefano Zampini   PetscBool      new_nearnullspace_provided,ismatis;
13760c7d97c5SJed Brown 
13770c7d97c5SJed Brown   PetscFunctionBegin;
13785e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
13795e8657edSStefano Zampini   if (!ismatis) {
13805e8657edSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
13815e8657edSStefano Zampini   }
138291e8d312SStefano Zampini   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
138391e8d312SStefano Zampini   if (nrows != ncols) {
138491e8d312SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
138591e8d312SStefano Zampini   }
13865e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1387f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
13883b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1389fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1390f4ddd8eeSStefano Zampini   /* split work */
1391674ae819SStefano Zampini   if (pc->setupcalled) {
1392d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1393674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1394674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1395674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1396674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1397674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1398674ae819SStefano Zampini     }
1399674ae819SStefano Zampini   } else {
1400674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1401674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1402674ae819SStefano Zampini   }
1403fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1404fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1405fb180af4SStefano Zampini   }
14068de1fae6SStefano Zampini   computeconstraintsmatrix = PETSC_FALSE;
1407b087196eSStefano Zampini 
1408b087196eSStefano Zampini   /* check parameters' compatibility */
14095a95e1ceSStefano Zampini   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
14105a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
14115a95e1ceSStefano Zampini   }
14125a95e1ceSStefano Zampini   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1413862806e4SStefano Zampini   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1414862806e4SStefano Zampini 
14155a95e1ceSStefano Zampini   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
14166816873aSStefano Zampini   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) {
14176816873aSStefano Zampini     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");
14186816873aSStefano Zampini   }
14192d3346b4SStefano Zampini 
14202d3346b4SStefano Zampini   /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */
14212d3346b4SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr);
142206a4e24aSStefano Zampini   if (pcbddc->benign_saddle_point && !ismatis) {
14232d3346b4SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS");
14242d3346b4SStefano Zampini   }
14252d3346b4SStefano Zampini 
1426f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
142770cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
142870cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
142958a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1430f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1431f4ddd8eeSStefano Zampini     }
143258a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1433f4ddd8eeSStefano Zampini   }
1434f4ddd8eeSStefano Zampini 
14355e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
14365e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
14375e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
14385e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
14395e8657edSStefano Zampini   } else {
1440b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
14415e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
14425e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
1443d16cbb6bSStefano Zampini   }
1444d16cbb6bSStefano Zampini 
14454f1b2e48SStefano Zampini   /* detect local disconnected subdomains if requested and not done before */
14464f1b2e48SStefano Zampini   if (pcbddc->detect_disconnected && !pcbddc->n_local_subs) {
14474f1b2e48SStefano Zampini     ierr = MatDetectDisconnectedComponents(pcbddc->local_mat,PETSC_FALSE,&pcbddc->n_local_subs,&pcbddc->local_subs);CHKERRQ(ierr);
14484f1b2e48SStefano Zampini   }
14494f1b2e48SStefano Zampini 
14504f1b2e48SStefano Zampini   /*
14514f1b2e48SStefano Zampini      change basis on local pressures (aka zerodiag dofs)
14524f1b2e48SStefano Zampini      This should come earlier then PCISSetUp for extracting the correct subdomain matrices
14534f1b2e48SStefano Zampini   */
1454d16cbb6bSStefano Zampini   if (pcbddc->benign_saddle_point) {
14559f47a83aSStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
14569f47a83aSStefano Zampini 
1457339f8db1SStefano Zampini     /* detect local saddle point and change the basis in pcbddc->local_mat */
1458339f8db1SStefano Zampini     ierr = PCBDDCBenignDetectSaddlePoint(pc,&zerodiag);CHKERRQ(ierr);
1459c263805aSStefano Zampini     /* pop B0 mat from pcbddc->local_mat */
1460c263805aSStefano Zampini     ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
14619f47a83aSStefano Zampini 
14629f47a83aSStefano Zampini     /* set flag in pcis to not reuse submatrices during PCISCreate */
14639f47a83aSStefano Zampini     pcis->reusesubmatrices = PETSC_FALSE;
1464674ae819SStefano Zampini   }
1465f4ddd8eeSStefano Zampini 
1466f1a72664SStefano Zampini   /* propagate relevant information */
14674f1b2e48SStefano Zampini #if !defined(PETSC_USE_COMPLEX) /* workaround for reals */
14683301b35fSStefano Zampini   if (matis->A->symmetric_set) {
14693301b35fSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
14703301b35fSStefano Zampini   }
1471e496cd5dSStefano Zampini #endif
147206a4e24aSStefano Zampini   if (matis->A->symmetric_set) {
147306a4e24aSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr);
147406a4e24aSStefano Zampini   }
147506a4e24aSStefano Zampini   if (matis->A->spd_set) {
147606a4e24aSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr);
147706a4e24aSStefano Zampini   }
1478e496cd5dSStefano Zampini 
14795e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
14805e8657edSStefano Zampini   {
14815e8657edSStefano Zampini     Mat temp_mat;
14825e8657edSStefano Zampini 
14835e8657edSStefano Zampini     temp_mat = matis->A;
14845e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
14855e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
14865e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
14875e8657edSStefano Zampini     matis->A = temp_mat;
14885e8657edSStefano Zampini   }
1489684f6988SStefano Zampini 
149081d14e9dSStefano Zampini   /* Analyze interface */
1491674ae819SStefano Zampini   if (computetopography) {
1492674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
14938de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1494674ae819SStefano Zampini   }
1495fb8d54d4SStefano Zampini 
14965408967cSStefano Zampini   /* check existence of a divergence free extension, i.e.
14975408967cSStefano Zampini      b(v_I,p_0) = 0 for all v_I (raise error if not).
14985408967cSStefano Zampini      Also, check that PCBDDCBenignGetOrSetP0 works */
14995408967cSStefano Zampini #if defined(PETSC_USE_DEBUG)
150009f581a4SStefano Zampini   if (pcbddc->benign_saddle_point) {
15015408967cSStefano Zampini     ierr = PCBDDCBenignCheck(pc,zerodiag);CHKERRQ(ierr);
150209f581a4SStefano Zampini   }
15035408967cSStefano Zampini #endif
15044f1b2e48SStefano Zampini   ierr = ISDestroy(&zerodiag);CHKERRQ(ierr);
150506f24817SStefano Zampini 
1506b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1507684f6988SStefano Zampini   if (computesolvers) {
1508d5574798SStefano Zampini     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1509d5574798SStefano Zampini 
1510d5574798SStefano Zampini     if (computesubschurs && computetopography) {
151108122e43SStefano Zampini       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1512b1b3d7a2SStefano Zampini     }
1513d5574798SStefano Zampini     if (sub_schurs->use_mumps) {
15142070dbb6SStefano Zampini       if (computesubschurs) {
151508122e43SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
15162070dbb6SStefano Zampini       }
1517d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1518d5574798SStefano Zampini     } else {
1519d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
15202070dbb6SStefano Zampini       if (computesubschurs) {
1521d5574798SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1522d5574798SStefano Zampini       }
15232070dbb6SStefano Zampini     }
152408122e43SStefano Zampini     if (pcbddc->adaptive_selection) {
152508122e43SStefano Zampini       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
15268de1fae6SStefano Zampini       computeconstraintsmatrix = PETSC_TRUE;
1527b7eb3628SStefano Zampini     }
1528b1b3d7a2SStefano Zampini   }
1529684f6988SStefano Zampini 
1530f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1531fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1532f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1533f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1534f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1535f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1536f4ddd8eeSStefano Zampini     } else {
1537f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1538f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1539f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1540165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1541f4ddd8eeSStefano Zampini         PetscInt         i;
1542165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1543165b64e2SStefano Zampini         PetscObjectState state;
1544165b64e2SStefano Zampini         PetscInt         nnsp_size;
1545165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1546f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1547f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1548165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1549f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1550f4ddd8eeSStefano Zampini             break;
1551f4ddd8eeSStefano Zampini           }
1552f4ddd8eeSStefano Zampini         }
1553f4ddd8eeSStefano Zampini       }
1554f4ddd8eeSStefano Zampini     }
1555f4ddd8eeSStefano Zampini   } else {
1556f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1557f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1558f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1559f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1560f4ddd8eeSStefano Zampini     }
1561f4ddd8eeSStefano Zampini   }
1562f4ddd8eeSStefano Zampini 
1563f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1564727cdba6SStefano Zampini   /* reset primal space flags */
1565f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1566727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
15678de1fae6SStefano Zampini   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1568727cdba6SStefano Zampini     /* It also sets the primal space flags */
1569674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1570e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1571f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
15723975b054SStefano Zampini   }
15735e8657edSStefano Zampini 
15743975b054SStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
15755e8657edSStefano Zampini     if (pcbddc->use_change_of_basis) {
15765e8657edSStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1577c263805aSStefano Zampini       Mat   temp_mat = NULL;
15785e8657edSStefano Zampini 
15794f1b2e48SStefano Zampini       if (pcbddc->benign_change) {
1580c263805aSStefano Zampini         /* insert B0 in pcbddc->local_mat */
1581c263805aSStefano Zampini         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr);
1582f1a72664SStefano Zampini         /* swap local matrices */
1583f1a72664SStefano Zampini         ierr = MatISGetLocalMat(pc->pmat,&temp_mat);CHKERRQ(ierr);
1584f1a72664SStefano Zampini         ierr = PetscObjectReference((PetscObject)temp_mat);CHKERRQ(ierr);
1585f1a72664SStefano Zampini         ierr = MatISSetLocalMat(pc->pmat,pcbddc->local_mat);CHKERRQ(ierr);
1586f1a72664SStefano Zampini         ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
1587c263805aSStefano Zampini       }
15885e8657edSStefano Zampini       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
15894f1b2e48SStefano Zampini       if (pcbddc->benign_change) {
1590c263805aSStefano Zampini         /* restore original matrix */
1591f1a72664SStefano Zampini         ierr = MatISSetLocalMat(pc->pmat,temp_mat);CHKERRQ(ierr);
1592f1a72664SStefano Zampini         ierr = PetscObjectDereference((PetscObject)temp_mat);CHKERRQ(ierr);
1593c263805aSStefano Zampini         /* pop B0 from pcbddc->local_mat */
1594c263805aSStefano Zampini         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1595c263805aSStefano Zampini       }
15965e8657edSStefano Zampini       /* get submatrices */
15975e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
15985e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
15995e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
16005e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
16015e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
16025e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
16033975b054SStefano Zampini       /* set flag in pcis to not reuse submatrices during PCISCreate */
16043975b054SStefano Zampini       pcis->reusesubmatrices = PETSC_FALSE;
16059c6a02ceSStefano Zampini     } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) {
1606b96c3477SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
16075e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
16085e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
16095e8657edSStefano Zampini     }
1610b96c3477SStefano Zampini     /* SetUp coarse and local Neumann solvers */
161199cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1612b96c3477SStefano Zampini     /* SetUp Scaling operator */
1613674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
16140c7d97c5SJed Brown   }
16150369aaf7SStefano Zampini 
161658a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
161758a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
16182b510759SStefano Zampini   }
16190c7d97c5SJed Brown   PetscFunctionReturn(0);
16200c7d97c5SJed Brown }
16210c7d97c5SJed Brown 
16220c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
16230c7d97c5SJed Brown /*
162450efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
16250c7d97c5SJed Brown 
16260c7d97c5SJed Brown    Input Parameters:
16270f202f7eSStefano Zampini +  pc - the preconditioner context
16280f202f7eSStefano Zampini -  r - input vector (global)
16290c7d97c5SJed Brown 
16300c7d97c5SJed Brown    Output Parameter:
16310c7d97c5SJed Brown .  z - output vector (global)
16320c7d97c5SJed Brown 
16330c7d97c5SJed Brown    Application Interface Routine: PCApply()
16340c7d97c5SJed Brown  */
16350c7d97c5SJed Brown #undef __FUNCT__
16360c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
163753cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
16380c7d97c5SJed Brown {
16390c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
16400c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1641b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
16420c7d97c5SJed Brown   PetscErrorCode    ierr;
16433b03a366Sstefano_zampini   const PetscScalar one = 1.0;
16443b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
16452617d88aSStefano Zampini   const PetscScalar zero = 0.0;
16460c7d97c5SJed Brown 
16470c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
16480c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
1649b097fa66SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
16500c7d97c5SJed Brown 
16510c7d97c5SJed Brown   PetscFunctionBegin;
1652015636ebSStefano Zampini   if (pcbddc->benign_saddle_point) { /* get p0 from r */
1653015636ebSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1654efc2fbd9SStefano Zampini   }
165585c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1656b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
16570c7d97c5SJed Brown     /* First Dirichlet solve */
16580c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16590c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16600c7d97c5SJed Brown     /*
16610c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1662b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1663674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
16640c7d97c5SJed Brown     */
1665b097fa66SStefano Zampini     if (n_D) {
1666b097fa66SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
16670c7d97c5SJed Brown       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
16688eeda7d8SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1669b097fa66SStefano Zampini       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1670b097fa66SStefano Zampini     } else {
1671b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1672b097fa66SStefano Zampini     }
16730c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16740c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1675674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1676b76ba322SStefano Zampini   } else {
1677b097fa66SStefano Zampini     if (pcbddc->switch_static) {
16780bdf917eSStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1679b097fa66SStefano Zampini     }
1680674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1681b76ba322SStefano Zampini   }
1682b76ba322SStefano Zampini 
16832617d88aSStefano Zampini   /* Apply interface preconditioner
16842617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1685dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
16862617d88aSStefano Zampini 
1687674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1688674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
16890c7d97c5SJed Brown 
16903b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
16910c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16920c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1693b097fa66SStefano Zampini   if (n_B) {
16940c7d97c5SJed Brown     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
16958eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1696b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1697b097fa66SStefano Zampini     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1698b097fa66SStefano Zampini   }
1699df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1700efc2fbd9SStefano Zampini 
1701b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1702b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1703b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1704b097fa66SStefano Zampini     } else {
1705b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1706b097fa66SStefano Zampini     }
17070c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17080c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1709b097fa66SStefano Zampini   } else {
1710b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1711b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1712b097fa66SStefano Zampini     } else {
1713b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1714b097fa66SStefano Zampini     }
1715b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1716b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1717b097fa66SStefano Zampini   }
1718efc2fbd9SStefano Zampini 
1719015636ebSStefano Zampini   if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */
1720015636ebSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
1721efc2fbd9SStefano Zampini   }
17220c7d97c5SJed Brown   PetscFunctionReturn(0);
17230c7d97c5SJed Brown }
172450efa1b5SStefano Zampini 
172550efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
172650efa1b5SStefano Zampini /*
172750efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
172850efa1b5SStefano Zampini 
172950efa1b5SStefano Zampini    Input Parameters:
17300f202f7eSStefano Zampini +  pc - the preconditioner context
17310f202f7eSStefano Zampini -  r - input vector (global)
173250efa1b5SStefano Zampini 
173350efa1b5SStefano Zampini    Output Parameter:
173450efa1b5SStefano Zampini .  z - output vector (global)
173550efa1b5SStefano Zampini 
173650efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
173750efa1b5SStefano Zampini  */
173850efa1b5SStefano Zampini #undef __FUNCT__
173950efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
174050efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
174150efa1b5SStefano Zampini {
174250efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
174350efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1744b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
174550efa1b5SStefano Zampini   PetscErrorCode    ierr;
174650efa1b5SStefano Zampini   const PetscScalar one = 1.0;
174750efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
174850efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
174950efa1b5SStefano Zampini 
175050efa1b5SStefano Zampini   PetscFunctionBegin;
1751537c1cdfSStefano Zampini   if (pcbddc->benign_saddle_point) { /* get p0 from r */
1752537c1cdfSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1753537c1cdfSStefano Zampini   }
175450efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1755b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
175650efa1b5SStefano Zampini     /* First Dirichlet solve */
175750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
175850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
175950efa1b5SStefano Zampini     /*
176050efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
1761b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
176250efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
176350efa1b5SStefano Zampini     */
1764b097fa66SStefano Zampini     if (n_D) {
1765b097fa66SStefano Zampini       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
176650efa1b5SStefano Zampini       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
176750efa1b5SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1768b097fa66SStefano Zampini       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1769b097fa66SStefano Zampini     } else {
1770b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1771b097fa66SStefano Zampini     }
177250efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
177350efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
177450efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
177550efa1b5SStefano Zampini   } else {
1776b097fa66SStefano Zampini     if (pcbddc->switch_static) {
177750efa1b5SStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1778b097fa66SStefano Zampini     }
177950efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
178050efa1b5SStefano Zampini   }
178150efa1b5SStefano Zampini 
178250efa1b5SStefano Zampini   /* Apply interface preconditioner
178350efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1784dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
178550efa1b5SStefano Zampini 
178650efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
178750efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
178850efa1b5SStefano Zampini 
178950efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
179050efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
179150efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1792b097fa66SStefano Zampini   if (n_B) {
179350efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
179450efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1795b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1796b097fa66SStefano Zampini     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1797b097fa66SStefano Zampini   }
1798b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1799b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1800b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1801b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1802b097fa66SStefano Zampini     } else {
1803b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1804b097fa66SStefano Zampini     }
180550efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
180650efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1807b097fa66SStefano Zampini   } else {
1808b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1809b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1810b097fa66SStefano Zampini     } else {
1811b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1812b097fa66SStefano Zampini     }
1813b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1814b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1815b097fa66SStefano Zampini   }
1816537c1cdfSStefano Zampini   if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */
1817537c1cdfSStefano Zampini     ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
1818537c1cdfSStefano Zampini   }
181950efa1b5SStefano Zampini   PetscFunctionReturn(0);
182050efa1b5SStefano Zampini }
1821da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1822674ae819SStefano Zampini 
1823da1bb401SStefano Zampini #undef __FUNCT__
1824da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1825da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1826da1bb401SStefano Zampini {
1827da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1828da1bb401SStefano Zampini   PetscErrorCode ierr;
1829da1bb401SStefano Zampini 
1830da1bb401SStefano Zampini   PetscFunctionBegin;
1831da1bb401SStefano Zampini   /* free data created by PCIS */
1832da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1833674ae819SStefano Zampini   /* free BDDC custom data  */
1834674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1835674ae819SStefano Zampini   /* destroy objects related to topography */
1836674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1837674ae819SStefano Zampini   /* free allocated graph structure */
1838da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1839b96c3477SStefano Zampini   /* free allocated sub schurs structure */
1840b96c3477SStefano Zampini   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
184134a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1842674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
184334a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1844674ae819SStefano Zampini   /* free solvers stuff */
1845674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
184662a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
184762a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
184862a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1849906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1850906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1851906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1852906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1853906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1854906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1855906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1856906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1857906d46d4SStefano Zampini   }
1858906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
18593425bc38SStefano Zampini   /* remove functions */
1860906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1861674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1862*30368db7SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",NULL);CHKERRQ(ierr);
1863bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
18642b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1865b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
18662b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1867bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1868bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
186982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1870bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
187182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1872bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
187382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1874bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1875785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1876bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
187763602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1878bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1879bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1880bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1881bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1882674ae819SStefano Zampini   /* Free the private data structure */
1883674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1884da1bb401SStefano Zampini   PetscFunctionReturn(0);
1885da1bb401SStefano Zampini }
18863425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
18871e6b0712SBarry Smith 
18883425bc38SStefano Zampini #undef __FUNCT__
18893425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
18903425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
18913425bc38SStefano Zampini {
1892674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
1893c08af4c6SStefano Zampini   Vec            copy_standard_rhs;
18943425bc38SStefano Zampini   PC_IS*         pcis;
18953425bc38SStefano Zampini   PC_BDDC*       pcbddc;
18963425bc38SStefano Zampini   PetscErrorCode ierr;
18970c7d97c5SJed Brown 
18983425bc38SStefano Zampini   PetscFunctionBegin;
18993425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
19003425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
19013425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
19023425bc38SStefano Zampini 
1903c08af4c6SStefano Zampini   /*
1904c08af4c6SStefano Zampini      change of basis for physical rhs if needed
1905c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
1906c08af4c6SStefano Zampini      TODO: better management when FETIDP will have its own class
1907c08af4c6SStefano Zampini   */
1908c08af4c6SStefano Zampini   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1909c08af4c6SStefano Zampini   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1910c08af4c6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
19113425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
1912c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1913c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1914fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1915fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
1916c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1917c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1918674ae819SStefano Zampini   /* Apply partition of unity */
19193425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1920c08af4c6SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
19218eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
19223425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
19233425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
19243425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
19253425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1926c08af4c6SStefano Zampini     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1927c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1928c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1929c08af4c6SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1930c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1931c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19323425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
19333425bc38SStefano Zampini   }
1934c08af4c6SStefano Zampini   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
19353425bc38SStefano Zampini   /* BDDC rhs */
19363425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
19378eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
19383425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
19393425bc38SStefano Zampini   }
19403425bc38SStefano Zampini   /* apply BDDC */
1941dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
19423425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
19433425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
19443425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
19453425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19463425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19473425bc38SStefano Zampini   PetscFunctionReturn(0);
19483425bc38SStefano Zampini }
19491e6b0712SBarry Smith 
19503425bc38SStefano Zampini #undef __FUNCT__
19513425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
19523425bc38SStefano Zampini /*@
19530f202f7eSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
19543425bc38SStefano Zampini 
19553425bc38SStefano Zampini    Collective
19563425bc38SStefano Zampini 
19573425bc38SStefano Zampini    Input Parameters:
19580f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
19590f202f7eSStefano Zampini -  standard_rhs    - the right-hand side of the original linear system
19603425bc38SStefano Zampini 
19613425bc38SStefano Zampini    Output Parameters:
19620f202f7eSStefano Zampini .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
19633425bc38SStefano Zampini 
19643425bc38SStefano Zampini    Level: developer
19653425bc38SStefano Zampini 
19663425bc38SStefano Zampini    Notes:
19673425bc38SStefano Zampini 
19680f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
19693425bc38SStefano Zampini @*/
19703425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
19713425bc38SStefano Zampini {
1972674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
19733425bc38SStefano Zampini   PetscErrorCode ierr;
19743425bc38SStefano Zampini 
19753425bc38SStefano Zampini   PetscFunctionBegin;
19763425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
19773425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
19783425bc38SStefano Zampini   PetscFunctionReturn(0);
19793425bc38SStefano Zampini }
19803425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
19811e6b0712SBarry Smith 
19823425bc38SStefano Zampini #undef __FUNCT__
19833425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
19843425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
19853425bc38SStefano Zampini {
1986674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
19873425bc38SStefano Zampini   PC_IS*         pcis;
19883425bc38SStefano Zampini   PC_BDDC*       pcbddc;
19893425bc38SStefano Zampini   PetscErrorCode ierr;
19903425bc38SStefano Zampini 
19913425bc38SStefano Zampini   PetscFunctionBegin;
19923425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
19933425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
19943425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
19953425bc38SStefano Zampini 
19963425bc38SStefano Zampini   /* apply B_delta^T */
19973425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19983425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19993425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
20003425bc38SStefano Zampini   /* compute rhs for BDDC application */
20013425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
20028eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
20033425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
20043425bc38SStefano Zampini   }
20053425bc38SStefano Zampini   /* apply BDDC */
2006dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
20073425bc38SStefano Zampini   /* put values into standard global vector */
20083425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20093425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20108eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
20113425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
20123425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
20133425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
20143425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
20153425bc38SStefano Zampini   }
20163425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20173425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20183425bc38SStefano Zampini   /* final change of basis if needed
20193425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
20203308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
20213425bc38SStefano Zampini   PetscFunctionReturn(0);
20223425bc38SStefano Zampini }
20231e6b0712SBarry Smith 
20243425bc38SStefano Zampini #undef __FUNCT__
20253425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
20263425bc38SStefano Zampini /*@
20270f202f7eSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
20283425bc38SStefano Zampini 
20293425bc38SStefano Zampini    Collective
20303425bc38SStefano Zampini 
20313425bc38SStefano Zampini    Input Parameters:
20320f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
20330f202f7eSStefano Zampini -  fetidp_flux_sol - the solution of the FETI-DP linear system
20343425bc38SStefano Zampini 
20353425bc38SStefano Zampini    Output Parameters:
20360f202f7eSStefano Zampini .  standard_sol    - the solution defined on the physical domain
20373425bc38SStefano Zampini 
20383425bc38SStefano Zampini    Level: developer
20393425bc38SStefano Zampini 
20403425bc38SStefano Zampini    Notes:
20413425bc38SStefano Zampini 
20420f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
20433425bc38SStefano Zampini @*/
20443425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
20453425bc38SStefano Zampini {
2046674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
20473425bc38SStefano Zampini   PetscErrorCode ierr;
20483425bc38SStefano Zampini 
20493425bc38SStefano Zampini   PetscFunctionBegin;
20503425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
20513425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
20523425bc38SStefano Zampini   PetscFunctionReturn(0);
20533425bc38SStefano Zampini }
20543425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
20551e6b0712SBarry Smith 
2056f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
2057edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
2058f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
2059f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
2060edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
2061f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
2062674ae819SStefano Zampini 
20633425bc38SStefano Zampini #undef __FUNCT__
20643425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
20653425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
20663425bc38SStefano Zampini {
2067674ae819SStefano Zampini 
2068674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
20693425bc38SStefano Zampini   Mat            newmat;
2070674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
20713425bc38SStefano Zampini   PC             newpc;
2072ce94432eSBarry Smith   MPI_Comm       comm;
20733425bc38SStefano Zampini   PetscErrorCode ierr;
20743425bc38SStefano Zampini 
20753425bc38SStefano Zampini   PetscFunctionBegin;
2076ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
20773425bc38SStefano Zampini   /* FETIDP linear matrix */
20783425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
20793425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
20803425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
20813425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
2082edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
20833425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
20843425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
20853425bc38SStefano Zampini   /* FETIDP preconditioner */
20863425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
20873425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
20883425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
20893425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
20903425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
20913425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
2092edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
20933425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
209423ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
20953425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
20963425bc38SStefano Zampini   /* return pointers for objects created */
20973425bc38SStefano Zampini   *fetidp_mat=newmat;
20983425bc38SStefano Zampini   *fetidp_pc=newpc;
20993425bc38SStefano Zampini   PetscFunctionReturn(0);
21003425bc38SStefano Zampini }
21011e6b0712SBarry Smith 
21023425bc38SStefano Zampini #undef __FUNCT__
21033425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
21043425bc38SStefano Zampini /*@
21050f202f7eSStefano Zampini  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
21063425bc38SStefano Zampini 
21073425bc38SStefano Zampini    Collective
21083425bc38SStefano Zampini 
21093425bc38SStefano Zampini    Input Parameters:
21100f202f7eSStefano Zampini .  pc - the BDDC preconditioning context (setup should have been called before)
211128509bceSStefano Zampini 
211228509bceSStefano Zampini    Output Parameters:
21130f202f7eSStefano Zampini +  fetidp_mat - shell FETI-DP matrix object
21140f202f7eSStefano Zampini -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
211528509bceSStefano Zampini 
211628509bceSStefano Zampini    Options Database Keys:
21170f202f7eSStefano Zampini .    -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
21183425bc38SStefano Zampini 
21193425bc38SStefano Zampini    Level: developer
21203425bc38SStefano Zampini 
21213425bc38SStefano Zampini    Notes:
21220f202f7eSStefano Zampini      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
21233425bc38SStefano Zampini 
21240f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
21253425bc38SStefano Zampini @*/
21263425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
21273425bc38SStefano Zampini {
21283425bc38SStefano Zampini   PetscErrorCode ierr;
21293425bc38SStefano Zampini 
21303425bc38SStefano Zampini   PetscFunctionBegin;
21313425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
21323425bc38SStefano Zampini   if (pc->setupcalled) {
2133516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
2134f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
21353425bc38SStefano Zampini   PetscFunctionReturn(0);
21363425bc38SStefano Zampini }
21370c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
2138da1bb401SStefano Zampini /*MC
2139da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
21400c7d97c5SJed Brown 
214128509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
214228509bceSStefano Zampini 
214328509bceSStefano Zampini .vb
214428509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
214528509bceSStefano 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
214628509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
21470f202f7eSStefano 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
214828509bceSStefano Zampini .ve
214928509bceSStefano Zampini 
215028509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
215128509bceSStefano Zampini 
21520f202f7eSStefano Zampini    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
215328509bceSStefano Zampini 
215428509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
215528509bceSStefano Zampini 
2156b6fdb6dfSStefano 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.
2157b6fdb6dfSStefano Zampini 
21580f202f7eSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
215928509bceSStefano Zampini 
21600f202f7eSStefano 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()
2161*30368db7SStefano Zampini    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesIS() and their local counterparts.
216228509bceSStefano Zampini 
21630f202f7eSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
216428509bceSStefano Zampini 
21650f202f7eSStefano 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.
21660f202f7eSStefano Zampini    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
216728509bceSStefano Zampini 
21680f202f7eSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
216928509bceSStefano Zampini 
21700f202f7eSStefano 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.
217128509bceSStefano Zampini 
21720f202f7eSStefano 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.
21730f202f7eSStefano Zampini    Deluxe scaling is not supported yet for FETI-DP.
21740f202f7eSStefano Zampini 
21750f202f7eSStefano Zampini    Options Database Keys (some of them, run with -h for a complete list):
21760f202f7eSStefano Zampini 
21770f202f7eSStefano Zampini .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
21780f202f7eSStefano Zampini .    -pc_bddc_use_edges <true> - use or not edges in primal space
21790f202f7eSStefano Zampini .    -pc_bddc_use_faces <false> - use or not faces in primal space
21800f202f7eSStefano Zampini .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
21810f202f7eSStefano Zampini .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
21820f202f7eSStefano Zampini .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
21830f202f7eSStefano Zampini .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
218428509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
21850f202f7eSStefano 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)
21860f202f7eSStefano 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)
21870f202f7eSStefano Zampini .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
21880f202f7eSStefano 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)
21890f202f7eSStefano 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)
219028509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
219128509bceSStefano Zampini 
219228509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
219328509bceSStefano Zampini .vb
219428509bceSStefano Zampini       -pc_bddc_dirichlet_
219528509bceSStefano Zampini       -pc_bddc_neumann_
219628509bceSStefano Zampini       -pc_bddc_coarse_
219728509bceSStefano Zampini .ve
21980f202f7eSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
219928509bceSStefano Zampini 
22000f202f7eSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
220128509bceSStefano Zampini .vb
2202312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
2203312be037SStefano Zampini       -pc_bddc_neumann_lN_
2204312be037SStefano Zampini       -pc_bddc_coarse_lN_
220528509bceSStefano Zampini .ve
22060f202f7eSStefano Zampini    Note that level number ranges from the finest (0) to the coarsest (N).
22070f202f7eSStefano 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.
22080f202f7eSStefano Zampini .vb
22090f202f7eSStefano Zampini      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
22100f202f7eSStefano Zampini .ve
22110f202f7eSStefano 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
2212da1bb401SStefano Zampini 
2213da1bb401SStefano Zampini    Level: intermediate
2214da1bb401SStefano Zampini 
2215b6fdb6dfSStefano Zampini    Developer notes:
2216da1bb401SStefano Zampini 
2217da1bb401SStefano Zampini    Contributed by Stefano Zampini
2218da1bb401SStefano Zampini 
2219da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2220da1bb401SStefano Zampini M*/
2221b2573a8aSBarry Smith 
2222da1bb401SStefano Zampini #undef __FUNCT__
2223da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
22248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2225da1bb401SStefano Zampini {
2226da1bb401SStefano Zampini   PetscErrorCode      ierr;
2227da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
2228da1bb401SStefano Zampini 
2229da1bb401SStefano Zampini   PetscFunctionBegin;
2230da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
2231b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
2232da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
2233da1bb401SStefano Zampini 
2234da1bb401SStefano Zampini   /* create PCIS data structure */
2235da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
2236da1bb401SStefano Zampini 
223747d04d0dSStefano Zampini   /* BDDC customization */
223808a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
223947d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
224047d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
224147d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
224247d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
224347d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
224447d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
2245fa434743SStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE;
2246fa434743SStefano Zampini   pcbddc->use_qr_single       = PETSC_FALSE;
22473301b35fSStefano Zampini   pcbddc->symmetric_primal    = PETSC_TRUE;
224806a4e24aSStefano Zampini   pcbddc->benign_saddle_point = PETSC_FALSE;
224947d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
2250b9d89cd5SStefano Zampini   /* private */
2251727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
22520e6343abSStefano Zampini   pcbddc->local_primal_size_cc       = 0;
22530e6343abSStefano Zampini   pcbddc->local_primal_ref_node      = 0;
22540e6343abSStefano Zampini   pcbddc->local_primal_ref_mult      = 0;
2255e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
2256727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
2257fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
225868457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
2259f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
2260727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
2261f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
2262f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
2263f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
2264674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
2265*30368db7SStefano Zampini   pcbddc->user_primal_vertices_local = 0;
22660bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
22673972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
2268534831adSStefano Zampini   pcbddc->original_rhs               = 0;
2269534831adSStefano Zampini   pcbddc->local_mat                  = 0;
2270534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
2271b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
2272906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
2273da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
2274da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
2275da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
2276da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
227715aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
227815aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
2279da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
2280da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
2281da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
2282da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
2283da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
2284da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
2285da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
2286da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
2287da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
2288da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
2289785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
2290785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
2291785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
229260d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
229360d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
229463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
2295da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
229663602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
2297da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
229885c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
229947d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
230047d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
2301b0c7d250SStefano Zampini   pcbddc->coarse_adj_red             = 0;
23024fad6a16SStefano Zampini   pcbddc->current_level              = 0;
23032b510759SStefano Zampini   pcbddc->max_levels                 = 0;
2304323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
230574e2c79eSStefano Zampini   pcbddc->redistribute_coarse        = 0;
2306f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
2307323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
23084f1b2e48SStefano Zampini   pcbddc->detect_disconnected        = PETSC_FALSE;
23094f1b2e48SStefano Zampini   pcbddc->n_local_subs               = 0;
23104f1b2e48SStefano Zampini   pcbddc->local_subs                 = NULL;
231181d14e9dSStefano Zampini 
231281d14e9dSStefano Zampini   /* benign subspace trick */
231381d14e9dSStefano Zampini   pcbddc->benign_change              = 0;
23140369aaf7SStefano Zampini   pcbddc->benign_vec                 = 0;
23150369aaf7SStefano Zampini   pcbddc->benign_original_mat        = 0;
23160369aaf7SStefano Zampini   pcbddc->benign_sf                  = 0;
23174f1b2e48SStefano Zampini   pcbddc->benign_B0                  = 0;
23184f1b2e48SStefano Zampini   pcbddc->benign_n                   = 0;
23194f1b2e48SStefano Zampini   pcbddc->benign_p0                  = NULL;
23204f1b2e48SStefano Zampini   pcbddc->benign_p0_lidx             = NULL;
23214f1b2e48SStefano Zampini   pcbddc->benign_p0_gidx             = NULL;
2322b0f5fe93SStefano Zampini   pcbddc->benign_null                = PETSC_FALSE;
232381d14e9dSStefano Zampini 
2324674ae819SStefano Zampini   /* create local graph structure */
2325674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2326674ae819SStefano Zampini 
2327674ae819SStefano Zampini   /* scaling */
2328674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
232934a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2330ac632422SStefano Zampini   pcbddc->faster_deluxe         = PETSC_FALSE;
2331b96c3477SStefano Zampini 
2332b96c3477SStefano Zampini   /* create sub schurs structure */
2333b96c3477SStefano Zampini   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2334b96c3477SStefano Zampini   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2335b96c3477SStefano Zampini   pcbddc->sub_schurs_layers      = -1;
2336b96c3477SStefano Zampini   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2337b96c3477SStefano Zampini 
2338b96c3477SStefano Zampini   pcbddc->computed_rowadj = PETSC_FALSE;
2339da1bb401SStefano Zampini 
2340b7eb3628SStefano Zampini   /* adaptivity */
2341f6f667cfSStefano Zampini   pcbddc->adaptive_threshold      = 0.0;
234208122e43SStefano Zampini   pcbddc->adaptive_nmax           = 0;
2343f6f667cfSStefano Zampini   pcbddc->adaptive_nmin           = 0;
2344b7eb3628SStefano Zampini 
2345da1bb401SStefano Zampini   /* function pointers */
2346da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
234793bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2348da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
2349da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
2350da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2351da1bb401SStefano Zampini   pc->ops->view                = 0;
2352da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
2353da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
2354da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
2355534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
2356534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
2357da1bb401SStefano Zampini 
2358da1bb401SStefano Zampini   /* composing function */
2359906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2360674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2361*30368db7SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",PCBDDCSetPrimalVerticesIS_BDDC);CHKERRQ(ierr);
2362bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
23632b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2364b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
23652b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2366bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2367bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
236882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2369bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
237082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2371bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
237282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2373bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
237482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2375bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
237663602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2377bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2378bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2379bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2380bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2381da1bb401SStefano Zampini   PetscFunctionReturn(0);
2382da1bb401SStefano Zampini }
23833425bc38SStefano Zampini 
2384