xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision efc2fbd988ffb32ec2d0f2ce2df54086b651f77c)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    Solvers
4a0d3c3abSStefano Zampini    - Add support for cholesky for coarse solver (similar to local solvers)
5eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
6eb97c9d2SStefano Zampini 
7eb97c9d2SStefano Zampini    User interface
80f202f7eSStefano Zampini    - ** DM attached to pc?
9eb97c9d2SStefano Zampini 
10eb97c9d2SStefano Zampini    Debugging output
11b9b85e73SStefano Zampini    - * Better management of verbosity levels of debugging output
12eb97c9d2SStefano Zampini 
13eb97c9d2SStefano Zampini    Extra
14b9b85e73SStefano Zampini    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15eb97c9d2SStefano Zampini    - BDDC with MG framework?
16eb97c9d2SStefano Zampini 
17eb97c9d2SStefano Zampini    FETIDP
18eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
19eb97c9d2SStefano Zampini 
20eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
21eb97c9d2SStefano Zampini    - Provide general case for subassembling
22eb97c9d2SStefano Zampini 
2353cdbc3dSStefano Zampini */
240c7d97c5SJed Brown 
25ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
26ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
273b03a366Sstefano_zampini #include <petscblaslapack.h>
28674ae819SStefano Zampini 
290c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
300c7d97c5SJed Brown #undef __FUNCT__
310c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
328c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc)
330c7d97c5SJed Brown {
340c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
350c7d97c5SJed Brown   PetscErrorCode ierr;
360c7d97c5SJed Brown 
370c7d97c5SJed Brown   PetscFunctionBegin;
38e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr);
398eeda7d8SStefano Zampini   /* Verbose debugging */
408eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
418eeda7d8SStefano Zampini   /* Primal space cumstomization */
4208a5cf49SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr);
438eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr);
448eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr);
458eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr);
466661aa1dSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);CHKERRQ(ierr);
47ac632422SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr);
488eeda7d8SStefano Zampini   /* Change of basis */
49b9b85e73SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr);
50b9b85e73SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr);
51674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
52674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
53674ae819SStefano Zampini   }
548eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
558eeda7d8SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);CHKERRQ(ierr);
5674e2c79eSStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);CHKERRQ(ierr);
570298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
582b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
59323d395dSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr);
60674ae819SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr);
61b96c3477SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);CHKERRQ(ierr);
62b96c3477SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);CHKERRQ(ierr);
63b96c3477SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);CHKERRQ(ierr);
64ac632422SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);CHKERRQ(ierr);
654c6709b3SStefano Zampini   ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr);
6608122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
6708122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
683301b35fSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr);
69b0c7d250SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);CHKERRQ(ierr);
7006a4e24aSStefano 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);
710c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
720c7d97c5SJed Brown   PetscFunctionReturn(0);
730c7d97c5SJed Brown }
740c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
75674ae819SStefano Zampini #undef __FUNCT__
76906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
77906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
78b9b85e73SStefano Zampini {
79b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
80b9b85e73SStefano Zampini   PetscErrorCode ierr;
81b9b85e73SStefano Zampini 
82b9b85e73SStefano Zampini   PetscFunctionBegin;
83b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
84b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
85b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
86b9b85e73SStefano Zampini   PetscFunctionReturn(0);
87b9b85e73SStefano Zampini }
88b9b85e73SStefano Zampini #undef __FUNCT__
89906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
90b9b85e73SStefano Zampini /*@
91906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
92b9b85e73SStefano Zampini 
93b9b85e73SStefano Zampini    Collective on PC
94b9b85e73SStefano Zampini 
95b9b85e73SStefano Zampini    Input Parameters:
96b9b85e73SStefano Zampini +  pc - the preconditioning context
97906d46d4SStefano Zampini -  change - the change of basis matrix
98b9b85e73SStefano Zampini 
99b9b85e73SStefano Zampini    Level: intermediate
100b9b85e73SStefano Zampini 
101b9b85e73SStefano Zampini    Notes:
102b9b85e73SStefano Zampini 
103b9b85e73SStefano Zampini .seealso: PCBDDC
104b9b85e73SStefano Zampini @*/
105906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
106b9b85e73SStefano Zampini {
107b9b85e73SStefano Zampini   PetscErrorCode ierr;
108b9b85e73SStefano Zampini 
109b9b85e73SStefano Zampini   PetscFunctionBegin;
110b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
111b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
112906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
113906d46d4SStefano Zampini   if (pc->mat) {
114906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
115906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
116906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
117906d46d4SStefano Zampini     if (rows_c != rows) {
118906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
119906d46d4SStefano Zampini     }
120906d46d4SStefano Zampini     if (cols_c != cols) {
121906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
122906d46d4SStefano Zampini     }
123906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
124906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
125906d46d4SStefano Zampini     if (rows_c != rows) {
126906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
127906d46d4SStefano Zampini     }
128906d46d4SStefano Zampini     if (cols_c != cols) {
129906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
130906d46d4SStefano Zampini     }
131906d46d4SStefano Zampini   }
132906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
133b9b85e73SStefano Zampini   PetscFunctionReturn(0);
134b9b85e73SStefano Zampini }
135b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
136b9b85e73SStefano Zampini #undef __FUNCT__
137674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
138674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
139674ae819SStefano Zampini {
140674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
141674ae819SStefano Zampini   PetscErrorCode ierr;
1421e6b0712SBarry Smith 
143674ae819SStefano Zampini   PetscFunctionBegin;
144674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
145674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
146674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
147674ae819SStefano Zampini   PetscFunctionReturn(0);
148674ae819SStefano Zampini }
149674ae819SStefano Zampini #undef __FUNCT__
150674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
151674ae819SStefano Zampini /*@
15228509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
153674ae819SStefano Zampini 
15417eb1463SStefano Zampini    Collective
155674ae819SStefano Zampini 
156674ae819SStefano Zampini    Input Parameters:
157674ae819SStefano Zampini +  pc - the preconditioning context
15817eb1463SStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering (can be empty)
159674ae819SStefano Zampini 
160674ae819SStefano Zampini    Level: intermediate
161674ae819SStefano Zampini 
162674ae819SStefano Zampini    Notes:
163674ae819SStefano Zampini 
164674ae819SStefano Zampini .seealso: PCBDDC
165674ae819SStefano Zampini @*/
166674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
167674ae819SStefano Zampini {
168674ae819SStefano Zampini   PetscErrorCode ierr;
169674ae819SStefano Zampini 
170674ae819SStefano Zampini   PetscFunctionBegin;
171674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
172674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
17317eb1463SStefano Zampini   PetscCheckSameComm(pc,1,PrimalVertices,2);
174674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
175674ae819SStefano Zampini   PetscFunctionReturn(0);
176674ae819SStefano Zampini }
177674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1780c7d97c5SJed Brown #undef __FUNCT__
1794fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1804fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1814fad6a16SStefano Zampini {
1824fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1834fad6a16SStefano Zampini 
1844fad6a16SStefano Zampini   PetscFunctionBegin;
1854fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1864fad6a16SStefano Zampini   PetscFunctionReturn(0);
1874fad6a16SStefano Zampini }
1881e6b0712SBarry Smith 
1894fad6a16SStefano Zampini #undef __FUNCT__
1904fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1914fad6a16SStefano Zampini /*@
19228509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1934fad6a16SStefano Zampini 
1944fad6a16SStefano Zampini    Logically collective on PC
1954fad6a16SStefano Zampini 
1964fad6a16SStefano Zampini    Input Parameters:
1974fad6a16SStefano Zampini +  pc - the preconditioning context
19828509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
1994fad6a16SStefano Zampini 
2000f202f7eSStefano Zampini    Options Database Keys:
2010f202f7eSStefano Zampini .    -pc_bddc_coarsening_ratio
2024fad6a16SStefano Zampini 
2034fad6a16SStefano Zampini    Level: intermediate
2044fad6a16SStefano Zampini 
2054fad6a16SStefano Zampini    Notes:
2060f202f7eSStefano Zampini      Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2074fad6a16SStefano Zampini 
2080f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetLevels()
2094fad6a16SStefano Zampini @*/
2104fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2114fad6a16SStefano Zampini {
2124fad6a16SStefano Zampini   PetscErrorCode ierr;
2134fad6a16SStefano Zampini 
2144fad6a16SStefano Zampini   PetscFunctionBegin;
2154fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2162b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2174fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2184fad6a16SStefano Zampini   PetscFunctionReturn(0);
2194fad6a16SStefano Zampini }
2202b510759SStefano Zampini 
221b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2222b510759SStefano Zampini #undef __FUNCT__
223b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
224b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
225b8ffe317SStefano Zampini {
226b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
227b8ffe317SStefano Zampini 
228b8ffe317SStefano Zampini   PetscFunctionBegin;
22985c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
230b8ffe317SStefano Zampini   PetscFunctionReturn(0);
231b8ffe317SStefano Zampini }
232b8ffe317SStefano Zampini 
233b8ffe317SStefano Zampini #undef __FUNCT__
234b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
235b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2362b510759SStefano Zampini {
2372b510759SStefano Zampini   PetscErrorCode ierr;
2382b510759SStefano Zampini 
2392b510759SStefano Zampini   PetscFunctionBegin;
2402b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
241b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
242b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2432b510759SStefano Zampini   PetscFunctionReturn(0);
2442b510759SStefano Zampini }
2451e6b0712SBarry Smith 
2464fad6a16SStefano Zampini #undef __FUNCT__
2472b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2482b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2494fad6a16SStefano Zampini {
2504fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2514fad6a16SStefano Zampini 
2524fad6a16SStefano Zampini   PetscFunctionBegin;
2532b510759SStefano Zampini   pcbddc->current_level = level;
2544fad6a16SStefano Zampini   PetscFunctionReturn(0);
2554fad6a16SStefano Zampini }
2561e6b0712SBarry Smith 
2574fad6a16SStefano Zampini #undef __FUNCT__
258b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
259b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
260b8ffe317SStefano Zampini {
261b8ffe317SStefano Zampini   PetscErrorCode ierr;
262b8ffe317SStefano Zampini 
263b8ffe317SStefano Zampini   PetscFunctionBegin;
264b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
265b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
266b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
267b8ffe317SStefano Zampini   PetscFunctionReturn(0);
268b8ffe317SStefano Zampini }
269b8ffe317SStefano Zampini 
270b8ffe317SStefano Zampini #undef __FUNCT__
2712b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2722b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2732b510759SStefano Zampini {
2742b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2752b510759SStefano Zampini 
2762b510759SStefano Zampini   PetscFunctionBegin;
2772b510759SStefano Zampini   pcbddc->max_levels = levels;
2782b510759SStefano Zampini   PetscFunctionReturn(0);
2792b510759SStefano Zampini }
2802b510759SStefano Zampini 
2812b510759SStefano Zampini #undef __FUNCT__
2822b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2834fad6a16SStefano Zampini /*@
28428509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2854fad6a16SStefano Zampini 
2864fad6a16SStefano Zampini    Logically collective on PC
2874fad6a16SStefano Zampini 
2884fad6a16SStefano Zampini    Input Parameters:
2894fad6a16SStefano Zampini +  pc - the preconditioning context
29028509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2914fad6a16SStefano Zampini 
2920f202f7eSStefano Zampini    Options Database Keys:
2930f202f7eSStefano Zampini .    -pc_bddc_levels
2944fad6a16SStefano Zampini 
2954fad6a16SStefano Zampini    Level: intermediate
2964fad6a16SStefano Zampini 
2974fad6a16SStefano Zampini    Notes:
2980f202f7eSStefano Zampini      Default value is 0, i.e. traditional one-level BDDC
2994fad6a16SStefano Zampini 
3000f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
3014fad6a16SStefano Zampini @*/
3022b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3034fad6a16SStefano Zampini {
3044fad6a16SStefano Zampini   PetscErrorCode ierr;
3054fad6a16SStefano Zampini 
3064fad6a16SStefano Zampini   PetscFunctionBegin;
3074fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3082b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
309312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3102b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3114fad6a16SStefano Zampini   PetscFunctionReturn(0);
3124fad6a16SStefano Zampini }
3134fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3141e6b0712SBarry Smith 
3154fad6a16SStefano Zampini #undef __FUNCT__
3160bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3170bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3180bdf917eSStefano Zampini {
3190bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3200bdf917eSStefano Zampini   PetscErrorCode ierr;
3210bdf917eSStefano Zampini 
3220bdf917eSStefano Zampini   PetscFunctionBegin;
3230bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3240bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3250bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3260bdf917eSStefano Zampini   PetscFunctionReturn(0);
3270bdf917eSStefano Zampini }
3281e6b0712SBarry Smith 
3290bdf917eSStefano Zampini #undef __FUNCT__
3300bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3310bdf917eSStefano Zampini /*@
33228509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3330bdf917eSStefano Zampini 
3340bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3350bdf917eSStefano Zampini 
3360bdf917eSStefano Zampini    Input Parameters:
3370bdf917eSStefano Zampini +  pc - the preconditioning context
33828509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3390bdf917eSStefano Zampini 
3400bdf917eSStefano Zampini    Level: intermediate
3410bdf917eSStefano Zampini 
3420bdf917eSStefano Zampini    Notes:
3430bdf917eSStefano Zampini 
3440bdf917eSStefano Zampini .seealso: PCBDDC
3450bdf917eSStefano Zampini @*/
3460bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3470bdf917eSStefano Zampini {
3480bdf917eSStefano Zampini   PetscErrorCode ierr;
3490bdf917eSStefano Zampini 
3500bdf917eSStefano Zampini   PetscFunctionBegin;
3510bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
352674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
3532b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
3540bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
3550bdf917eSStefano Zampini   PetscFunctionReturn(0);
3560bdf917eSStefano Zampini }
3570bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3581e6b0712SBarry Smith 
3590bdf917eSStefano Zampini #undef __FUNCT__
3603b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3613b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3623b03a366Sstefano_zampini {
3633b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3643b03a366Sstefano_zampini   PetscErrorCode ierr;
3653b03a366Sstefano_zampini 
3663b03a366Sstefano_zampini   PetscFunctionBegin;
367785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
368785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3693b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
37036e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
37136e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
372fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3733b03a366Sstefano_zampini   PetscFunctionReturn(0);
3743b03a366Sstefano_zampini }
3751e6b0712SBarry Smith 
3763b03a366Sstefano_zampini #undef __FUNCT__
3773b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3783b03a366Sstefano_zampini /*@
37928509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3803b03a366Sstefano_zampini 
381785d1243SStefano Zampini    Collective
3823b03a366Sstefano_zampini 
3833b03a366Sstefano_zampini    Input Parameters:
3843b03a366Sstefano_zampini +  pc - the preconditioning context
385785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3863b03a366Sstefano_zampini 
3873b03a366Sstefano_zampini    Level: intermediate
3883b03a366Sstefano_zampini 
3890f202f7eSStefano Zampini    Notes:
3900f202f7eSStefano Zampini      Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
3913b03a366Sstefano_zampini 
3920f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
3933b03a366Sstefano_zampini @*/
3943b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3953b03a366Sstefano_zampini {
3963b03a366Sstefano_zampini   PetscErrorCode ierr;
3973b03a366Sstefano_zampini 
3983b03a366Sstefano_zampini   PetscFunctionBegin;
3993b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
400674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
401785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4023b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4033b03a366Sstefano_zampini   PetscFunctionReturn(0);
4043b03a366Sstefano_zampini }
4053b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4061e6b0712SBarry Smith 
4073b03a366Sstefano_zampini #undef __FUNCT__
40882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
40982ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4100c7d97c5SJed Brown {
4110c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4120c7d97c5SJed Brown   PetscErrorCode ierr;
4130c7d97c5SJed Brown 
4140c7d97c5SJed Brown   PetscFunctionBegin;
415785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
416785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4170c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4180c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
419785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4207d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4210c7d97c5SJed Brown   PetscFunctionReturn(0);
4220c7d97c5SJed Brown }
4230c7d97c5SJed Brown 
4240c7d97c5SJed Brown #undef __FUNCT__
42582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4269c0446d6SStefano Zampini /*@
42782ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4289c0446d6SStefano Zampini 
429785d1243SStefano Zampini    Collective
43053cdbc3dSStefano Zampini 
43153cdbc3dSStefano Zampini    Input Parameters:
43253cdbc3dSStefano Zampini +  pc - the preconditioning context
43382ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
43453cdbc3dSStefano Zampini 
43553cdbc3dSStefano Zampini    Level: intermediate
43653cdbc3dSStefano Zampini 
4379c0446d6SStefano Zampini    Notes:
43853cdbc3dSStefano Zampini 
4390f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
44053cdbc3dSStefano Zampini @*/
44182ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4420c7d97c5SJed Brown {
4430c7d97c5SJed Brown   PetscErrorCode ierr;
4440c7d97c5SJed Brown 
4450c7d97c5SJed Brown   PetscFunctionBegin;
4460c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4470c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
44882ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
44982ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4500c7d97c5SJed Brown   PetscFunctionReturn(0);
4510c7d97c5SJed Brown }
4520c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4530c7d97c5SJed Brown 
4540c7d97c5SJed Brown #undef __FUNCT__
4550c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
45653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4570c7d97c5SJed Brown {
4580c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
45953cdbc3dSStefano Zampini   PetscErrorCode ierr;
4600c7d97c5SJed Brown 
4610c7d97c5SJed Brown   PetscFunctionBegin;
462785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
463785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
46453cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
46536e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
46636e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
467fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4680c7d97c5SJed Brown   PetscFunctionReturn(0);
4690c7d97c5SJed Brown }
4701e6b0712SBarry Smith 
4710c7d97c5SJed Brown #undef __FUNCT__
4720c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
47357527edcSJed Brown /*@
47428509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
47557527edcSJed Brown 
476785d1243SStefano Zampini    Collective
47757527edcSJed Brown 
47857527edcSJed Brown    Input Parameters:
47957527edcSJed Brown +  pc - the preconditioning context
480785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
48157527edcSJed Brown 
48257527edcSJed Brown    Level: intermediate
48357527edcSJed Brown 
4840f202f7eSStefano Zampini    Notes:
4850f202f7eSStefano Zampini      Any process can list any global node
48657527edcSJed Brown 
4870f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
48857527edcSJed Brown @*/
48953cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4900c7d97c5SJed Brown {
4910c7d97c5SJed Brown   PetscErrorCode ierr;
4920c7d97c5SJed Brown 
4930c7d97c5SJed Brown   PetscFunctionBegin;
4940c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
495674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
496785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
49753cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
49853cdbc3dSStefano Zampini   PetscFunctionReturn(0);
49953cdbc3dSStefano Zampini }
50053cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5011e6b0712SBarry Smith 
50253cdbc3dSStefano Zampini #undef __FUNCT__
50382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
50482ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5050c7d97c5SJed Brown {
5060c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5070c7d97c5SJed Brown   PetscErrorCode ierr;
5080c7d97c5SJed Brown 
5090c7d97c5SJed Brown   PetscFunctionBegin;
510785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
511785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5120c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5130c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
514785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5157d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5160c7d97c5SJed Brown   PetscFunctionReturn(0);
5170c7d97c5SJed Brown }
5180c7d97c5SJed Brown 
5190c7d97c5SJed Brown #undef __FUNCT__
52082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5210c7d97c5SJed Brown /*@
52282ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5230c7d97c5SJed Brown 
524785d1243SStefano Zampini    Collective
5250c7d97c5SJed Brown 
5260c7d97c5SJed Brown    Input Parameters:
5270c7d97c5SJed Brown +  pc - the preconditioning context
52882ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5290c7d97c5SJed Brown 
5300c7d97c5SJed Brown    Level: intermediate
5310c7d97c5SJed Brown 
5320c7d97c5SJed Brown    Notes:
5330c7d97c5SJed Brown 
5340f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
5350c7d97c5SJed Brown @*/
53682ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5370c7d97c5SJed Brown {
5380c7d97c5SJed Brown   PetscErrorCode ierr;
5390c7d97c5SJed Brown 
5400c7d97c5SJed Brown   PetscFunctionBegin;
5410c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5420c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
54382ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54482ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
54553cdbc3dSStefano Zampini   PetscFunctionReturn(0);
54653cdbc3dSStefano Zampini }
54753cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
54853cdbc3dSStefano Zampini 
54953cdbc3dSStefano Zampini #undef __FUNCT__
550da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
551da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
552da1bb401SStefano Zampini {
553da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
554da1bb401SStefano Zampini 
555da1bb401SStefano Zampini   PetscFunctionBegin;
556da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
557da1bb401SStefano Zampini   PetscFunctionReturn(0);
558da1bb401SStefano Zampini }
5591e6b0712SBarry Smith 
560da1bb401SStefano Zampini #undef __FUNCT__
561da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
562da1bb401SStefano Zampini /*@
563785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
564da1bb401SStefano Zampini 
565785d1243SStefano Zampini    Collective
566785d1243SStefano Zampini 
567785d1243SStefano Zampini    Input Parameters:
568785d1243SStefano Zampini .  pc - the preconditioning context
569785d1243SStefano Zampini 
570785d1243SStefano Zampini    Output Parameters:
571785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
572785d1243SStefano Zampini 
573785d1243SStefano Zampini    Level: intermediate
574785d1243SStefano Zampini 
5750f202f7eSStefano Zampini    Notes:
5760f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
577785d1243SStefano Zampini 
578785d1243SStefano Zampini .seealso: PCBDDC
579785d1243SStefano Zampini @*/
580785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
581785d1243SStefano Zampini {
582785d1243SStefano Zampini   PetscErrorCode ierr;
583785d1243SStefano Zampini 
584785d1243SStefano Zampini   PetscFunctionBegin;
585785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
586785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
587785d1243SStefano Zampini   PetscFunctionReturn(0);
588785d1243SStefano Zampini }
589785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
590785d1243SStefano Zampini 
591785d1243SStefano Zampini #undef __FUNCT__
592785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
593785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
594785d1243SStefano Zampini {
595785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
596785d1243SStefano Zampini 
597785d1243SStefano Zampini   PetscFunctionBegin;
598785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
599785d1243SStefano Zampini   PetscFunctionReturn(0);
600785d1243SStefano Zampini }
601785d1243SStefano Zampini 
602785d1243SStefano Zampini #undef __FUNCT__
60382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
604da1bb401SStefano Zampini /*@
60582ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
606da1bb401SStefano Zampini 
607785d1243SStefano Zampini    Collective
608da1bb401SStefano Zampini 
609da1bb401SStefano Zampini    Input Parameters:
61028509bceSStefano Zampini .  pc - the preconditioning context
611da1bb401SStefano Zampini 
612da1bb401SStefano Zampini    Output Parameters:
61328509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
614da1bb401SStefano Zampini 
615da1bb401SStefano Zampini    Level: intermediate
616da1bb401SStefano Zampini 
617da1bb401SStefano Zampini    Notes:
6180f202f7eSStefano 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).
6190f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
620da1bb401SStefano Zampini 
621da1bb401SStefano Zampini .seealso: PCBDDC
622da1bb401SStefano Zampini @*/
62382ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
624da1bb401SStefano Zampini {
625da1bb401SStefano Zampini   PetscErrorCode ierr;
626da1bb401SStefano Zampini 
627da1bb401SStefano Zampini   PetscFunctionBegin;
628da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
62982ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
630da1bb401SStefano Zampini   PetscFunctionReturn(0);
631da1bb401SStefano Zampini }
632da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6331e6b0712SBarry Smith 
634da1bb401SStefano Zampini #undef __FUNCT__
63553cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
63653cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
63753cdbc3dSStefano Zampini {
63853cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
63953cdbc3dSStefano Zampini 
64053cdbc3dSStefano Zampini   PetscFunctionBegin;
64153cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
64253cdbc3dSStefano Zampini   PetscFunctionReturn(0);
64353cdbc3dSStefano Zampini }
6441e6b0712SBarry Smith 
64553cdbc3dSStefano Zampini #undef __FUNCT__
64653cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
64753cdbc3dSStefano Zampini /*@
648785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
64953cdbc3dSStefano Zampini 
650785d1243SStefano Zampini    Collective
651785d1243SStefano Zampini 
652785d1243SStefano Zampini    Input Parameters:
653785d1243SStefano Zampini .  pc - the preconditioning context
654785d1243SStefano Zampini 
655785d1243SStefano Zampini    Output Parameters:
656785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
657785d1243SStefano Zampini 
658785d1243SStefano Zampini    Level: intermediate
659785d1243SStefano Zampini 
6600f202f7eSStefano Zampini    Notes:
6610f202f7eSStefano Zampini      The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
662785d1243SStefano Zampini 
663785d1243SStefano Zampini .seealso: PCBDDC
664785d1243SStefano Zampini @*/
665785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
666785d1243SStefano Zampini {
667785d1243SStefano Zampini   PetscErrorCode ierr;
668785d1243SStefano Zampini 
669785d1243SStefano Zampini   PetscFunctionBegin;
670785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
671785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
672785d1243SStefano Zampini   PetscFunctionReturn(0);
673785d1243SStefano Zampini }
674785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
675785d1243SStefano Zampini 
676785d1243SStefano Zampini #undef __FUNCT__
677785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
678785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
679785d1243SStefano Zampini {
680785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
681785d1243SStefano Zampini 
682785d1243SStefano Zampini   PetscFunctionBegin;
683785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
684785d1243SStefano Zampini   PetscFunctionReturn(0);
685785d1243SStefano Zampini }
686785d1243SStefano Zampini 
687785d1243SStefano Zampini #undef __FUNCT__
68882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
68953cdbc3dSStefano Zampini /*@
69082ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
69153cdbc3dSStefano Zampini 
692785d1243SStefano Zampini    Collective
69353cdbc3dSStefano Zampini 
69453cdbc3dSStefano Zampini    Input Parameters:
69528509bceSStefano Zampini .  pc - the preconditioning context
69653cdbc3dSStefano Zampini 
69753cdbc3dSStefano Zampini    Output Parameters:
69828509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
69953cdbc3dSStefano Zampini 
70053cdbc3dSStefano Zampini    Level: intermediate
70153cdbc3dSStefano Zampini 
70253cdbc3dSStefano Zampini    Notes:
7030f202f7eSStefano 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).
7040f202f7eSStefano Zampini           In the latter case, the IS will be available after PCSetUp.
70553cdbc3dSStefano Zampini 
70653cdbc3dSStefano Zampini .seealso: PCBDDC
70753cdbc3dSStefano Zampini @*/
70882ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
70953cdbc3dSStefano Zampini {
71053cdbc3dSStefano Zampini   PetscErrorCode ierr;
71153cdbc3dSStefano Zampini 
71253cdbc3dSStefano Zampini   PetscFunctionBegin;
71353cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
71482ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7150c7d97c5SJed Brown   PetscFunctionReturn(0);
7160c7d97c5SJed Brown }
71736e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7181e6b0712SBarry Smith 
71936e030ebSStefano Zampini #undef __FUNCT__
720da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7211a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
72236e030ebSStefano Zampini {
72336e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
724da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
725da1bb401SStefano Zampini   PetscErrorCode ierr;
72636e030ebSStefano Zampini 
72736e030ebSStefano Zampini   PetscFunctionBegin;
728674ae819SStefano Zampini   /* free old CSR */
729674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
730fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
731674ae819SStefano Zampini   /* get CSR into graph structure */
732da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
733854ce69bSBarry Smith     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
734785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
735674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
736674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
737da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7381a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7391a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
740674ae819SStefano Zampini   }
741575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
74236e030ebSStefano Zampini   PetscFunctionReturn(0);
74336e030ebSStefano Zampini }
7441e6b0712SBarry Smith 
74536e030ebSStefano Zampini #undef __FUNCT__
746da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
74736e030ebSStefano Zampini /*@
7480f202f7eSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix
74936e030ebSStefano Zampini 
75036e030ebSStefano Zampini    Not collective
75136e030ebSStefano Zampini 
75236e030ebSStefano Zampini    Input Parameters:
75336e030ebSStefano Zampini +  pc - the preconditioning context
7540f202f7eSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the size of the local problem)
75528509bceSStefano Zampini .  xadj, adjncy - the CSR graph
75628509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
75736e030ebSStefano Zampini 
75836e030ebSStefano Zampini    Level: intermediate
75936e030ebSStefano Zampini 
76036e030ebSStefano Zampini    Notes:
76136e030ebSStefano Zampini 
76228509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
76336e030ebSStefano Zampini @*/
7641a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
76536e030ebSStefano Zampini {
766575ad6abSStefano Zampini   void (*f)(void) = 0;
76736e030ebSStefano Zampini   PetscErrorCode ierr;
76836e030ebSStefano Zampini 
76936e030ebSStefano Zampini   PetscFunctionBegin;
77036e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
771674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
772d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
773674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
774674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
775da1bb401SStefano Zampini   }
77636e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
777575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
778575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
779575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
780575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
781575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
78236e030ebSStefano Zampini   }
78336e030ebSStefano Zampini   PetscFunctionReturn(0);
78436e030ebSStefano Zampini }
7859c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7861e6b0712SBarry Smith 
7879c0446d6SStefano Zampini #undef __FUNCT__
78863602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
78963602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
79063602bcaSStefano Zampini {
79163602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
79263602bcaSStefano Zampini   PetscInt i;
79363602bcaSStefano Zampini   PetscErrorCode ierr;
79463602bcaSStefano Zampini 
79563602bcaSStefano Zampini   PetscFunctionBegin;
79663602bcaSStefano Zampini   /* Destroy ISes if they were already set */
79763602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
79863602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
79963602bcaSStefano Zampini   }
80063602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
80163602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
80263602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
80363602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
80463602bcaSStefano Zampini   }
80563602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
80663602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
80763602bcaSStefano Zampini   /* allocate space then set */
808d02579f5SStefano Zampini   if (n_is) {
809d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
810d02579f5SStefano Zampini   }
81163602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
81263602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
81363602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
81463602bcaSStefano Zampini   }
81563602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
81663602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8171a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
81863602bcaSStefano Zampini   PetscFunctionReturn(0);
81963602bcaSStefano Zampini }
82063602bcaSStefano Zampini 
82163602bcaSStefano Zampini #undef __FUNCT__
82263602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
82363602bcaSStefano Zampini /*@
82463602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
82563602bcaSStefano Zampini 
82663602bcaSStefano Zampini    Collective
82763602bcaSStefano Zampini 
82863602bcaSStefano Zampini    Input Parameters:
82963602bcaSStefano Zampini +  pc - the preconditioning context
8300f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
8310f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in local ordering
83263602bcaSStefano Zampini 
83363602bcaSStefano Zampini    Level: intermediate
83463602bcaSStefano Zampini 
8350f202f7eSStefano Zampini    Notes:
8360f202f7eSStefano Zampini      n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
83763602bcaSStefano Zampini 
83863602bcaSStefano Zampini .seealso: PCBDDC
83963602bcaSStefano Zampini @*/
84063602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
84163602bcaSStefano Zampini {
84263602bcaSStefano Zampini   PetscInt       i;
84363602bcaSStefano Zampini   PetscErrorCode ierr;
84463602bcaSStefano Zampini 
84563602bcaSStefano Zampini   PetscFunctionBegin;
84663602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
84763602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
84863602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
84963602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
85063602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
85163602bcaSStefano Zampini   }
852e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
85363602bcaSStefano Zampini   PetscFunctionReturn(0);
85463602bcaSStefano Zampini }
85563602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
85663602bcaSStefano Zampini 
85763602bcaSStefano Zampini #undef __FUNCT__
8589c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8599c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8609c0446d6SStefano Zampini {
8619c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8629c0446d6SStefano Zampini   PetscInt i;
8639c0446d6SStefano Zampini   PetscErrorCode ierr;
8649c0446d6SStefano Zampini 
8659c0446d6SStefano Zampini   PetscFunctionBegin;
866da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8679c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8689c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8699c0446d6SStefano Zampini   }
870d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
87163602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
87263602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
87363602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
87463602bcaSStefano Zampini   }
87563602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
87663602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
877da1bb401SStefano Zampini   /* allocate space then set */
878d02579f5SStefano Zampini   if (n_is) {
879785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
880d02579f5SStefano Zampini   }
8819c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
882da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
883da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8849c0446d6SStefano Zampini   }
8859c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
88663602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8871a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8889c0446d6SStefano Zampini   PetscFunctionReturn(0);
8899c0446d6SStefano Zampini }
8901e6b0712SBarry Smith 
8919c0446d6SStefano Zampini #undef __FUNCT__
8929c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8939c0446d6SStefano Zampini /*@
89463602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8959c0446d6SStefano Zampini 
89663602bcaSStefano Zampini    Collective
8979c0446d6SStefano Zampini 
8989c0446d6SStefano Zampini    Input Parameters:
8999c0446d6SStefano Zampini +  pc - the preconditioning context
9000f202f7eSStefano Zampini .  n_is - number of index sets defining the fields
9010f202f7eSStefano Zampini -  ISForDofs - array of IS describing the fields in global ordering
9029c0446d6SStefano Zampini 
9039c0446d6SStefano Zampini    Level: intermediate
9049c0446d6SStefano Zampini 
9050f202f7eSStefano Zampini    Notes:
9060f202f7eSStefano Zampini      Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
9079c0446d6SStefano Zampini 
9089c0446d6SStefano Zampini .seealso: PCBDDC
9099c0446d6SStefano Zampini @*/
9109c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9119c0446d6SStefano Zampini {
9122b510759SStefano Zampini   PetscInt       i;
9139c0446d6SStefano Zampini   PetscErrorCode ierr;
9149c0446d6SStefano Zampini 
9159c0446d6SStefano Zampini   PetscFunctionBegin;
9169c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
91763602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9182b510759SStefano Zampini   for (i=0;i<n_is;i++) {
91963602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
92063602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9212b510759SStefano Zampini   }
9229c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9239c0446d6SStefano Zampini   PetscFunctionReturn(0);
9249c0446d6SStefano Zampini }
925906d46d4SStefano Zampini 
926da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
927534831adSStefano Zampini #undef __FUNCT__
928534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
929534831adSStefano Zampini /* -------------------------------------------------------------------------- */
930534831adSStefano Zampini /*
931534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
932534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9339c0446d6SStefano Zampini 
934534831adSStefano Zampini    Input Parameter:
935534831adSStefano Zampini +  pc - the preconditioner contex
936534831adSStefano Zampini 
937534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
938534831adSStefano Zampini 
939534831adSStefano Zampini    Notes:
940534831adSStefano Zampini      The interface routine PCPreSolve() is not usually called directly by
941534831adSStefano Zampini    the user, but instead is called by KSPSolve().
942534831adSStefano Zampini */
943534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
944534831adSStefano Zampini {
945534831adSStefano Zampini   PetscErrorCode ierr;
946534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
947534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9483972b0daSStefano Zampini   Vec            used_vec;
9498d00608fSStefano Zampini   PetscBool      copy_rhs = PETSC_TRUE;
950534831adSStefano Zampini 
951534831adSStefano Zampini   PetscFunctionBegin;
95285c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
95385c4d303SStefano Zampini   if (ksp) {
95485c4d303SStefano Zampini     PetscBool iscg;
95585c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
95685c4d303SStefano Zampini     if (!iscg) {
95785c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
95885c4d303SStefano Zampini     }
95985c4d303SStefano Zampini   }
96085c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
96162a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
96262a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
96362a6ff1dSStefano Zampini   }
96462a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
96562a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
96662a6ff1dSStefano Zampini   }
9678d00608fSStefano Zampini 
9683972b0daSStefano Zampini   if (x) {
9693972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9703972b0daSStefano Zampini     used_vec = x;
9718d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
9723972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9733972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9743972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9753972b0daSStefano Zampini   }
9768efcfb23SStefano Zampini 
9778efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
9783972b0daSStefano Zampini   if (ksp) {
979a0cb1b98SStefano Zampini     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
9808efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
9818efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
9823972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9833972b0daSStefano Zampini     }
9843972b0daSStefano Zampini   }
9853308cffdSStefano Zampini 
9868d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
9873972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
988a07ea27aSStefano Zampini   if (rhs) {
9893975b054SStefano Zampini     IS dirIS = NULL;
9903975b054SStefano Zampini 
991a07ea27aSStefano Zampini     /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
9923975b054SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
9933975b054SStefano Zampini     if (dirIS) {
994906d46d4SStefano Zampini       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
995785d1243SStefano Zampini       PetscInt          dirsize,i,*is_indices;
9962b095fd8SStefano Zampini       PetscScalar       *array_x;
9972b095fd8SStefano Zampini       const PetscScalar *array_diagonal;
998785d1243SStefano Zampini 
9993972b0daSStefano Zampini       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
10003972b0daSStefano Zampini       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
1001e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1002e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1003e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1004e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
100582ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
10063972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10072b095fd8SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10083972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10092fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
10103972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10112b095fd8SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10123972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
1013e176bc59SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1014e176bc59SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10158d00608fSStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
10161b968477SStefano Zampini       ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
10178efcfb23SStefano Zampini     }
1018a07ea27aSStefano Zampini   }
1019b76ba322SStefano Zampini 
10208efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
10218d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
10228d00608fSStefano Zampini     /* store the original rhs */
10238d00608fSStefano Zampini     if (copy_rhs) {
10248d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10258d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10268d00608fSStefano Zampini     }
10278d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10283972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10293972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10303972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10318efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10327acc28cbSStefano Zampini     if (ksp) {
10337acc28cbSStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
10347acc28cbSStefano Zampini     }
10353308cffdSStefano Zampini   }
10368efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1037b76ba322SStefano Zampini 
103806f24817SStefano Zampini   /* If using the benign trick, change rhs on pressures (iteration matrix is surely of type MATIS) */
103906f24817SStefano Zampini   if (pcbddc->benign_saddle_point) {
104006f24817SStefano Zampini     Mat_IS *matis = (Mat_IS*)(pc->mat->data);
104106f24817SStefano Zampini 
104206f24817SStefano Zampini     /* store the original rhs */
104306f24817SStefano Zampini     if (copy_rhs) {
104406f24817SStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
104506f24817SStefano Zampini       copy_rhs = PETSC_FALSE;
104606f24817SStefano Zampini     }
104706f24817SStefano Zampini 
104806f24817SStefano Zampini     /* now change (locally) the basis */
104906f24817SStefano Zampini     ierr = VecScatterBegin(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
105006f24817SStefano Zampini     ierr = VecScatterEnd(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
105106f24817SStefano Zampini     if (pcbddc->benign_change) {
105206f24817SStefano Zampini       ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
105306f24817SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
105406f24817SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
105506f24817SStefano Zampini       /* change local iteration matrix */
105606f24817SStefano Zampini       ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_REUSE_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
105706f24817SStefano Zampini       ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr);
105806f24817SStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
105906f24817SStefano Zampini       pcbddc->benign_original_mat = matis->A;
106006f24817SStefano Zampini       ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
106106f24817SStefano Zampini       ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr);
106206f24817SStefano Zampini       matis->A = pcbddc->local_mat;
106306f24817SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
106406f24817SStefano Zampini     } else {
106506f24817SStefano Zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
106606f24817SStefano Zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
106706f24817SStefano Zampini     }
106806f24817SStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
106906f24817SStefano Zampini   }
107006f24817SStefano Zampini 
1071b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
1072b097fa66SStefano Zampini   if (x && pcbddc->use_exact_dirichlet_trick) {
10738efcfb23SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1074b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1075b76ba322SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1076b76ba322SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
10778efcfb23SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10788efcfb23SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1079a0cb1b98SStefano Zampini     if (ksp) {
1080b76ba322SStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1081b76ba322SStefano Zampini     }
1082b76ba322SStefano Zampini   }
1083b76ba322SStefano Zampini 
1084b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1085906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1086906d46d4SStefano Zampini 
1087906d46d4SStefano Zampini     /* get change ctx */
1088906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1089906d46d4SStefano Zampini 
1090906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1091906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1092906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1093906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1094906d46d4SStefano Zampini 
1095906d46d4SStefano Zampini     /* change iteration matrix */
1096906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1097906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1098906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1099906d46d4SStefano Zampini 
11008d00608fSStefano Zampini     /* store the original rhs */
11018d00608fSStefano Zampini     if (copy_rhs) {
11028d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
11038d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
11048d00608fSStefano Zampini     }
11058d00608fSStefano Zampini 
1106906d46d4SStefano Zampini     /* change rhs */
1107906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1108906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
11098d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
111092e3dcfbSStefano Zampini   }
111192e3dcfbSStefano Zampini 
111292e3dcfbSStefano Zampini   /* remove nullspace if present */
11138efcfb23SStefano Zampini   if (ksp && x && pcbddc->NullSpace) {
11148efcfb23SStefano Zampini     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
11158d00608fSStefano Zampini     /* store the original rhs */
11168d00608fSStefano Zampini     if (copy_rhs) {
11178d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
11188d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
11198d00608fSStefano Zampini     }
11208d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
1121d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1122b76ba322SStefano Zampini   }
1123534831adSStefano Zampini   PetscFunctionReturn(0);
1124534831adSStefano Zampini }
1125906d46d4SStefano Zampini 
1126534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1127534831adSStefano Zampini #undef __FUNCT__
1128534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1129534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1130534831adSStefano Zampini /*
1131534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1132534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1133534831adSStefano Zampini 
1134534831adSStefano Zampini    Input Parameter:
1135534831adSStefano Zampini +  pc - the preconditioner contex
1136534831adSStefano Zampini 
1137534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1138534831adSStefano Zampini 
1139534831adSStefano Zampini    Notes:
1140534831adSStefano Zampini      The interface routine PCPostSolve() is not usually called directly by
1141534831adSStefano Zampini      the user, but instead is called by KSPSolve().
1142534831adSStefano Zampini */
1143534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1144534831adSStefano Zampini {
1145534831adSStefano Zampini   PetscErrorCode ierr;
1146534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1147534831adSStefano Zampini 
1148534831adSStefano Zampini   PetscFunctionBegin;
1149b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1150906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1151906d46d4SStefano Zampini 
1152906d46d4SStefano Zampini     /* get change ctx */
1153906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1154906d46d4SStefano Zampini 
1155906d46d4SStefano Zampini     /* restore iteration matrix */
1156906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1157906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1158906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1159906d46d4SStefano Zampini 
1160906d46d4SStefano Zampini     /* get solution in original basis */
1161906d46d4SStefano Zampini     if (x) {
1162906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
116306f24817SStefano Zampini 
1164906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
11652d3346b4SStefano Zampini       /* restore solution on pressures */
116606a4e24aSStefano Zampini       if (pcbddc->benign_saddle_point) {
11672d3346b4SStefano Zampini         Mat_IS *matis = (Mat_IS*)pc->mat->data;
11682d3346b4SStefano Zampini 
116906f24817SStefano Zampini         ierr = VecScatterBegin(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
117006f24817SStefano Zampini         ierr = VecScatterEnd(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11712d3346b4SStefano Zampini         if (pcbddc->benign_change) {
11722d3346b4SStefano Zampini 
11732d3346b4SStefano Zampini           ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr);
117406f24817SStefano Zampini           ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
117506f24817SStefano Zampini           ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11762d3346b4SStefano Zampini         } else {
117706f24817SStefano Zampini           ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
117806f24817SStefano Zampini           ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11792d3346b4SStefano Zampini         }
11802d3346b4SStefano Zampini       }
1181906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
11823425bc38SStefano Zampini     }
1183534831adSStefano Zampini   }
1184906d46d4SStefano Zampini 
11853972b0daSStefano Zampini   /* add solution removed in presolve */
11866bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
11873425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
11883425bc38SStefano Zampini   }
1189906d46d4SStefano Zampini 
1190fb223d50SStefano Zampini   /* restore rhs to its original state */
11918d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
1192fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1193fb223d50SStefano Zampini   }
11948d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
11958efcfb23SStefano Zampini 
11968efcfb23SStefano Zampini   /* restore ksp guess state */
11978efcfb23SStefano Zampini   if (ksp) {
11988efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
11998efcfb23SStefano Zampini   }
1200534831adSStefano Zampini   PetscFunctionReturn(0);
1201534831adSStefano Zampini }
1202534831adSStefano Zampini /* -------------------------------------------------------------------------- */
120353cdbc3dSStefano Zampini #undef __FUNCT__
120453cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
12050c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
12060c7d97c5SJed Brown /*
12070c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
12080c7d97c5SJed Brown                   by setting data structures and options.
12090c7d97c5SJed Brown 
12100c7d97c5SJed Brown    Input Parameter:
121153cdbc3dSStefano Zampini +  pc - the preconditioner context
12120c7d97c5SJed Brown 
12130c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
12140c7d97c5SJed Brown 
12150c7d97c5SJed Brown    Notes:
12160c7d97c5SJed Brown      The interface routine PCSetUp() is not usually called directly by
12170c7d97c5SJed Brown      the user, but instead is called by PCApply() if necessary.
12180c7d97c5SJed Brown */
121953cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
12200c7d97c5SJed Brown {
12210c7d97c5SJed Brown   PetscErrorCode ierr;
12220c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
12235e8657edSStefano Zampini   Mat_IS*        matis;
122408122e43SStefano Zampini   MatNullSpace   nearnullspace;
122591e8d312SStefano Zampini   PetscInt       nrows,ncols;
122608122e43SStefano Zampini   PetscBool      computetopography,computesolvers,computesubschurs;
12278de1fae6SStefano Zampini   PetscBool      computeconstraintsmatrix;
12285e8657edSStefano Zampini   PetscBool      new_nearnullspace_provided,ismatis;
12290c7d97c5SJed Brown 
12300c7d97c5SJed Brown   PetscFunctionBegin;
12315e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
12325e8657edSStefano Zampini   if (!ismatis) {
12335e8657edSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
12345e8657edSStefano Zampini   }
123591e8d312SStefano Zampini   ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr);
123691e8d312SStefano Zampini   if (nrows != ncols) {
123791e8d312SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
123891e8d312SStefano Zampini   }
12395e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1240f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
12413b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1242fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1243f4ddd8eeSStefano Zampini   /* split work */
1244674ae819SStefano Zampini   if (pc->setupcalled) {
1245d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1246674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1247674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1248674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1249674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1250674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1251674ae819SStefano Zampini     }
1252674ae819SStefano Zampini   } else {
1253674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1254674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1255674ae819SStefano Zampini   }
1256fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1257fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1258fb180af4SStefano Zampini   }
12598de1fae6SStefano Zampini   computeconstraintsmatrix = PETSC_FALSE;
12605a95e1ceSStefano Zampini   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
12615a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
12625a95e1ceSStefano Zampini   }
12635a95e1ceSStefano Zampini   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
1264862806e4SStefano Zampini   if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1265862806e4SStefano Zampini 
12665a95e1ceSStefano Zampini   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
12676816873aSStefano Zampini   if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) {
12686816873aSStefano 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");
12696816873aSStefano Zampini   }
12702d3346b4SStefano Zampini 
12712d3346b4SStefano Zampini   /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */
12722d3346b4SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr);
127306a4e24aSStefano Zampini   if (pcbddc->benign_saddle_point && !ismatis) {
12742d3346b4SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS");
12752d3346b4SStefano Zampini   }
12762d3346b4SStefano Zampini 
1277f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
127870cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
127970cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
128058a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1281f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1282f4ddd8eeSStefano Zampini     }
128358a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1284f4ddd8eeSStefano Zampini   }
1285f4ddd8eeSStefano Zampini 
12865e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
12875e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
12885e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
12895e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
12905e8657edSStefano Zampini   } else {
1291b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
129206a4e24aSStefano Zampini     if (!pcbddc->benign_saddle_point) {
12935e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
12945e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
1295c263805aSStefano Zampini     } else { /* TODO: handle user change of basis */
129681d14e9dSStefano Zampini       PetscInt  nz;
129781d14e9dSStefano Zampini       PetscBool sorted;
1298f5fd1fbfSStefano Zampini 
1299f5fd1fbfSStefano Zampini       ierr = ISDestroy(&pcbddc->zerodiag);CHKERRQ(ierr);
130081d14e9dSStefano Zampini       /* should I also include nonzero pressures? */
1301f5fd1fbfSStefano Zampini       ierr = MatFindZeroDiagonals(matis->A,&pcbddc->zerodiag);CHKERRQ(ierr);
130281d14e9dSStefano Zampini       ierr = ISSorted(pcbddc->zerodiag,&sorted);CHKERRQ(ierr);
130381d14e9dSStefano Zampini       if (!sorted) {
130481d14e9dSStefano Zampini         ierr = ISSort(pcbddc->zerodiag);CHKERRQ(ierr);
1305f5fd1fbfSStefano Zampini       }
130606a4e24aSStefano Zampini       //if (!PetscGlobalRank) printf("ZERODIAG\n");
130706a4e24aSStefano Zampini       //if (!PetscGlobalRank) ISView(pcbddc->zerodiag,NULL);
130881d14e9dSStefano Zampini       ierr = ISGetLocalSize(pcbddc->zerodiag,&nz);CHKERRQ(ierr);
130981d14e9dSStefano Zampini       if (nz) {
131081d14e9dSStefano Zampini         IS                zerodiagc;
131181d14e9dSStefano Zampini         PetscScalar       *array;
1312c263805aSStefano Zampini         const PetscInt    *idxs,*idxsc;
1313c263805aSStefano Zampini         PetscInt          i,n,*nnz;
131481d14e9dSStefano Zampini 
131581d14e9dSStefano Zampini         /* TODO: add check for shared dofs */
131681d14e9dSStefano Zampini         pcbddc->use_change_of_basis = PETSC_TRUE;
131781d14e9dSStefano Zampini         ierr = MatGetLocalSize(matis->A,&n,NULL);CHKERRQ(ierr);
131881d14e9dSStefano Zampini         ierr = ISComplement(pcbddc->zerodiag,0,n,&zerodiagc);CHKERRQ(ierr);
131981d14e9dSStefano Zampini         ierr = ISGetIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
132081d14e9dSStefano Zampini         ierr = ISGetIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
132181d14e9dSStefano Zampini         /* local change of basis for pressures */
132281d14e9dSStefano Zampini         ierr = MatDestroy(&pcbddc->benign_change);CHKERRQ(ierr);
132381d14e9dSStefano Zampini         ierr = MatCreate(PetscObjectComm((PetscObject)matis->A),&pcbddc->benign_change);CHKERRQ(ierr);
132481d14e9dSStefano Zampini         ierr = MatSetType(pcbddc->benign_change,MATAIJ);CHKERRQ(ierr);
132581d14e9dSStefano Zampini         ierr = MatSetSizes(pcbddc->benign_change,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
132681d14e9dSStefano Zampini         ierr = PetscMalloc1(n,&nnz);CHKERRQ(ierr);
132781d14e9dSStefano Zampini         for (i=0;i<n-nz;i++) nnz[idxsc[i]] = 1; /* identity on velocities */
132881d14e9dSStefano Zampini         for (i=0;i<nz-1;i++) nnz[idxs[i]] = 2; /* change on pressures */
132981d14e9dSStefano Zampini         nnz[idxs[nz-1]] = nz; /* last local pressure dof: _0 set */
133081d14e9dSStefano Zampini         ierr = MatSeqAIJSetPreallocation(pcbddc->benign_change,0,nnz);CHKERRQ(ierr);
133181d14e9dSStefano Zampini         ierr = PetscFree(nnz);CHKERRQ(ierr);
133281d14e9dSStefano Zampini         /* set identity on velocities */
133381d14e9dSStefano Zampini         for (i=0;i<n-nz;i++) {
133481d14e9dSStefano Zampini           ierr = MatSetValue(pcbddc->benign_change,idxsc[i],idxsc[i],1.,INSERT_VALUES);CHKERRQ(ierr);
133581d14e9dSStefano Zampini         }
133681d14e9dSStefano Zampini         /* set change on pressures */
133781d14e9dSStefano Zampini         for (i=0;i<nz-1;i++) {
133881d14e9dSStefano Zampini           PetscScalar vals[2];
133981d14e9dSStefano Zampini           PetscInt    cols[2];
134081d14e9dSStefano Zampini 
134181d14e9dSStefano Zampini           /* TODO: add quadrature */
134281d14e9dSStefano Zampini           cols[0] = idxs[i];
134381d14e9dSStefano Zampini           cols[1] = idxs[nz-1];
134481d14e9dSStefano Zampini           vals[0] = 1.;
134581d14e9dSStefano Zampini           vals[1] = 1./nz;
134681d14e9dSStefano Zampini           ierr = MatSetValues(pcbddc->benign_change,1,idxs+i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
134781d14e9dSStefano Zampini         }
134881d14e9dSStefano Zampini         ierr = PetscMalloc1(nz,&array);CHKERRQ(ierr);
134981d14e9dSStefano Zampini         for (i=0;i<nz-1;i++) array[i] = -1.;
135081d14e9dSStefano Zampini         array[nz-1] = 1./nz;
135181d14e9dSStefano Zampini         ierr = MatSetValues(pcbddc->benign_change,1,idxs+nz-1,nz,idxs,array,INSERT_VALUES);CHKERRQ(ierr);
135281d14e9dSStefano Zampini         ierr = PetscFree(array);CHKERRQ(ierr);
135381d14e9dSStefano Zampini         ierr = MatAssemblyBegin(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135481d14e9dSStefano Zampini         ierr = MatAssemblyEnd(pcbddc->benign_change,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135581d14e9dSStefano Zampini         /* TODO: need optimization? */
135681d14e9dSStefano Zampini         ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_INITIAL_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
1357*efc2fbd9SStefano Zampini         /* store local and global idxs for p0 */
1358*efc2fbd9SStefano Zampini         pcbddc->benign_p0_lidx = idxs[nz-1];
1359*efc2fbd9SStefano Zampini         ierr = ISLocalToGlobalMappingApply(pc->pmat->rmap->mapping,1,&idxs[nz-1],&pcbddc->benign_p0_gidx);CHKERRQ(ierr);
1360c263805aSStefano Zampini         ierr = ISRestoreIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
1361c263805aSStefano Zampini         ierr = ISRestoreIndices(zerodiagc,&idxsc);CHKERRQ(ierr);
1362c263805aSStefano Zampini         ierr = ISDestroy(&zerodiagc);CHKERRQ(ierr);
1363c263805aSStefano Zampini         /* pop B0 mat from pcbddc->local_mat */
1364c263805aSStefano Zampini         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1365019a44ceSStefano Zampini       } else { /* this is unlikely to happen but, just in case, destroy the empty IS */
1366019a44ceSStefano Zampini         ierr = ISDestroy(&pcbddc->zerodiag);CHKERRQ(ierr);
136781d14e9dSStefano Zampini         ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
136881d14e9dSStefano Zampini         pcbddc->local_mat = matis->A;
136981d14e9dSStefano Zampini       }
1370f5fd1fbfSStefano Zampini     }
1371674ae819SStefano Zampini   }
1372f4ddd8eeSStefano Zampini 
1373e496cd5dSStefano Zampini   /* workaround for reals */
1374e496cd5dSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
13753301b35fSStefano Zampini   if (matis->A->symmetric_set) {
13763301b35fSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr);
13773301b35fSStefano Zampini   }
1378e496cd5dSStefano Zampini #endif
137906a4e24aSStefano Zampini   if (matis->A->symmetric_set) {
138006a4e24aSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr);
138106a4e24aSStefano Zampini   }
138206a4e24aSStefano Zampini   if (matis->A->spd_set) {
138306a4e24aSStefano Zampini     ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr);
138406a4e24aSStefano Zampini   }
1385e496cd5dSStefano Zampini 
13865e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
13875e8657edSStefano Zampini   {
13885e8657edSStefano Zampini     Mat temp_mat;
13895e8657edSStefano Zampini 
13905e8657edSStefano Zampini     temp_mat = matis->A;
13915e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
13925e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
13935e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
13945e8657edSStefano Zampini     matis->A = temp_mat;
13955e8657edSStefano Zampini   }
1396684f6988SStefano Zampini 
139781d14e9dSStefano Zampini   /* Analyze interface */
1398674ae819SStefano Zampini   if (computetopography) {
1399674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
14008de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1401674ae819SStefano Zampini   }
1402fb8d54d4SStefano Zampini 
140306f24817SStefano Zampini   if (pcbddc->dbg_flag && pcbddc->zerodiag) {
140406f24817SStefano Zampini     PC_IS          *pcis = (PC_IS*)(pc->data);
140506f24817SStefano Zampini     IS             dirIS = NULL;
140606f24817SStefano Zampini     PetscScalar    *vals;
140706f24817SStefano Zampini     const PetscInt *idxs;
140806f24817SStefano Zampini     PetscInt       i,nz;
140906f24817SStefano Zampini 
141006f24817SStefano Zampini     /* p0 */
141106f24817SStefano Zampini     ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
141206f24817SStefano Zampini     ierr = PetscMalloc1(pcis->n,&vals);CHKERRQ(ierr);
141306f24817SStefano Zampini     ierr = ISGetLocalSize(pcbddc->zerodiag,&nz);CHKERRQ(ierr);
141406f24817SStefano Zampini     ierr = ISGetIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
141506f24817SStefano Zampini     for (i=0;i<nz;i++) vals[i] = 1.;
141606f24817SStefano Zampini     ierr = VecSetValues(pcis->vec1_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
141706f24817SStefano Zampini     ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
141806f24817SStefano Zampini     ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
141906f24817SStefano Zampini     /* v_I */
142006f24817SStefano Zampini     ierr = VecSetRandom(pcis->vec2_N,NULL);CHKERRQ(ierr);
142106f24817SStefano Zampini     for (i=0;i<nz;i++) vals[i] = 0.;
142206f24817SStefano Zampini     ierr = VecSetValues(pcis->vec2_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
142306f24817SStefano Zampini     ierr = ISRestoreIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
142406f24817SStefano Zampini     ierr = ISGetIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
142506f24817SStefano Zampini     for (i=0;i<pcis->n_B;i++) vals[i] = 0.;
142606f24817SStefano Zampini     ierr = VecSetValues(pcis->vec2_N,pcis->n_B,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
142706f24817SStefano Zampini     ierr = ISRestoreIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
142806f24817SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
142906f24817SStefano Zampini     if (dirIS) {
143006f24817SStefano Zampini       PetscInt n;
143106f24817SStefano Zampini 
143206f24817SStefano Zampini       ierr = ISGetLocalSize(dirIS,&n);CHKERRQ(ierr);
143306f24817SStefano Zampini       ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
143406f24817SStefano Zampini       for (i=0;i<n;i++) vals[i] = 0.;
143506f24817SStefano Zampini       ierr = VecSetValues(pcis->vec2_N,n,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
143606f24817SStefano Zampini       ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
143706f24817SStefano Zampini     }
143806f24817SStefano Zampini     ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
143906f24817SStefano Zampini     ierr = VecAssemblyBegin(pcis->vec2_N);CHKERRQ(ierr);
144006f24817SStefano Zampini     ierr = VecAssemblyEnd(pcis->vec2_N);CHKERRQ(ierr);
144106f24817SStefano Zampini     ierr = VecSet(matis->x,0.);CHKERRQ(ierr);
144206f24817SStefano Zampini     ierr = MatMult(matis->A,pcis->vec1_N,matis->x);CHKERRQ(ierr);
144306f24817SStefano Zampini     ierr = VecDot(matis->x,pcis->vec2_N,&vals[0]);CHKERRQ(ierr);
144406f24817SStefano Zampini     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] check original mat: %1.4e\n",PetscGlobalRank,vals[0]);CHKERRQ(ierr);
144506f24817SStefano Zampini     ierr = PetscFree(vals);CHKERRQ(ierr);
1446*efc2fbd9SStefano Zampini 
1447*efc2fbd9SStefano Zampini     /* check PCBDDCBenignPopOrPush */
1448*efc2fbd9SStefano Zampini     ierr = VecSetRandom(pcis->vec1_global,NULL);CHKERRQ(ierr);
1449*efc2fbd9SStefano Zampini     pcbddc->benign_p0 = -PetscGlobalRank;
1450*efc2fbd9SStefano Zampini     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr);
1451*efc2fbd9SStefano Zampini     pcbddc->benign_p0 = 1;
1452*efc2fbd9SStefano Zampini     ierr = PCBDDCBenignPopOrPushP0(pc,pcis->vec1_global,PETSC_TRUE);CHKERRQ(ierr);
1453*efc2fbd9SStefano Zampini     if (pcbddc->benign_p0 != -PetscGlobalRank) {
1454*efc2fbd9SStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error testing PCBDDCBenignPopOrPushP0! Found %1.4e instead of %1.4e\n",pcbddc->benign_p0,-PetscGlobalRank);CHKERRQ(ierr);
1455*efc2fbd9SStefano Zampini     }
145606f24817SStefano Zampini   }
145706f24817SStefano Zampini 
1458b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1459684f6988SStefano Zampini   if (computesolvers) {
1460d5574798SStefano Zampini     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1461d5574798SStefano Zampini 
1462d5574798SStefano Zampini     if (computesubschurs && computetopography) {
146308122e43SStefano Zampini       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1464b1b3d7a2SStefano Zampini     }
1465d5574798SStefano Zampini     if (sub_schurs->use_mumps) {
14662070dbb6SStefano Zampini       if (computesubschurs) {
146708122e43SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
14682070dbb6SStefano Zampini       }
1469d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1470d5574798SStefano Zampini     } else {
1471d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
14722070dbb6SStefano Zampini       if (computesubschurs) {
1473d5574798SStefano Zampini         ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1474d5574798SStefano Zampini       }
14752070dbb6SStefano Zampini     }
147608122e43SStefano Zampini     if (pcbddc->adaptive_selection) {
147708122e43SStefano Zampini       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
14788de1fae6SStefano Zampini       computeconstraintsmatrix = PETSC_TRUE;
1479b7eb3628SStefano Zampini     }
1480b1b3d7a2SStefano Zampini   }
1481684f6988SStefano Zampini 
1482f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1483fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1484f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1485f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1486f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1487f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1488f4ddd8eeSStefano Zampini     } else {
1489f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1490f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1491f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1492165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1493f4ddd8eeSStefano Zampini         PetscInt         i;
1494165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1495165b64e2SStefano Zampini         PetscObjectState state;
1496165b64e2SStefano Zampini         PetscInt         nnsp_size;
1497165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1498f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1499f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1500165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1501f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1502f4ddd8eeSStefano Zampini             break;
1503f4ddd8eeSStefano Zampini           }
1504f4ddd8eeSStefano Zampini         }
1505f4ddd8eeSStefano Zampini       }
1506f4ddd8eeSStefano Zampini     }
1507f4ddd8eeSStefano Zampini   } else {
1508f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1509f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1510f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1511f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1512f4ddd8eeSStefano Zampini     }
1513f4ddd8eeSStefano Zampini   }
1514f4ddd8eeSStefano Zampini 
1515f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1516727cdba6SStefano Zampini   /* reset primal space flags */
1517f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1518727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
15198de1fae6SStefano Zampini   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1520727cdba6SStefano Zampini     /* It also sets the primal space flags */
1521674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1522e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1523f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
15243975b054SStefano Zampini   }
15255e8657edSStefano Zampini 
15263975b054SStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
15275e8657edSStefano Zampini     if (pcbddc->use_change_of_basis) {
15285e8657edSStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1529c263805aSStefano Zampini       Mat   temp_mat = NULL;
15305e8657edSStefano Zampini 
1531c263805aSStefano Zampini       if (pcbddc->zerodiag) {
1532c263805aSStefano Zampini         /* insert B0 in pcbddc->local_mat */
1533c263805aSStefano Zampini         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr);
153406f24817SStefano Zampini         if (pcbddc->dbg_flag) {
153506f24817SStefano Zampini           PC_IS          *pcis = (PC_IS*)(pc->data);
153606f24817SStefano Zampini           IS             dirIS = NULL;
153706f24817SStefano Zampini           PetscScalar    *vals;
153806f24817SStefano Zampini           const PetscInt *idxs;
153906f24817SStefano Zampini           PetscInt       i,nz;
154006f24817SStefano Zampini 
154106f24817SStefano Zampini           /* p0 */
154206f24817SStefano Zampini           ierr = VecSet(pcis->vec1_N,0.);CHKERRQ(ierr);
154306f24817SStefano Zampini           ierr = PetscMalloc1(pcis->n,&vals);CHKERRQ(ierr);
154406f24817SStefano Zampini           ierr = ISGetLocalSize(pcbddc->zerodiag,&nz);CHKERRQ(ierr);
154506f24817SStefano Zampini           ierr = ISGetIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
154606f24817SStefano Zampini           vals[0] = 1.;
154706f24817SStefano Zampini           ierr = VecSetValues(pcis->vec1_N,1,&idxs[nz-1],vals,INSERT_VALUES);CHKERRQ(ierr);
154806f24817SStefano Zampini           ierr = VecAssemblyBegin(pcis->vec1_N);CHKERRQ(ierr);
154906f24817SStefano Zampini           ierr = VecAssemblyEnd(pcis->vec1_N);CHKERRQ(ierr);
155006f24817SStefano Zampini           /* v_I */
155106f24817SStefano Zampini           ierr = VecSetRandom(pcis->vec2_N,NULL);CHKERRQ(ierr);
155206f24817SStefano Zampini           for (i=0;i<nz;i++) vals[i] = 0.;
155306f24817SStefano Zampini           ierr = VecSetValues(pcis->vec2_N,nz,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
155406f24817SStefano Zampini           ierr = ISRestoreIndices(pcbddc->zerodiag,&idxs);CHKERRQ(ierr);
155506f24817SStefano Zampini           ierr = ISGetIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
155606f24817SStefano Zampini           for (i=0;i<pcis->n_B;i++) vals[i] = 0.;
155706f24817SStefano Zampini           ierr = VecSetValues(pcis->vec2_N,pcis->n_B,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
155806f24817SStefano Zampini           ierr = ISRestoreIndices(pcis->is_B_local,&idxs);CHKERRQ(ierr);
155906f24817SStefano Zampini           ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
156006f24817SStefano Zampini           if (dirIS) {
156106f24817SStefano Zampini             PetscInt n;
156206f24817SStefano Zampini 
156306f24817SStefano Zampini             ierr = ISGetLocalSize(dirIS,&n);CHKERRQ(ierr);
156406f24817SStefano Zampini             ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
156506f24817SStefano Zampini             for (i=0;i<n;i++) vals[i] = 0.;
156606f24817SStefano Zampini             ierr = VecSetValues(pcis->vec2_N,n,idxs,vals,INSERT_VALUES);CHKERRQ(ierr);
156706f24817SStefano Zampini             ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
156806f24817SStefano Zampini           }
156906f24817SStefano Zampini           ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
157006f24817SStefano Zampini           ierr = VecAssemblyBegin(pcis->vec2_N);CHKERRQ(ierr);
157106f24817SStefano Zampini           ierr = VecAssemblyEnd(pcis->vec2_N);CHKERRQ(ierr);
157206f24817SStefano Zampini           ierr = VecSet(matis->x,0.);CHKERRQ(ierr);
157306f24817SStefano Zampini           ierr = MatMult(pcbddc->local_mat,pcis->vec1_N,matis->x);CHKERRQ(ierr);
157406f24817SStefano Zampini           ierr = VecDot(matis->x,pcis->vec2_N,&vals[0]);CHKERRQ(ierr);
157506f24817SStefano Zampini           ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] check new mat: %1.4e\n",PetscGlobalRank,vals[0]);CHKERRQ(ierr);
157606f24817SStefano Zampini           ierr = PetscFree(vals);CHKERRQ(ierr);
157706f24817SStefano Zampini         }
1578c263805aSStefano Zampini         /* hack: swap pointers */
1579c263805aSStefano Zampini         temp_mat = matis->A;
1580c263805aSStefano Zampini         matis->A = pcbddc->local_mat;
1581c263805aSStefano Zampini         pcbddc->local_mat = NULL;
1582c263805aSStefano Zampini       }
15835e8657edSStefano Zampini       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
1584c263805aSStefano Zampini       if (pcbddc->zerodiag) {
1585c263805aSStefano Zampini         /* restore original matrix */
1586c263805aSStefano Zampini         ierr = MatDestroy(&matis->A);CHKERRQ(ierr);
1587c263805aSStefano Zampini         matis->A = temp_mat;
1588c263805aSStefano Zampini         /* pop B0 from pcbddc->local_mat */
1589c263805aSStefano Zampini         ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr);
1590c263805aSStefano Zampini       }
15915e8657edSStefano Zampini       /* get submatrices */
15925e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
15935e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
15945e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
15955e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
15965e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
15975e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
15983975b054SStefano Zampini       /* set flag in pcis to not reuse submatrices during PCISCreate */
15993975b054SStefano Zampini       pcis->reusesubmatrices = PETSC_FALSE;
16005e8657edSStefano Zampini     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1601b96c3477SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
16025e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
16035e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
16045e8657edSStefano Zampini     }
1605b96c3477SStefano Zampini     /* SetUp coarse and local Neumann solvers */
160699cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1607b96c3477SStefano Zampini     /* SetUp Scaling operator */
1608674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
16090c7d97c5SJed Brown   }
161070cf5478SStefano Zampini 
161158a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
161258a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
16132b510759SStefano Zampini   }
16140c7d97c5SJed Brown   PetscFunctionReturn(0);
16150c7d97c5SJed Brown }
16160c7d97c5SJed Brown 
16170c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
16180c7d97c5SJed Brown /*
161950efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
16200c7d97c5SJed Brown 
16210c7d97c5SJed Brown    Input Parameters:
16220f202f7eSStefano Zampini +  pc - the preconditioner context
16230f202f7eSStefano Zampini -  r - input vector (global)
16240c7d97c5SJed Brown 
16250c7d97c5SJed Brown    Output Parameter:
16260c7d97c5SJed Brown .  z - output vector (global)
16270c7d97c5SJed Brown 
16280c7d97c5SJed Brown    Application Interface Routine: PCApply()
16290c7d97c5SJed Brown  */
16300c7d97c5SJed Brown #undef __FUNCT__
16310c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
163253cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
16330c7d97c5SJed Brown {
16340c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
16350c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1636b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
16370c7d97c5SJed Brown   PetscErrorCode    ierr;
16383b03a366Sstefano_zampini   const PetscScalar one = 1.0;
16393b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
16402617d88aSStefano Zampini   const PetscScalar zero = 0.0;
16410c7d97c5SJed Brown 
16420c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
16430c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
1644b097fa66SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
16450c7d97c5SJed Brown 
16460c7d97c5SJed Brown   PetscFunctionBegin;
1647*efc2fbd9SStefano Zampini   if (pcbddc->benign_saddle_point) { /* extract p0 from r */
1648*efc2fbd9SStefano Zampini     ierr = PCBDDCBenignPopOrPushP0(pc,r,PETSC_TRUE);CHKERRQ(ierr);
1649*efc2fbd9SStefano Zampini   }
165085c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1651b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
16520c7d97c5SJed Brown     /* First Dirichlet solve */
16530c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16540c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16550c7d97c5SJed Brown     /*
16560c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1657b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1658674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
16590c7d97c5SJed Brown     */
1660b097fa66SStefano Zampini     if (n_D) {
1661b097fa66SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
16620c7d97c5SJed Brown       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
16638eeda7d8SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1664b097fa66SStefano Zampini       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1665b097fa66SStefano Zampini     } else {
1666b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1667b097fa66SStefano Zampini     }
16680c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16690c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1670674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1671b76ba322SStefano Zampini   } else {
1672b097fa66SStefano Zampini     if (pcbddc->switch_static) {
16730bdf917eSStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1674b097fa66SStefano Zampini     }
1675674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1676b76ba322SStefano Zampini   }
1677b76ba322SStefano Zampini 
16782617d88aSStefano Zampini   /* Apply interface preconditioner
16792617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1680dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
16812617d88aSStefano Zampini 
1682674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1683674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
16840c7d97c5SJed Brown 
16853b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
16860c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16870c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1688b097fa66SStefano Zampini   if (n_B) {
16890c7d97c5SJed Brown     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
16908eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1691b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1692b097fa66SStefano Zampini     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1693b097fa66SStefano Zampini   }
1694df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1695*efc2fbd9SStefano Zampini 
1696b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1697b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1698b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1699b097fa66SStefano Zampini     } else {
1700b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1701b097fa66SStefano Zampini     }
17020c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17030c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1704b097fa66SStefano Zampini   } else {
1705b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1706b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1707b097fa66SStefano Zampini     } else {
1708b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1709b097fa66SStefano Zampini     }
1710b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1711b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1712b097fa66SStefano Zampini   }
1713*efc2fbd9SStefano Zampini 
1714*efc2fbd9SStefano Zampini   if (pcbddc->benign_saddle_point) { /* push p0 (computed in PCBDDCApplyInterface) */
1715*efc2fbd9SStefano Zampini     ierr = PCBDDCBenignPopOrPushP0(pc,z,PETSC_FALSE);CHKERRQ(ierr);
1716*efc2fbd9SStefano Zampini   }
17170c7d97c5SJed Brown   PetscFunctionReturn(0);
17180c7d97c5SJed Brown }
171950efa1b5SStefano Zampini 
172050efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
172150efa1b5SStefano Zampini /*
172250efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
172350efa1b5SStefano Zampini 
172450efa1b5SStefano Zampini    Input Parameters:
17250f202f7eSStefano Zampini +  pc - the preconditioner context
17260f202f7eSStefano Zampini -  r - input vector (global)
172750efa1b5SStefano Zampini 
172850efa1b5SStefano Zampini    Output Parameter:
172950efa1b5SStefano Zampini .  z - output vector (global)
173050efa1b5SStefano Zampini 
173150efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
173250efa1b5SStefano Zampini  */
173350efa1b5SStefano Zampini #undef __FUNCT__
173450efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
173550efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
173650efa1b5SStefano Zampini {
173750efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
173850efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1739b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
174050efa1b5SStefano Zampini   PetscErrorCode    ierr;
174150efa1b5SStefano Zampini   const PetscScalar one = 1.0;
174250efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
174350efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
174450efa1b5SStefano Zampini 
174550efa1b5SStefano Zampini   PetscFunctionBegin;
174650efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1747b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
174850efa1b5SStefano Zampini     /* First Dirichlet solve */
174950efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
175050efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
175150efa1b5SStefano Zampini     /*
175250efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
1753b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
175450efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
175550efa1b5SStefano Zampini     */
1756b097fa66SStefano Zampini     if (n_D) {
1757b097fa66SStefano Zampini       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
175850efa1b5SStefano Zampini       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
175950efa1b5SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1760b097fa66SStefano Zampini       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1761b097fa66SStefano Zampini     } else {
1762b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1763b097fa66SStefano Zampini     }
176450efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
176550efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
176650efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
176750efa1b5SStefano Zampini   } else {
1768b097fa66SStefano Zampini     if (pcbddc->switch_static) {
176950efa1b5SStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1770b097fa66SStefano Zampini     }
177150efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
177250efa1b5SStefano Zampini   }
177350efa1b5SStefano Zampini 
177450efa1b5SStefano Zampini   /* Apply interface preconditioner
177550efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1776dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
177750efa1b5SStefano Zampini 
177850efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
177950efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
178050efa1b5SStefano Zampini 
178150efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
178250efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178350efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1784b097fa66SStefano Zampini   if (n_B) {
178550efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
178650efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1787b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1788b097fa66SStefano Zampini     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1789b097fa66SStefano Zampini   }
1790b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1791b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1792b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1793b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1794b097fa66SStefano Zampini     } else {
1795b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1796b097fa66SStefano Zampini     }
179750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
179850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1799b097fa66SStefano Zampini   } else {
1800b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1801b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1802b097fa66SStefano Zampini     } else {
1803b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1804b097fa66SStefano Zampini     }
1805b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1806b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1807b097fa66SStefano Zampini   }
180850efa1b5SStefano Zampini   PetscFunctionReturn(0);
180950efa1b5SStefano Zampini }
1810da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1811674ae819SStefano Zampini 
1812da1bb401SStefano Zampini #undef __FUNCT__
1813da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1814da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1815da1bb401SStefano Zampini {
1816da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1817da1bb401SStefano Zampini   PetscErrorCode ierr;
1818da1bb401SStefano Zampini 
1819da1bb401SStefano Zampini   PetscFunctionBegin;
1820da1bb401SStefano Zampini   /* free data created by PCIS */
1821da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1822674ae819SStefano Zampini   /* free BDDC custom data  */
1823674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1824674ae819SStefano Zampini   /* destroy objects related to topography */
1825674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1826674ae819SStefano Zampini   /* free allocated graph structure */
1827da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1828b96c3477SStefano Zampini   /* free allocated sub schurs structure */
1829b96c3477SStefano Zampini   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
183034a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1831674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
183234a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1833674ae819SStefano Zampini   /* free solvers stuff */
1834674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
183562a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
183662a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
183762a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1838906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1839906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1840906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1841906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1842906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1843906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1844906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1845906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1846906d46d4SStefano Zampini   }
1847906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
18483425bc38SStefano Zampini   /* remove functions */
1849906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1850674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1851bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
18522b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1853b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
18542b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1855bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1856bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
185782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1858bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
185982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1860bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
186182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1862bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1863785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1864bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
186563602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1866bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1867bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1868bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1869bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1870674ae819SStefano Zampini   /* Free the private data structure */
1871674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1872da1bb401SStefano Zampini   PetscFunctionReturn(0);
1873da1bb401SStefano Zampini }
18743425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
18751e6b0712SBarry Smith 
18763425bc38SStefano Zampini #undef __FUNCT__
18773425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
18783425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
18793425bc38SStefano Zampini {
1880674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
1881c08af4c6SStefano Zampini   Vec            copy_standard_rhs;
18823425bc38SStefano Zampini   PC_IS*         pcis;
18833425bc38SStefano Zampini   PC_BDDC*       pcbddc;
18843425bc38SStefano Zampini   PetscErrorCode ierr;
18850c7d97c5SJed Brown 
18863425bc38SStefano Zampini   PetscFunctionBegin;
18873425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
18883425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
18893425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
18903425bc38SStefano Zampini 
1891c08af4c6SStefano Zampini   /*
1892c08af4c6SStefano Zampini      change of basis for physical rhs if needed
1893c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
1894c08af4c6SStefano Zampini      TODO: better management when FETIDP will have its own class
1895c08af4c6SStefano Zampini   */
1896c08af4c6SStefano Zampini   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1897c08af4c6SStefano Zampini   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1898c08af4c6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
18993425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
1900c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1901c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1902fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1903fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
1904c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1905c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1906674ae819SStefano Zampini   /* Apply partition of unity */
19073425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1908c08af4c6SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
19098eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
19103425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
19113425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
19123425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
19133425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1914c08af4c6SStefano Zampini     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1915c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1916c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1917c08af4c6SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1918c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1919c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19203425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
19213425bc38SStefano Zampini   }
1922c08af4c6SStefano Zampini   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
19233425bc38SStefano Zampini   /* BDDC rhs */
19243425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
19258eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
19263425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
19273425bc38SStefano Zampini   }
19283425bc38SStefano Zampini   /* apply BDDC */
1929dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
19303425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
19313425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
19323425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
19333425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19343425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19353425bc38SStefano Zampini   PetscFunctionReturn(0);
19363425bc38SStefano Zampini }
19371e6b0712SBarry Smith 
19383425bc38SStefano Zampini #undef __FUNCT__
19393425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
19403425bc38SStefano Zampini /*@
19410f202f7eSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
19423425bc38SStefano Zampini 
19433425bc38SStefano Zampini    Collective
19443425bc38SStefano Zampini 
19453425bc38SStefano Zampini    Input Parameters:
19460f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
19470f202f7eSStefano Zampini -  standard_rhs    - the right-hand side of the original linear system
19483425bc38SStefano Zampini 
19493425bc38SStefano Zampini    Output Parameters:
19500f202f7eSStefano Zampini .  fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
19513425bc38SStefano Zampini 
19523425bc38SStefano Zampini    Level: developer
19533425bc38SStefano Zampini 
19543425bc38SStefano Zampini    Notes:
19553425bc38SStefano Zampini 
19560f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
19573425bc38SStefano Zampini @*/
19583425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
19593425bc38SStefano Zampini {
1960674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
19613425bc38SStefano Zampini   PetscErrorCode ierr;
19623425bc38SStefano Zampini 
19633425bc38SStefano Zampini   PetscFunctionBegin;
19643425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
19653425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
19663425bc38SStefano Zampini   PetscFunctionReturn(0);
19673425bc38SStefano Zampini }
19683425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
19691e6b0712SBarry Smith 
19703425bc38SStefano Zampini #undef __FUNCT__
19713425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
19723425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
19733425bc38SStefano Zampini {
1974674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
19753425bc38SStefano Zampini   PC_IS*         pcis;
19763425bc38SStefano Zampini   PC_BDDC*       pcbddc;
19773425bc38SStefano Zampini   PetscErrorCode ierr;
19783425bc38SStefano Zampini 
19793425bc38SStefano Zampini   PetscFunctionBegin;
19803425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
19813425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
19823425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
19833425bc38SStefano Zampini 
19843425bc38SStefano Zampini   /* apply B_delta^T */
19853425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19863425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19873425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
19883425bc38SStefano Zampini   /* compute rhs for BDDC application */
19893425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
19908eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
19913425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
19923425bc38SStefano Zampini   }
19933425bc38SStefano Zampini   /* apply BDDC */
1994dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
19953425bc38SStefano Zampini   /* put values into standard global vector */
19963425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19973425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19988eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
19993425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
20003425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
20013425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
20023425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
20033425bc38SStefano Zampini   }
20043425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20053425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20063425bc38SStefano Zampini   /* final change of basis if needed
20073425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
20083308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
20093425bc38SStefano Zampini   PetscFunctionReturn(0);
20103425bc38SStefano Zampini }
20111e6b0712SBarry Smith 
20123425bc38SStefano Zampini #undef __FUNCT__
20133425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
20143425bc38SStefano Zampini /*@
20150f202f7eSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
20163425bc38SStefano Zampini 
20173425bc38SStefano Zampini    Collective
20183425bc38SStefano Zampini 
20193425bc38SStefano Zampini    Input Parameters:
20200f202f7eSStefano Zampini +  fetidp_mat      - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
20210f202f7eSStefano Zampini -  fetidp_flux_sol - the solution of the FETI-DP linear system
20223425bc38SStefano Zampini 
20233425bc38SStefano Zampini    Output Parameters:
20240f202f7eSStefano Zampini .  standard_sol    - the solution defined on the physical domain
20253425bc38SStefano Zampini 
20263425bc38SStefano Zampini    Level: developer
20273425bc38SStefano Zampini 
20283425bc38SStefano Zampini    Notes:
20293425bc38SStefano Zampini 
20300f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
20313425bc38SStefano Zampini @*/
20323425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
20333425bc38SStefano Zampini {
2034674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
20353425bc38SStefano Zampini   PetscErrorCode ierr;
20363425bc38SStefano Zampini 
20373425bc38SStefano Zampini   PetscFunctionBegin;
20383425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
20393425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
20403425bc38SStefano Zampini   PetscFunctionReturn(0);
20413425bc38SStefano Zampini }
20423425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
20431e6b0712SBarry Smith 
2044f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
2045edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
2046f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
2047f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
2048edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
2049f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
2050674ae819SStefano Zampini 
20513425bc38SStefano Zampini #undef __FUNCT__
20523425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
20533425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
20543425bc38SStefano Zampini {
2055674ae819SStefano Zampini 
2056674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
20573425bc38SStefano Zampini   Mat            newmat;
2058674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
20593425bc38SStefano Zampini   PC             newpc;
2060ce94432eSBarry Smith   MPI_Comm       comm;
20613425bc38SStefano Zampini   PetscErrorCode ierr;
20623425bc38SStefano Zampini 
20633425bc38SStefano Zampini   PetscFunctionBegin;
2064ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
20653425bc38SStefano Zampini   /* FETIDP linear matrix */
20663425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
20673425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
20683425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
20693425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
2070edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
20713425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
20723425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
20733425bc38SStefano Zampini   /* FETIDP preconditioner */
20743425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
20753425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
20763425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
20773425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
20783425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
20793425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
2080edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
20813425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
208223ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
20833425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
20843425bc38SStefano Zampini   /* return pointers for objects created */
20853425bc38SStefano Zampini   *fetidp_mat=newmat;
20863425bc38SStefano Zampini   *fetidp_pc=newpc;
20873425bc38SStefano Zampini   PetscFunctionReturn(0);
20883425bc38SStefano Zampini }
20891e6b0712SBarry Smith 
20903425bc38SStefano Zampini #undef __FUNCT__
20913425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
20923425bc38SStefano Zampini /*@
20930f202f7eSStefano Zampini  PCBDDCCreateFETIDPOperators - Create FETI-DP operators
20943425bc38SStefano Zampini 
20953425bc38SStefano Zampini    Collective
20963425bc38SStefano Zampini 
20973425bc38SStefano Zampini    Input Parameters:
20980f202f7eSStefano Zampini .  pc - the BDDC preconditioning context (setup should have been called before)
209928509bceSStefano Zampini 
210028509bceSStefano Zampini    Output Parameters:
21010f202f7eSStefano Zampini +  fetidp_mat - shell FETI-DP matrix object
21020f202f7eSStefano Zampini -  fetidp_pc  - shell Dirichlet preconditioner for FETI-DP matrix
210328509bceSStefano Zampini 
210428509bceSStefano Zampini    Options Database Keys:
21050f202f7eSStefano Zampini .    -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers
21063425bc38SStefano Zampini 
21073425bc38SStefano Zampini    Level: developer
21083425bc38SStefano Zampini 
21093425bc38SStefano Zampini    Notes:
21100f202f7eSStefano Zampini      Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
21113425bc38SStefano Zampini 
21120f202f7eSStefano Zampini .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
21133425bc38SStefano Zampini @*/
21143425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
21153425bc38SStefano Zampini {
21163425bc38SStefano Zampini   PetscErrorCode ierr;
21173425bc38SStefano Zampini 
21183425bc38SStefano Zampini   PetscFunctionBegin;
21193425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
21203425bc38SStefano Zampini   if (pc->setupcalled) {
2121516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
2122f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
21233425bc38SStefano Zampini   PetscFunctionReturn(0);
21243425bc38SStefano Zampini }
21250c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
2126da1bb401SStefano Zampini /*MC
2127da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
21280c7d97c5SJed Brown 
212928509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
213028509bceSStefano Zampini 
213128509bceSStefano Zampini .vb
213228509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
213328509bceSStefano 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
213428509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
21350f202f7eSStefano 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
213628509bceSStefano Zampini .ve
213728509bceSStefano Zampini 
213828509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
213928509bceSStefano Zampini 
21400f202f7eSStefano Zampini    Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
214128509bceSStefano Zampini 
214228509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
214328509bceSStefano Zampini 
2144b6fdb6dfSStefano 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.
2145b6fdb6dfSStefano Zampini 
21460f202f7eSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace()
214728509bceSStefano Zampini 
21480f202f7eSStefano 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()
21490f202f7eSStefano Zampini    Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS()
215028509bceSStefano Zampini 
21510f202f7eSStefano Zampini    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
215228509bceSStefano Zampini 
21530f202f7eSStefano 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.
21540f202f7eSStefano Zampini    User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
215528509bceSStefano Zampini 
21560f202f7eSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
215728509bceSStefano Zampini 
21580f202f7eSStefano 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.
215928509bceSStefano Zampini 
21600f202f7eSStefano 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.
21610f202f7eSStefano Zampini    Deluxe scaling is not supported yet for FETI-DP.
21620f202f7eSStefano Zampini 
21630f202f7eSStefano Zampini    Options Database Keys (some of them, run with -h for a complete list):
21640f202f7eSStefano Zampini 
21650f202f7eSStefano Zampini .    -pc_bddc_use_vertices <true> - use or not vertices in primal space
21660f202f7eSStefano Zampini .    -pc_bddc_use_edges <true> - use or not edges in primal space
21670f202f7eSStefano Zampini .    -pc_bddc_use_faces <false> - use or not faces in primal space
21680f202f7eSStefano Zampini .    -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
21690f202f7eSStefano Zampini .    -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
21700f202f7eSStefano Zampini .    -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
21710f202f7eSStefano Zampini .    -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
217228509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
21730f202f7eSStefano 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)
21740f202f7eSStefano 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)
21750f202f7eSStefano Zampini .    -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
21760f202f7eSStefano 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)
21770f202f7eSStefano 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)
217828509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
217928509bceSStefano Zampini 
218028509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
218128509bceSStefano Zampini .vb
218228509bceSStefano Zampini       -pc_bddc_dirichlet_
218328509bceSStefano Zampini       -pc_bddc_neumann_
218428509bceSStefano Zampini       -pc_bddc_coarse_
218528509bceSStefano Zampini .ve
21860f202f7eSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
218728509bceSStefano Zampini 
21880f202f7eSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
218928509bceSStefano Zampini .vb
2190312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
2191312be037SStefano Zampini       -pc_bddc_neumann_lN_
2192312be037SStefano Zampini       -pc_bddc_coarse_lN_
219328509bceSStefano Zampini .ve
21940f202f7eSStefano Zampini    Note that level number ranges from the finest (0) to the coarsest (N).
21950f202f7eSStefano 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.
21960f202f7eSStefano Zampini .vb
21970f202f7eSStefano Zampini      -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
21980f202f7eSStefano Zampini .ve
21990f202f7eSStefano 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
2200da1bb401SStefano Zampini 
2201da1bb401SStefano Zampini    Level: intermediate
2202da1bb401SStefano Zampini 
2203b6fdb6dfSStefano Zampini    Developer notes:
2204da1bb401SStefano Zampini 
2205da1bb401SStefano Zampini    Contributed by Stefano Zampini
2206da1bb401SStefano Zampini 
2207da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
2208da1bb401SStefano Zampini M*/
2209b2573a8aSBarry Smith 
2210da1bb401SStefano Zampini #undef __FUNCT__
2211da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
22128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2213da1bb401SStefano Zampini {
2214da1bb401SStefano Zampini   PetscErrorCode      ierr;
2215da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
2216da1bb401SStefano Zampini 
2217da1bb401SStefano Zampini   PetscFunctionBegin;
2218da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
2219b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
2220da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
2221da1bb401SStefano Zampini 
2222da1bb401SStefano Zampini   /* create PCIS data structure */
2223da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
2224da1bb401SStefano Zampini 
222547d04d0dSStefano Zampini   /* BDDC customization */
222608a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
222747d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
222847d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
222947d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
223047d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
223147d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
223247d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
2233fa434743SStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE;
2234fa434743SStefano Zampini   pcbddc->use_qr_single       = PETSC_FALSE;
22353301b35fSStefano Zampini   pcbddc->symmetric_primal    = PETSC_TRUE;
223606a4e24aSStefano Zampini   pcbddc->benign_saddle_point = PETSC_FALSE;
223747d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
2238b9d89cd5SStefano Zampini   /* private */
2239727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
22400e6343abSStefano Zampini   pcbddc->local_primal_size_cc       = 0;
22410e6343abSStefano Zampini   pcbddc->local_primal_ref_node      = 0;
22420e6343abSStefano Zampini   pcbddc->local_primal_ref_mult      = 0;
2243e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
2244727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
2245fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
224668457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
2247f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
2248727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
2249f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
2250f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
2251f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
2252674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
22530bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
22543972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
2255534831adSStefano Zampini   pcbddc->original_rhs               = 0;
2256534831adSStefano Zampini   pcbddc->local_mat                  = 0;
2257534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
2258b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
2259906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
2260da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
2261da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
2262da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
2263da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
226415aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
226515aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
2266da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
2267da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
2268da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
2269da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
2270da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
2271da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
2272da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
2273da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
2274da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
2275da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
2276785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
2277785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
2278785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
227960d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
228060d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
228163602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
2282da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
228363602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
2284da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
228585c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
228647d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
228747d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
2288b0c7d250SStefano Zampini   pcbddc->coarse_adj_red             = 0;
22894fad6a16SStefano Zampini   pcbddc->current_level              = 0;
22902b510759SStefano Zampini   pcbddc->max_levels                 = 0;
2291323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
229274e2c79eSStefano Zampini   pcbddc->redistribute_coarse        = 0;
2293f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
2294323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
229581d14e9dSStefano Zampini 
229681d14e9dSStefano Zampini   /* benign subspace trick */
2297f5fd1fbfSStefano Zampini   pcbddc->zerodiag                   = 0;
229881d14e9dSStefano Zampini   pcbddc->B0_ncol                    = 0;
229981d14e9dSStefano Zampini   pcbddc->B0_cols                    = NULL;
230081d14e9dSStefano Zampini   pcbddc->B0_vals                    = NULL;
230181d14e9dSStefano Zampini   pcbddc->benign_change              = 0;
230281d14e9dSStefano Zampini 
2303674ae819SStefano Zampini   /* create local graph structure */
2304674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2305674ae819SStefano Zampini 
2306674ae819SStefano Zampini   /* scaling */
2307674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
230834a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2309ac632422SStefano Zampini   pcbddc->faster_deluxe         = PETSC_FALSE;
2310b96c3477SStefano Zampini 
2311b96c3477SStefano Zampini   /* create sub schurs structure */
2312b96c3477SStefano Zampini   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2313b96c3477SStefano Zampini   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2314b96c3477SStefano Zampini   pcbddc->sub_schurs_layers      = -1;
2315b96c3477SStefano Zampini   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2316b96c3477SStefano Zampini 
2317b96c3477SStefano Zampini   pcbddc->computed_rowadj = PETSC_FALSE;
2318da1bb401SStefano Zampini 
2319b7eb3628SStefano Zampini   /* adaptivity */
2320f6f667cfSStefano Zampini   pcbddc->adaptive_threshold      = 0.0;
232108122e43SStefano Zampini   pcbddc->adaptive_nmax           = 0;
2322f6f667cfSStefano Zampini   pcbddc->adaptive_nmin           = 0;
2323b7eb3628SStefano Zampini 
2324da1bb401SStefano Zampini   /* function pointers */
2325da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
232693bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2327da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
2328da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
2329da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2330da1bb401SStefano Zampini   pc->ops->view                = 0;
2331da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
2332da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
2333da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
2334534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
2335534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
2336da1bb401SStefano Zampini 
2337da1bb401SStefano Zampini   /* composing function */
2338906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2339674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2340bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
23412b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2342b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
23432b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2344bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2345bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
234682ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2347bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
234882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2349bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
235082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2351bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
235282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2353bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
235463602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2355bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2356bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2357bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2358bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2359da1bb401SStefano Zampini   PetscFunctionReturn(0);
2360da1bb401SStefano Zampini }
23613425bc38SStefano Zampini 
2362