xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision b7eb36285fd32db6a2205f2794656e34b279d97c)
153cdbc3dSStefano Zampini /* TODOLIST
2eb97c9d2SStefano Zampini 
3eb97c9d2SStefano Zampini    ConstraintsSetup
4eb97c9d2SStefano Zampini    - tolerances for constraints as an option (take care of single precision!)
5eb97c9d2SStefano Zampini 
6eb97c9d2SStefano Zampini    Solvers
7a0d3c3abSStefano Zampini    - Add support for cholesky for coarse solver (similar to local solvers)
8eb97c9d2SStefano Zampini    - Propagate ksp prefixes for solvers to mat objects?
9eb97c9d2SStefano Zampini    - Propagate nearnullspace info among levels
10eb97c9d2SStefano Zampini 
11eb97c9d2SStefano Zampini    User interface
12b9b85e73SStefano Zampini    - ** DofSplitting and DM attached to pc?
13eb97c9d2SStefano Zampini 
14eb97c9d2SStefano Zampini    Debugging output
15b9b85e73SStefano Zampini    - * Better management of verbosity levels of debugging output
16eb97c9d2SStefano Zampini 
17eb97c9d2SStefano Zampini    Build
18eb97c9d2SStefano Zampini    - make runexe59
19eb97c9d2SStefano Zampini 
20eb97c9d2SStefano Zampini    Extra
21b9b85e73SStefano Zampini    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
22c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
23eb97c9d2SStefano Zampini    - BDDC with MG framework?
24eb97c9d2SStefano Zampini 
25eb97c9d2SStefano Zampini    FETIDP
26eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
27eb97c9d2SStefano Zampini 
28eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
29eb97c9d2SStefano Zampini    - Provide general case for subassembling
30b9b85e73SStefano Zampini    - *** Preallocation routines in MatISGetMPIAXAIJ
31eb97c9d2SStefano Zampini 
3253cdbc3dSStefano Zampini */
330c7d97c5SJed Brown 
3453cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
350c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
360c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
3753cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
3853cdbc3dSStefano Zampini 
39ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
40ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
413b03a366Sstefano_zampini #include <petscblaslapack.h>
42674ae819SStefano Zampini 
430c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
440c7d97c5SJed Brown #undef __FUNCT__
450c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
460c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
470c7d97c5SJed Brown {
480c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
490c7d97c5SJed Brown   PetscErrorCode ierr;
500c7d97c5SJed Brown 
510c7d97c5SJed Brown   PetscFunctionBegin;
520c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
538eeda7d8SStefano Zampini   /* Verbose debugging */
548eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
558eeda7d8SStefano Zampini   /* Primal space cumstomization */
5608a5cf49SStefano 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);
578eeda7d8SStefano 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);
588eeda7d8SStefano 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);
598eeda7d8SStefano 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);
606661aa1dSStefano 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);
61fa434743SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is always used for multiple constraints on the same cc)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr);
628eeda7d8SStefano Zampini   /* Change of basis */
63b9b85e73SStefano 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);
64b9b85e73SStefano 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);
65674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
66674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
67674ae819SStefano Zampini   }
688eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
698eeda7d8SStefano 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);
7074e2c79eSStefano 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);
710298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr);
722b510759SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr);
73323d395dSStefano 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);
74674ae819SStefano 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);
75b96c3477SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_schur_threshold","Schur principal minors smaller than this value are explicilty computed (-1 computes all)","none",pcbddc->sub_schurs_threshold,&pcbddc->sub_schurs_threshold,NULL);CHKERRQ(ierr);
76b96c3477SStefano 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);
77b96c3477SStefano 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);
78b96c3477SStefano 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);
79*b7eb3628SStefano Zampini   ierr = PetscOptionsScalar("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr);
800c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
810c7d97c5SJed Brown   PetscFunctionReturn(0);
820c7d97c5SJed Brown }
830c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
84674ae819SStefano Zampini #undef __FUNCT__
85906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
86906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
87b9b85e73SStefano Zampini {
88b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
89b9b85e73SStefano Zampini   PetscErrorCode ierr;
90b9b85e73SStefano Zampini 
91b9b85e73SStefano Zampini   PetscFunctionBegin;
92b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
93b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
94b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
95b9b85e73SStefano Zampini   PetscFunctionReturn(0);
96b9b85e73SStefano Zampini }
97b9b85e73SStefano Zampini #undef __FUNCT__
98906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
99b9b85e73SStefano Zampini /*@
100906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
101b9b85e73SStefano Zampini 
102b9b85e73SStefano Zampini    Collective on PC
103b9b85e73SStefano Zampini 
104b9b85e73SStefano Zampini    Input Parameters:
105b9b85e73SStefano Zampini +  pc - the preconditioning context
106906d46d4SStefano Zampini -  change - the change of basis matrix
107b9b85e73SStefano Zampini 
108b9b85e73SStefano Zampini    Level: intermediate
109b9b85e73SStefano Zampini 
110b9b85e73SStefano Zampini    Notes:
111b9b85e73SStefano Zampini 
112b9b85e73SStefano Zampini .seealso: PCBDDC
113b9b85e73SStefano Zampini @*/
114906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
115b9b85e73SStefano Zampini {
116b9b85e73SStefano Zampini   PetscErrorCode ierr;
117b9b85e73SStefano Zampini 
118b9b85e73SStefano Zampini   PetscFunctionBegin;
119b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
120b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
121906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
122906d46d4SStefano Zampini   if (pc->mat) {
123906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
124906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
125906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
126906d46d4SStefano Zampini     if (rows_c != rows) {
127906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
128906d46d4SStefano Zampini     }
129906d46d4SStefano Zampini     if (cols_c != cols) {
130906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
131906d46d4SStefano Zampini     }
132906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
133906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
134906d46d4SStefano Zampini     if (rows_c != rows) {
135906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
136906d46d4SStefano Zampini     }
137906d46d4SStefano Zampini     if (cols_c != cols) {
138906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
139906d46d4SStefano Zampini     }
140906d46d4SStefano Zampini   }
141906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
142b9b85e73SStefano Zampini   PetscFunctionReturn(0);
143b9b85e73SStefano Zampini }
144b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
145b9b85e73SStefano Zampini #undef __FUNCT__
146674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
147674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
148674ae819SStefano Zampini {
149674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
150674ae819SStefano Zampini   PetscErrorCode ierr;
1511e6b0712SBarry Smith 
152674ae819SStefano Zampini   PetscFunctionBegin;
153674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
154674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
155674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
156674ae819SStefano Zampini   PetscFunctionReturn(0);
157674ae819SStefano Zampini }
158674ae819SStefano Zampini #undef __FUNCT__
159674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
160674ae819SStefano Zampini /*@
16128509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
162674ae819SStefano Zampini 
163674ae819SStefano Zampini    Not collective
164674ae819SStefano Zampini 
165674ae819SStefano Zampini    Input Parameters:
166674ae819SStefano Zampini +  pc - the preconditioning context
16728509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
168674ae819SStefano Zampini 
169674ae819SStefano Zampini    Level: intermediate
170674ae819SStefano Zampini 
171674ae819SStefano Zampini    Notes:
172674ae819SStefano Zampini 
173674ae819SStefano Zampini .seealso: PCBDDC
174674ae819SStefano Zampini @*/
175674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
176674ae819SStefano Zampini {
177674ae819SStefano Zampini   PetscErrorCode ierr;
178674ae819SStefano Zampini 
179674ae819SStefano Zampini   PetscFunctionBegin;
180674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
182674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
183674ae819SStefano Zampini   PetscFunctionReturn(0);
184674ae819SStefano Zampini }
185674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1860c7d97c5SJed Brown #undef __FUNCT__
1874fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1884fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1894fad6a16SStefano Zampini {
1904fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1914fad6a16SStefano Zampini 
1924fad6a16SStefano Zampini   PetscFunctionBegin;
1934fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1944fad6a16SStefano Zampini   PetscFunctionReturn(0);
1954fad6a16SStefano Zampini }
1961e6b0712SBarry Smith 
1974fad6a16SStefano Zampini #undef __FUNCT__
1984fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1994fad6a16SStefano Zampini /*@
20028509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
2014fad6a16SStefano Zampini 
2024fad6a16SStefano Zampini    Logically collective on PC
2034fad6a16SStefano Zampini 
2044fad6a16SStefano Zampini    Input Parameters:
2054fad6a16SStefano Zampini +  pc - the preconditioning context
20628509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2074fad6a16SStefano Zampini 
20828509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2094fad6a16SStefano Zampini 
2104fad6a16SStefano Zampini    Level: intermediate
2114fad6a16SStefano Zampini 
2124fad6a16SStefano Zampini    Notes:
2134fad6a16SStefano Zampini 
2144fad6a16SStefano Zampini .seealso: PCBDDC
2154fad6a16SStefano Zampini @*/
2164fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2174fad6a16SStefano Zampini {
2184fad6a16SStefano Zampini   PetscErrorCode ierr;
2194fad6a16SStefano Zampini 
2204fad6a16SStefano Zampini   PetscFunctionBegin;
2214fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2222b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2234fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2244fad6a16SStefano Zampini   PetscFunctionReturn(0);
2254fad6a16SStefano Zampini }
2262b510759SStefano Zampini 
227b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2282b510759SStefano Zampini #undef __FUNCT__
229b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
230b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
231b8ffe317SStefano Zampini {
232b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
233b8ffe317SStefano Zampini 
234b8ffe317SStefano Zampini   PetscFunctionBegin;
23585c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
236b8ffe317SStefano Zampini   PetscFunctionReturn(0);
237b8ffe317SStefano Zampini }
238b8ffe317SStefano Zampini 
239b8ffe317SStefano Zampini #undef __FUNCT__
240b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
241b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2422b510759SStefano Zampini {
2432b510759SStefano Zampini   PetscErrorCode ierr;
2442b510759SStefano Zampini 
2452b510759SStefano Zampini   PetscFunctionBegin;
2462b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
247b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
248b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2492b510759SStefano Zampini   PetscFunctionReturn(0);
2502b510759SStefano Zampini }
2511e6b0712SBarry Smith 
2524fad6a16SStefano Zampini #undef __FUNCT__
2532b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2542b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2554fad6a16SStefano Zampini {
2564fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2574fad6a16SStefano Zampini 
2584fad6a16SStefano Zampini   PetscFunctionBegin;
2592b510759SStefano Zampini   pcbddc->current_level = level;
2604fad6a16SStefano Zampini   PetscFunctionReturn(0);
2614fad6a16SStefano Zampini }
2621e6b0712SBarry Smith 
2634fad6a16SStefano Zampini #undef __FUNCT__
264b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
265b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
266b8ffe317SStefano Zampini {
267b8ffe317SStefano Zampini   PetscErrorCode ierr;
268b8ffe317SStefano Zampini 
269b8ffe317SStefano Zampini   PetscFunctionBegin;
270b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
271b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
272b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
273b8ffe317SStefano Zampini   PetscFunctionReturn(0);
274b8ffe317SStefano Zampini }
275b8ffe317SStefano Zampini 
276b8ffe317SStefano Zampini #undef __FUNCT__
2772b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2782b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2792b510759SStefano Zampini {
2802b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2812b510759SStefano Zampini 
2822b510759SStefano Zampini   PetscFunctionBegin;
2832b510759SStefano Zampini   pcbddc->max_levels = levels;
2842b510759SStefano Zampini   PetscFunctionReturn(0);
2852b510759SStefano Zampini }
2862b510759SStefano Zampini 
2872b510759SStefano Zampini #undef __FUNCT__
2882b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2894fad6a16SStefano Zampini /*@
29028509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2914fad6a16SStefano Zampini 
2924fad6a16SStefano Zampini    Logically collective on PC
2934fad6a16SStefano Zampini 
2944fad6a16SStefano Zampini    Input Parameters:
2954fad6a16SStefano Zampini +  pc - the preconditioning context
29628509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2974fad6a16SStefano Zampini 
29828509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2994fad6a16SStefano Zampini 
3004fad6a16SStefano Zampini    Level: intermediate
3014fad6a16SStefano Zampini 
3024fad6a16SStefano Zampini    Notes:
3034fad6a16SStefano Zampini 
3044fad6a16SStefano Zampini .seealso: PCBDDC
3054fad6a16SStefano Zampini @*/
3062b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3074fad6a16SStefano Zampini {
3084fad6a16SStefano Zampini   PetscErrorCode ierr;
3094fad6a16SStefano Zampini 
3104fad6a16SStefano Zampini   PetscFunctionBegin;
3114fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3122b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
313312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3142b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3154fad6a16SStefano Zampini   PetscFunctionReturn(0);
3164fad6a16SStefano Zampini }
3174fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3181e6b0712SBarry Smith 
3194fad6a16SStefano Zampini #undef __FUNCT__
3200bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3210bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3220bdf917eSStefano Zampini {
3230bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3240bdf917eSStefano Zampini   PetscErrorCode ierr;
3250bdf917eSStefano Zampini 
3260bdf917eSStefano Zampini   PetscFunctionBegin;
3270bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3280bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3290bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3300bdf917eSStefano Zampini   PetscFunctionReturn(0);
3310bdf917eSStefano Zampini }
3321e6b0712SBarry Smith 
3330bdf917eSStefano Zampini #undef __FUNCT__
3340bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3350bdf917eSStefano Zampini /*@
33628509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3370bdf917eSStefano Zampini 
3380bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3390bdf917eSStefano Zampini 
3400bdf917eSStefano Zampini    Input Parameters:
3410bdf917eSStefano Zampini +  pc - the preconditioning context
34228509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3430bdf917eSStefano Zampini 
3440bdf917eSStefano Zampini    Level: intermediate
3450bdf917eSStefano Zampini 
3460bdf917eSStefano Zampini    Notes:
3470bdf917eSStefano Zampini 
3480bdf917eSStefano Zampini .seealso: PCBDDC
3490bdf917eSStefano Zampini @*/
3500bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3510bdf917eSStefano Zampini {
3520bdf917eSStefano Zampini   PetscErrorCode ierr;
3530bdf917eSStefano Zampini 
3540bdf917eSStefano Zampini   PetscFunctionBegin;
3550bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
356674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
3572b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
3580bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
3590bdf917eSStefano Zampini   PetscFunctionReturn(0);
3600bdf917eSStefano Zampini }
3610bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3621e6b0712SBarry Smith 
3630bdf917eSStefano Zampini #undef __FUNCT__
3643b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3653b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3663b03a366Sstefano_zampini {
3673b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3683b03a366Sstefano_zampini   PetscErrorCode ierr;
3693b03a366Sstefano_zampini 
3703b03a366Sstefano_zampini   PetscFunctionBegin;
371785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
372785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3733b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
37436e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
37536e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
376fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3773b03a366Sstefano_zampini   PetscFunctionReturn(0);
3783b03a366Sstefano_zampini }
3791e6b0712SBarry Smith 
3803b03a366Sstefano_zampini #undef __FUNCT__
3813b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3823b03a366Sstefano_zampini /*@
38328509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3843b03a366Sstefano_zampini 
385785d1243SStefano Zampini    Collective
3863b03a366Sstefano_zampini 
3873b03a366Sstefano_zampini    Input Parameters:
3883b03a366Sstefano_zampini +  pc - the preconditioning context
389785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3903b03a366Sstefano_zampini 
3913b03a366Sstefano_zampini    Level: intermediate
3923b03a366Sstefano_zampini 
393785d1243SStefano Zampini    Notes: Any process can list any global node
3943b03a366Sstefano_zampini 
3953b03a366Sstefano_zampini .seealso: PCBDDC
3963b03a366Sstefano_zampini @*/
3973b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3983b03a366Sstefano_zampini {
3993b03a366Sstefano_zampini   PetscErrorCode ierr;
4003b03a366Sstefano_zampini 
4013b03a366Sstefano_zampini   PetscFunctionBegin;
4023b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
403674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
404785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4053b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4063b03a366Sstefano_zampini   PetscFunctionReturn(0);
4073b03a366Sstefano_zampini }
4083b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4091e6b0712SBarry Smith 
4103b03a366Sstefano_zampini #undef __FUNCT__
41182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
41282ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4130c7d97c5SJed Brown {
4140c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4150c7d97c5SJed Brown   PetscErrorCode ierr;
4160c7d97c5SJed Brown 
4170c7d97c5SJed Brown   PetscFunctionBegin;
418785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
419785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4200c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4210c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
422785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4237d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4240c7d97c5SJed Brown   PetscFunctionReturn(0);
4250c7d97c5SJed Brown }
4260c7d97c5SJed Brown 
4270c7d97c5SJed Brown #undef __FUNCT__
42882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4299c0446d6SStefano Zampini /*@
43082ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4319c0446d6SStefano Zampini 
432785d1243SStefano Zampini    Collective
43353cdbc3dSStefano Zampini 
43453cdbc3dSStefano Zampini    Input Parameters:
43553cdbc3dSStefano Zampini +  pc - the preconditioning context
43682ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
43753cdbc3dSStefano Zampini 
43853cdbc3dSStefano Zampini    Level: intermediate
43953cdbc3dSStefano Zampini 
4409c0446d6SStefano Zampini    Notes:
44153cdbc3dSStefano Zampini 
44253cdbc3dSStefano Zampini .seealso: PCBDDC
44353cdbc3dSStefano Zampini @*/
44482ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4450c7d97c5SJed Brown {
4460c7d97c5SJed Brown   PetscErrorCode ierr;
4470c7d97c5SJed Brown 
4480c7d97c5SJed Brown   PetscFunctionBegin;
4490c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4500c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
45182ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
45282ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4530c7d97c5SJed Brown   PetscFunctionReturn(0);
4540c7d97c5SJed Brown }
4550c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4560c7d97c5SJed Brown 
4570c7d97c5SJed Brown #undef __FUNCT__
4580c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
45953cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4600c7d97c5SJed Brown {
4610c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
46253cdbc3dSStefano Zampini   PetscErrorCode ierr;
4630c7d97c5SJed Brown 
4640c7d97c5SJed Brown   PetscFunctionBegin;
465785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
466785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
46753cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
46836e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
46936e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
470fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4710c7d97c5SJed Brown   PetscFunctionReturn(0);
4720c7d97c5SJed Brown }
4731e6b0712SBarry Smith 
4740c7d97c5SJed Brown #undef __FUNCT__
4750c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
47657527edcSJed Brown /*@
47728509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
47857527edcSJed Brown 
479785d1243SStefano Zampini    Collective
48057527edcSJed Brown 
48157527edcSJed Brown    Input Parameters:
48257527edcSJed Brown +  pc - the preconditioning context
483785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
48457527edcSJed Brown 
48557527edcSJed Brown    Level: intermediate
48657527edcSJed Brown 
487785d1243SStefano Zampini    Notes: Any process can list any global node
48857527edcSJed Brown 
48957527edcSJed Brown .seealso: PCBDDC
49057527edcSJed Brown @*/
49153cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4920c7d97c5SJed Brown {
4930c7d97c5SJed Brown   PetscErrorCode ierr;
4940c7d97c5SJed Brown 
4950c7d97c5SJed Brown   PetscFunctionBegin;
4960c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
497674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
498785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
49953cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
50053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
50153cdbc3dSStefano Zampini }
50253cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5031e6b0712SBarry Smith 
50453cdbc3dSStefano Zampini #undef __FUNCT__
50582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
50682ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5070c7d97c5SJed Brown {
5080c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5090c7d97c5SJed Brown   PetscErrorCode ierr;
5100c7d97c5SJed Brown 
5110c7d97c5SJed Brown   PetscFunctionBegin;
512785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
513785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5140c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5150c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
516785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5177d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5180c7d97c5SJed Brown   PetscFunctionReturn(0);
5190c7d97c5SJed Brown }
5200c7d97c5SJed Brown 
5210c7d97c5SJed Brown #undef __FUNCT__
52282ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5230c7d97c5SJed Brown /*@
52482ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5250c7d97c5SJed Brown 
526785d1243SStefano Zampini    Collective
5270c7d97c5SJed Brown 
5280c7d97c5SJed Brown    Input Parameters:
5290c7d97c5SJed Brown +  pc - the preconditioning context
53082ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5310c7d97c5SJed Brown 
5320c7d97c5SJed Brown    Level: intermediate
5330c7d97c5SJed Brown 
5340c7d97c5SJed Brown    Notes:
5350c7d97c5SJed Brown 
5360c7d97c5SJed Brown .seealso: PCBDDC
5370c7d97c5SJed Brown @*/
53882ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5390c7d97c5SJed Brown {
5400c7d97c5SJed Brown   PetscErrorCode ierr;
5410c7d97c5SJed Brown 
5420c7d97c5SJed Brown   PetscFunctionBegin;
5430c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5440c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
54582ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54682ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
54753cdbc3dSStefano Zampini   PetscFunctionReturn(0);
54853cdbc3dSStefano Zampini }
54953cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
55053cdbc3dSStefano Zampini 
55153cdbc3dSStefano Zampini #undef __FUNCT__
552da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
553da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
554da1bb401SStefano Zampini {
555da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
556da1bb401SStefano Zampini 
557da1bb401SStefano Zampini   PetscFunctionBegin;
558da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
559da1bb401SStefano Zampini   PetscFunctionReturn(0);
560da1bb401SStefano Zampini }
5611e6b0712SBarry Smith 
562da1bb401SStefano Zampini #undef __FUNCT__
563da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
564da1bb401SStefano Zampini /*@
565785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
566da1bb401SStefano Zampini 
567785d1243SStefano Zampini    Collective
568785d1243SStefano Zampini 
569785d1243SStefano Zampini    Input Parameters:
570785d1243SStefano Zampini .  pc - the preconditioning context
571785d1243SStefano Zampini 
572785d1243SStefano Zampini    Output Parameters:
573785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
574785d1243SStefano Zampini 
575785d1243SStefano Zampini    Level: intermediate
576785d1243SStefano Zampini 
577785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
578785d1243SStefano Zampini 
579785d1243SStefano Zampini .seealso: PCBDDC
580785d1243SStefano Zampini @*/
581785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
582785d1243SStefano Zampini {
583785d1243SStefano Zampini   PetscErrorCode ierr;
584785d1243SStefano Zampini 
585785d1243SStefano Zampini   PetscFunctionBegin;
586785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
587785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
588785d1243SStefano Zampini   PetscFunctionReturn(0);
589785d1243SStefano Zampini }
590785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
591785d1243SStefano Zampini 
592785d1243SStefano Zampini #undef __FUNCT__
593785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
594785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
595785d1243SStefano Zampini {
596785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
597785d1243SStefano Zampini 
598785d1243SStefano Zampini   PetscFunctionBegin;
599785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
600785d1243SStefano Zampini   PetscFunctionReturn(0);
601785d1243SStefano Zampini }
602785d1243SStefano Zampini 
603785d1243SStefano Zampini #undef __FUNCT__
60482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
605da1bb401SStefano Zampini /*@
60682ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
607da1bb401SStefano Zampini 
608785d1243SStefano Zampini    Collective
609da1bb401SStefano Zampini 
610da1bb401SStefano Zampini    Input Parameters:
61128509bceSStefano Zampini .  pc - the preconditioning context
612da1bb401SStefano Zampini 
613da1bb401SStefano Zampini    Output Parameters:
61428509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
615da1bb401SStefano Zampini 
616da1bb401SStefano Zampini    Level: intermediate
617da1bb401SStefano Zampini 
618da1bb401SStefano Zampini    Notes:
619da1bb401SStefano Zampini 
620da1bb401SStefano Zampini .seealso: PCBDDC
621da1bb401SStefano Zampini @*/
62282ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
623da1bb401SStefano Zampini {
624da1bb401SStefano Zampini   PetscErrorCode ierr;
625da1bb401SStefano Zampini 
626da1bb401SStefano Zampini   PetscFunctionBegin;
627da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
62882ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
629da1bb401SStefano Zampini   PetscFunctionReturn(0);
630da1bb401SStefano Zampini }
631da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6321e6b0712SBarry Smith 
633da1bb401SStefano Zampini #undef __FUNCT__
63453cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
63553cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
63653cdbc3dSStefano Zampini {
63753cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
63853cdbc3dSStefano Zampini 
63953cdbc3dSStefano Zampini   PetscFunctionBegin;
64053cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
64153cdbc3dSStefano Zampini   PetscFunctionReturn(0);
64253cdbc3dSStefano Zampini }
6431e6b0712SBarry Smith 
64453cdbc3dSStefano Zampini #undef __FUNCT__
64553cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
64653cdbc3dSStefano Zampini /*@
647785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
64853cdbc3dSStefano Zampini 
649785d1243SStefano Zampini    Collective
650785d1243SStefano Zampini 
651785d1243SStefano Zampini    Input Parameters:
652785d1243SStefano Zampini .  pc - the preconditioning context
653785d1243SStefano Zampini 
654785d1243SStefano Zampini    Output Parameters:
655785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
656785d1243SStefano Zampini 
657785d1243SStefano Zampini    Level: intermediate
658785d1243SStefano Zampini 
659785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
660785d1243SStefano Zampini 
661785d1243SStefano Zampini .seealso: PCBDDC
662785d1243SStefano Zampini @*/
663785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
664785d1243SStefano Zampini {
665785d1243SStefano Zampini   PetscErrorCode ierr;
666785d1243SStefano Zampini 
667785d1243SStefano Zampini   PetscFunctionBegin;
668785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
669785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
670785d1243SStefano Zampini   PetscFunctionReturn(0);
671785d1243SStefano Zampini }
672785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
673785d1243SStefano Zampini 
674785d1243SStefano Zampini #undef __FUNCT__
675785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
676785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
677785d1243SStefano Zampini {
678785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
679785d1243SStefano Zampini 
680785d1243SStefano Zampini   PetscFunctionBegin;
681785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
682785d1243SStefano Zampini   PetscFunctionReturn(0);
683785d1243SStefano Zampini }
684785d1243SStefano Zampini 
685785d1243SStefano Zampini #undef __FUNCT__
68682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
68753cdbc3dSStefano Zampini /*@
68882ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
68953cdbc3dSStefano Zampini 
690785d1243SStefano Zampini    Collective
69153cdbc3dSStefano Zampini 
69253cdbc3dSStefano Zampini    Input Parameters:
69328509bceSStefano Zampini .  pc - the preconditioning context
69453cdbc3dSStefano Zampini 
69553cdbc3dSStefano Zampini    Output Parameters:
69628509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
69753cdbc3dSStefano Zampini 
69853cdbc3dSStefano Zampini    Level: intermediate
69953cdbc3dSStefano Zampini 
70053cdbc3dSStefano Zampini    Notes:
70153cdbc3dSStefano Zampini 
70253cdbc3dSStefano Zampini .seealso: PCBDDC
70353cdbc3dSStefano Zampini @*/
70482ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
70553cdbc3dSStefano Zampini {
70653cdbc3dSStefano Zampini   PetscErrorCode ierr;
70753cdbc3dSStefano Zampini 
70853cdbc3dSStefano Zampini   PetscFunctionBegin;
70953cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
71082ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7110c7d97c5SJed Brown   PetscFunctionReturn(0);
7120c7d97c5SJed Brown }
71336e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7141e6b0712SBarry Smith 
71536e030ebSStefano Zampini #undef __FUNCT__
716da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7171a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
71836e030ebSStefano Zampini {
71936e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
720da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
721da1bb401SStefano Zampini   PetscErrorCode ierr;
72236e030ebSStefano Zampini 
72336e030ebSStefano Zampini   PetscFunctionBegin;
724674ae819SStefano Zampini   /* free old CSR */
725674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
726fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
727674ae819SStefano Zampini   /* get CSR into graph structure */
728da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
729854ce69bSBarry Smith     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
730785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
731674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
732674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
733da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7341a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7351a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
736674ae819SStefano Zampini   }
737575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
73836e030ebSStefano Zampini   PetscFunctionReturn(0);
73936e030ebSStefano Zampini }
7401e6b0712SBarry Smith 
74136e030ebSStefano Zampini #undef __FUNCT__
742da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
74336e030ebSStefano Zampini /*@
74428509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
74536e030ebSStefano Zampini 
74636e030ebSStefano Zampini    Not collective
74736e030ebSStefano Zampini 
74836e030ebSStefano Zampini    Input Parameters:
74936e030ebSStefano Zampini +  pc - the preconditioning context
75028509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
75128509bceSStefano Zampini .  xadj, adjncy - the CSR graph
75228509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
75336e030ebSStefano Zampini 
75436e030ebSStefano Zampini    Level: intermediate
75536e030ebSStefano Zampini 
75636e030ebSStefano Zampini    Notes:
75736e030ebSStefano Zampini 
75828509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
75936e030ebSStefano Zampini @*/
7601a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
76136e030ebSStefano Zampini {
762575ad6abSStefano Zampini   void (*f)(void) = 0;
76336e030ebSStefano Zampini   PetscErrorCode ierr;
76436e030ebSStefano Zampini 
76536e030ebSStefano Zampini   PetscFunctionBegin;
76636e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
767674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
768d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
769674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
770674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
771da1bb401SStefano Zampini   }
77236e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
773575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
774575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
775575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
776575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
777575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
77836e030ebSStefano Zampini   }
77936e030ebSStefano Zampini   PetscFunctionReturn(0);
78036e030ebSStefano Zampini }
7819c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7821e6b0712SBarry Smith 
7839c0446d6SStefano Zampini #undef __FUNCT__
78463602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
78563602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
78663602bcaSStefano Zampini {
78763602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
78863602bcaSStefano Zampini   PetscInt i;
78963602bcaSStefano Zampini   PetscErrorCode ierr;
79063602bcaSStefano Zampini 
79163602bcaSStefano Zampini   PetscFunctionBegin;
79263602bcaSStefano Zampini   /* Destroy ISes if they were already set */
79363602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
79463602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
79563602bcaSStefano Zampini   }
79663602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
79763602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
79863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
79963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
80063602bcaSStefano Zampini   }
80163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
80263602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
80363602bcaSStefano Zampini   /* allocate space then set */
804d02579f5SStefano Zampini   if (n_is) {
805d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
806d02579f5SStefano Zampini   }
80763602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
80863602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
80963602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
81063602bcaSStefano Zampini   }
81163602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
81263602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8131a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
81463602bcaSStefano Zampini   PetscFunctionReturn(0);
81563602bcaSStefano Zampini }
81663602bcaSStefano Zampini 
81763602bcaSStefano Zampini #undef __FUNCT__
81863602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
81963602bcaSStefano Zampini /*@
82063602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
82163602bcaSStefano Zampini 
82263602bcaSStefano Zampini    Collective
82363602bcaSStefano Zampini 
82463602bcaSStefano Zampini    Input Parameters:
82563602bcaSStefano Zampini +  pc - the preconditioning context
82663602bcaSStefano Zampini -  n_is - number of index sets defining the fields
82763602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in local ordering
82863602bcaSStefano Zampini 
82963602bcaSStefano Zampini    Level: intermediate
83063602bcaSStefano Zampini 
83163602bcaSStefano Zampini    Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field.
83263602bcaSStefano Zampini 
83363602bcaSStefano Zampini .seealso: PCBDDC
83463602bcaSStefano Zampini @*/
83563602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
83663602bcaSStefano Zampini {
83763602bcaSStefano Zampini   PetscInt       i;
83863602bcaSStefano Zampini   PetscErrorCode ierr;
83963602bcaSStefano Zampini 
84063602bcaSStefano Zampini   PetscFunctionBegin;
84163602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
84263602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
84363602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
84463602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
84563602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
84663602bcaSStefano Zampini   }
847e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
84863602bcaSStefano Zampini   PetscFunctionReturn(0);
84963602bcaSStefano Zampini }
85063602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
85163602bcaSStefano Zampini 
85263602bcaSStefano Zampini #undef __FUNCT__
8539c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8549c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8559c0446d6SStefano Zampini {
8569c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8579c0446d6SStefano Zampini   PetscInt i;
8589c0446d6SStefano Zampini   PetscErrorCode ierr;
8599c0446d6SStefano Zampini 
8609c0446d6SStefano Zampini   PetscFunctionBegin;
861da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8629c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8639c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8649c0446d6SStefano Zampini   }
865d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
86663602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
86763602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
86863602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
86963602bcaSStefano Zampini   }
87063602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
87163602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
872da1bb401SStefano Zampini   /* allocate space then set */
873d02579f5SStefano Zampini   if (n_is) {
874785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
875d02579f5SStefano Zampini   }
8769c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
877da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
878da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8799c0446d6SStefano Zampini   }
8809c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
88163602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8821a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8839c0446d6SStefano Zampini   PetscFunctionReturn(0);
8849c0446d6SStefano Zampini }
8851e6b0712SBarry Smith 
8869c0446d6SStefano Zampini #undef __FUNCT__
8879c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8889c0446d6SStefano Zampini /*@
88963602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8909c0446d6SStefano Zampini 
89163602bcaSStefano Zampini    Collective
8929c0446d6SStefano Zampini 
8939c0446d6SStefano Zampini    Input Parameters:
8949c0446d6SStefano Zampini +  pc - the preconditioning context
89528509bceSStefano Zampini -  n_is - number of index sets defining the fields
89663602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in global ordering
8979c0446d6SStefano Zampini 
8989c0446d6SStefano Zampini    Level: intermediate
8999c0446d6SStefano Zampini 
90063602bcaSStefano Zampini    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
9019c0446d6SStefano Zampini 
9029c0446d6SStefano Zampini .seealso: PCBDDC
9039c0446d6SStefano Zampini @*/
9049c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9059c0446d6SStefano Zampini {
9062b510759SStefano Zampini   PetscInt       i;
9079c0446d6SStefano Zampini   PetscErrorCode ierr;
9089c0446d6SStefano Zampini 
9099c0446d6SStefano Zampini   PetscFunctionBegin;
9109c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
91163602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9122b510759SStefano Zampini   for (i=0;i<n_is;i++) {
91363602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
91463602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9152b510759SStefano Zampini   }
9169c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9179c0446d6SStefano Zampini   PetscFunctionReturn(0);
9189c0446d6SStefano Zampini }
919906d46d4SStefano Zampini 
920da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
921534831adSStefano Zampini #undef __FUNCT__
922534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
923534831adSStefano Zampini /* -------------------------------------------------------------------------- */
924534831adSStefano Zampini /*
925534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
926534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9279c0446d6SStefano Zampini 
928534831adSStefano Zampini    Input Parameter:
929534831adSStefano Zampini +  pc - the preconditioner contex
930534831adSStefano Zampini 
931534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
932534831adSStefano Zampini 
933534831adSStefano Zampini    Notes:
934534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
935534831adSStefano Zampini    the user, but instead is called by KSPSolve().
936534831adSStefano Zampini */
937534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
938534831adSStefano Zampini {
939534831adSStefano Zampini   PetscErrorCode ierr;
940534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
941534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9423972b0daSStefano Zampini   IS             dirIS;
9433972b0daSStefano Zampini   Vec            used_vec;
9448d00608fSStefano Zampini   PetscBool      copy_rhs = PETSC_TRUE;
945534831adSStefano Zampini 
946534831adSStefano Zampini   PetscFunctionBegin;
94785c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
94885c4d303SStefano Zampini   if (ksp) {
94985c4d303SStefano Zampini     PetscBool iscg;
95085c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
95185c4d303SStefano Zampini     if (!iscg) {
95285c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
95385c4d303SStefano Zampini     }
95485c4d303SStefano Zampini   }
95585c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
95662a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
95762a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
95862a6ff1dSStefano Zampini   }
95962a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
96062a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
96162a6ff1dSStefano Zampini   }
9628d00608fSStefano Zampini 
9633972b0daSStefano Zampini   if (x) {
9643972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9653972b0daSStefano Zampini     used_vec = x;
9668d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
9673972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9683972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9693972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9703972b0daSStefano Zampini   }
9718efcfb23SStefano Zampini 
9728efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
9733972b0daSStefano Zampini   if (ksp) {
9748efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
9758efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
9763972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9773972b0daSStefano Zampini     }
9783972b0daSStefano Zampini   }
9793308cffdSStefano Zampini 
9808d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
9813972b0daSStefano Zampini 
9823972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
983785d1243SStefano Zampini   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
98482ba6b80SStefano Zampini   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
98582ba6b80SStefano Zampini   if (rhs && dirIS) {
986906d46d4SStefano Zampini     Mat_IS      *matis = (Mat_IS*)pc->pmat->data;
987785d1243SStefano Zampini     PetscInt    dirsize,i,*is_indices;
988785d1243SStefano Zampini     PetscScalar *array_x,*array_diagonal;
989785d1243SStefano Zampini 
9903972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
9913972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
9923972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9933972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9943972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9953972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
99682ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
9973972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9983972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9993972b0daSStefano Zampini     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10002fa5cd67SKarl Rupp     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
10013972b0daSStefano Zampini     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10023972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10033972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10043972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10053972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10068d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10078efcfb23SStefano Zampini   }
1008b76ba322SStefano Zampini 
10098efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
10108d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
10118d00608fSStefano Zampini     /* store the original rhs */
10128d00608fSStefano Zampini     if (copy_rhs) {
10138d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10148d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10158d00608fSStefano Zampini     }
10168d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10173972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10183972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10193972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10208efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10213308cffdSStefano Zampini   }
10228efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1023b76ba322SStefano Zampini 
1024b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
10253972b0daSStefano Zampini   if (x) {
10268efcfb23SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
102785c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
1028b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1029b76ba322SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1030b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
10318efcfb23SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10328efcfb23SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10338efcfb23SStefano Zampini       if (ksp && !pcbddc->ksp_guess_nonzero) {
1034b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1035b76ba322SStefano Zampini       }
1036b76ba322SStefano Zampini     }
10373972b0daSStefano Zampini   }
1038b76ba322SStefano Zampini 
1039b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1040906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1041906d46d4SStefano Zampini 
1042906d46d4SStefano Zampini     /* get change ctx */
1043906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1044906d46d4SStefano Zampini 
1045906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1046906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1047906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1048906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1049906d46d4SStefano Zampini 
1050906d46d4SStefano Zampini     /* change iteration matrix */
1051906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1052906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1053906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1054906d46d4SStefano Zampini 
10558d00608fSStefano Zampini     /* store the original rhs */
10568d00608fSStefano Zampini     if (copy_rhs) {
10578d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10588d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10598d00608fSStefano Zampini     }
10608d00608fSStefano Zampini 
1061906d46d4SStefano Zampini     /* change rhs */
1062906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1063906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
10648d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
106592e3dcfbSStefano Zampini   }
106692e3dcfbSStefano Zampini 
106792e3dcfbSStefano Zampini   /* remove nullspace if present */
10688efcfb23SStefano Zampini   if (ksp && x && pcbddc->NullSpace) {
10698efcfb23SStefano Zampini     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
10708d00608fSStefano Zampini     /* store the original rhs */
10718d00608fSStefano Zampini     if (copy_rhs) {
10728d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10738d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10748d00608fSStefano Zampini     }
10758d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
1076d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1077b76ba322SStefano Zampini   }
1078534831adSStefano Zampini   PetscFunctionReturn(0);
1079534831adSStefano Zampini }
1080906d46d4SStefano Zampini 
1081534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1082534831adSStefano Zampini #undef __FUNCT__
1083534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1084534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1085534831adSStefano Zampini /*
1086534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1087534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1088534831adSStefano Zampini 
1089534831adSStefano Zampini    Input Parameter:
1090534831adSStefano Zampini +  pc - the preconditioner contex
1091534831adSStefano Zampini 
1092534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1093534831adSStefano Zampini 
1094534831adSStefano Zampini    Notes:
1095534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
1096534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1097534831adSStefano Zampini */
1098534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1099534831adSStefano Zampini {
1100534831adSStefano Zampini   PetscErrorCode ierr;
1101534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1102534831adSStefano Zampini 
1103534831adSStefano Zampini   PetscFunctionBegin;
1104b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1105906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1106906d46d4SStefano Zampini 
1107906d46d4SStefano Zampini     /* get change ctx */
1108906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1109906d46d4SStefano Zampini 
1110906d46d4SStefano Zampini     /* restore iteration matrix */
1111906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1112906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1113906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1114906d46d4SStefano Zampini 
1115906d46d4SStefano Zampini     /* get solution in original basis */
1116906d46d4SStefano Zampini     if (x) {
1117906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1118906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1119906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
11203425bc38SStefano Zampini     }
1121534831adSStefano Zampini   }
1122906d46d4SStefano Zampini 
11233972b0daSStefano Zampini   /* add solution removed in presolve */
11246bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
11253425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
11263425bc38SStefano Zampini   }
1127906d46d4SStefano Zampini 
1128fb223d50SStefano Zampini   /* restore rhs to its original state */
11298d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
1130fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1131fb223d50SStefano Zampini   }
11328d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
11338efcfb23SStefano Zampini 
11348efcfb23SStefano Zampini   /* restore ksp guess state */
11358efcfb23SStefano Zampini   if (ksp) {
11368efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
11378efcfb23SStefano Zampini   }
1138534831adSStefano Zampini   PetscFunctionReturn(0);
1139534831adSStefano Zampini }
1140534831adSStefano Zampini /* -------------------------------------------------------------------------- */
114153cdbc3dSStefano Zampini #undef __FUNCT__
114253cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
11430c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
11440c7d97c5SJed Brown /*
11450c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
11460c7d97c5SJed Brown                   by setting data structures and options.
11470c7d97c5SJed Brown 
11480c7d97c5SJed Brown    Input Parameter:
114953cdbc3dSStefano Zampini +  pc - the preconditioner context
11500c7d97c5SJed Brown 
11510c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
11520c7d97c5SJed Brown 
11530c7d97c5SJed Brown    Notes:
11540c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
11550c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
11560c7d97c5SJed Brown */
115753cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
11580c7d97c5SJed Brown {
11590c7d97c5SJed Brown   PetscErrorCode ierr;
11600c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
11615e8657edSStefano Zampini   Mat_IS*        matis;
1162*b7eb3628SStefano Zampini   MatNullSpace   nearnullspace,ad_nearnullspace;
1163*b7eb3628SStefano Zampini   PetscBool      computetopography,computesolvers,computesubschurs,adaptive;
11645e8657edSStefano Zampini   PetscBool      new_nearnullspace_provided,ismatis;
11650c7d97c5SJed Brown 
11660c7d97c5SJed Brown   PetscFunctionBegin;
11675e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
11685e8657edSStefano Zampini   if (!ismatis) {
11695e8657edSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
11705e8657edSStefano Zampini   }
11715e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1172f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
11733b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1174fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1175f4ddd8eeSStefano Zampini   /* split work */
1176674ae819SStefano Zampini   if (pc->setupcalled) {
1177d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1178674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1179674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1180674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1181674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1182674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1183674ae819SStefano Zampini     }
1184674ae819SStefano Zampini   } else {
1185674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1186674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1187674ae819SStefano Zampini   }
1188fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1189fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1190fb180af4SStefano Zampini   }
1191*b7eb3628SStefano Zampini   adaptive = (PetscBool)(pcbddc->adaptive_threshold > 0.0);
1192*b7eb3628SStefano Zampini   computesubschurs = pcbddc->use_deluxe_scaling || adaptive;
1193f4ddd8eeSStefano Zampini 
1194f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
119570cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
119670cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
119758a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1198f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1199f4ddd8eeSStefano Zampini     }
120058a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1201f4ddd8eeSStefano Zampini   }
1202f4ddd8eeSStefano Zampini 
1203d65f70fdSStefano Zampini   /*ierr = MatIsSymmetric(pc->pmat,0.,&pcbddc->issym);CHKERRQ(ierr);*/
1204d65f70fdSStefano Zampini   { /* this is a temporary workaround since seqbaij matrices does not have support for symmetry checking */
1205d65f70fdSStefano Zampini     PetscBool setsym;
1206d65f70fdSStefano Zampini     ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&pcbddc->issym);CHKERRQ(ierr);
1207d65f70fdSStefano Zampini     if (!setsym) pcbddc->issym = PETSC_FALSE;
1208d65f70fdSStefano Zampini   }
1209d65f70fdSStefano Zampini 
12105e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
12115e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
12125e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
12135e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
12145e8657edSStefano Zampini   } else {
1215b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
12165e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
12175e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
1218674ae819SStefano Zampini   }
1219f4ddd8eeSStefano Zampini 
12205e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
12215e8657edSStefano Zampini   {
12225e8657edSStefano Zampini     Mat temp_mat;
12235e8657edSStefano Zampini 
12245e8657edSStefano Zampini     temp_mat = matis->A;
12255e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
12265e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
12275e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
12285e8657edSStefano Zampini     matis->A = temp_mat;
12295e8657edSStefano Zampini   }
1230684f6988SStefano Zampini 
1231b96c3477SStefano Zampini   /* Analyze interface and setup sub_schurs data */
1232674ae819SStefano Zampini   if (computetopography) {
1233674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1234674ae819SStefano Zampini   }
1235fb8d54d4SStefano Zampini 
1236b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1237*b7eb3628SStefano Zampini   ad_nearnullspace = NULL;
1238684f6988SStefano Zampini   if (computesolvers) {
1239684f6988SStefano Zampini     ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1240b96c3477SStefano Zampini     if (computesubschurs) {
1241b96c3477SStefano Zampini       if (computetopography) {
1242b96c3477SStefano Zampini         ierr = PCBDDCInitSubSchurs(pc,pcbddc->sub_schurs_rebuild,pcbddc->sub_schurs_threshold);CHKERRQ(ierr);
1243b1b3d7a2SStefano Zampini       }
1244*b7eb3628SStefano Zampini       ierr = PCBDDCSetUpSubSchurs(pc,pcbddc->sub_schurs_layers,pcbddc->sub_schurs_use_useradj,adaptive);CHKERRQ(ierr);
1245*b7eb3628SStefano Zampini #if 0
1246*b7eb3628SStefano Zampini       if (adaptive) {
1247*b7eb3628SStefano Zampini         ierr = PCBDDCAdaptiveSelection(pc,&ad_nearnullspace);CHKERRQ(ierr);
1248*b7eb3628SStefano Zampini       }
1249*b7eb3628SStefano Zampini #endif
1250b1b3d7a2SStefano Zampini     }
1251684f6988SStefano Zampini   }
1252684f6988SStefano Zampini 
1253f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1254fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1255f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1256f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1257f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1258f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1259f4ddd8eeSStefano Zampini     } else {
1260f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1261f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1262f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1263165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1264f4ddd8eeSStefano Zampini         PetscInt         i;
1265165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1266165b64e2SStefano Zampini         PetscObjectState state;
1267165b64e2SStefano Zampini         PetscInt         nnsp_size;
1268165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1269f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1270f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1271165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1272f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1273f4ddd8eeSStefano Zampini             break;
1274f4ddd8eeSStefano Zampini           }
1275f4ddd8eeSStefano Zampini         }
1276f4ddd8eeSStefano Zampini       }
1277f4ddd8eeSStefano Zampini     }
1278f4ddd8eeSStefano Zampini   } else {
1279f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1280f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1281f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1282f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1283f4ddd8eeSStefano Zampini     }
1284f4ddd8eeSStefano Zampini   }
1285*b7eb3628SStefano Zampini #if 0
1286*b7eb3628SStefano Zampini   if (adaptive) {
1287*b7eb3628SStefano Zampini     /* join nullspaces */
1288*b7eb3628SStefano Zampini   }
1289*b7eb3628SStefano Zampini   /* set null space in BDDC */
1290*b7eb3628SStefano Zampini #endif
1291f4ddd8eeSStefano Zampini 
1292f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1293727cdba6SStefano Zampini   /* reset primal space flags */
1294f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1295727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
1296fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
1297727cdba6SStefano Zampini     /* It also sets the primal space flags */
1298674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1299e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1300f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
13015e8657edSStefano Zampini 
13025e8657edSStefano Zampini     if (pcbddc->use_change_of_basis) {
13035e8657edSStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
13045e8657edSStefano Zampini 
13055e8657edSStefano Zampini       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
13065e8657edSStefano Zampini       /* get submatrices */
13075e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
13085e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
13095e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
13105e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
13115e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
13125e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
13135e8657edSStefano Zampini     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1314b96c3477SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
13155e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
13165e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
13175e8657edSStefano Zampini     }
1318674ae819SStefano Zampini   }
1319fb8d54d4SStefano Zampini 
1320f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
1321b96c3477SStefano Zampini     /* SetUp coarse and local Neumann solvers */
132299cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1323b96c3477SStefano Zampini     /* SetUp Scaling operator */
1324674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
13250c7d97c5SJed Brown   }
132670cf5478SStefano Zampini 
132758a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
132858a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
13292b510759SStefano Zampini   }
13300c7d97c5SJed Brown   PetscFunctionReturn(0);
13310c7d97c5SJed Brown }
13320c7d97c5SJed Brown 
13330c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
13340c7d97c5SJed Brown /*
133550efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
13360c7d97c5SJed Brown 
13370c7d97c5SJed Brown    Input Parameters:
13380c7d97c5SJed Brown .  pc - the preconditioner context
13390c7d97c5SJed Brown .  r - input vector (global)
13400c7d97c5SJed Brown 
13410c7d97c5SJed Brown    Output Parameter:
13420c7d97c5SJed Brown .  z - output vector (global)
13430c7d97c5SJed Brown 
13440c7d97c5SJed Brown    Application Interface Routine: PCApply()
13450c7d97c5SJed Brown  */
13460c7d97c5SJed Brown #undef __FUNCT__
13470c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
134853cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
13490c7d97c5SJed Brown {
13500c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
13510c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
13520c7d97c5SJed Brown   PetscErrorCode    ierr;
13533b03a366Sstefano_zampini   const PetscScalar one = 1.0;
13543b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
13552617d88aSStefano Zampini   const PetscScalar zero = 0.0;
13560c7d97c5SJed Brown 
13570c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
13580c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
13598eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
13600c7d97c5SJed Brown 
13610c7d97c5SJed Brown   PetscFunctionBegin;
136285c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
13630c7d97c5SJed Brown     /* First Dirichlet solve */
13640c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13650c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
136653cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
13670c7d97c5SJed Brown     /*
13680c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1369674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1370674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
13710c7d97c5SJed Brown     */
13720c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
13730c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
13748eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
13750c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
13760c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
13770c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13780c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1379674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1380b76ba322SStefano Zampini   } else {
13810bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1382b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1383674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1384b76ba322SStefano Zampini   }
1385b76ba322SStefano Zampini 
13862617d88aSStefano Zampini   /* Apply interface preconditioner
13872617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1388dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
13892617d88aSStefano Zampini 
1390674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1391674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
13920c7d97c5SJed Brown 
13933b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
13940c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13950c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13960c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
13978eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1398df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1399df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1400df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1401df187020SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
14020c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14030c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14040c7d97c5SJed Brown   PetscFunctionReturn(0);
14050c7d97c5SJed Brown }
140650efa1b5SStefano Zampini 
140750efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
140850efa1b5SStefano Zampini /*
140950efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
141050efa1b5SStefano Zampini 
141150efa1b5SStefano Zampini    Input Parameters:
141250efa1b5SStefano Zampini .  pc - the preconditioner context
141350efa1b5SStefano Zampini .  r - input vector (global)
141450efa1b5SStefano Zampini 
141550efa1b5SStefano Zampini    Output Parameter:
141650efa1b5SStefano Zampini .  z - output vector (global)
141750efa1b5SStefano Zampini 
141850efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
141950efa1b5SStefano Zampini  */
142050efa1b5SStefano Zampini #undef __FUNCT__
142150efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
142250efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
142350efa1b5SStefano Zampini {
142450efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
142550efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
142650efa1b5SStefano Zampini   PetscErrorCode    ierr;
142750efa1b5SStefano Zampini   const PetscScalar one = 1.0;
142850efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
142950efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
143050efa1b5SStefano Zampini 
143150efa1b5SStefano Zampini   PetscFunctionBegin;
143250efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
143350efa1b5SStefano Zampini     /* First Dirichlet solve */
143450efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
143550efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
143650efa1b5SStefano Zampini     ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
143750efa1b5SStefano Zampini     /*
143850efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
143950efa1b5SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
144050efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
144150efa1b5SStefano Zampini     */
144250efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
144350efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
144450efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
144550efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
144650efa1b5SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
144750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
144850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
144950efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
145050efa1b5SStefano Zampini   } else {
145150efa1b5SStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
145250efa1b5SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
145350efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
145450efa1b5SStefano Zampini   }
145550efa1b5SStefano Zampini 
145650efa1b5SStefano Zampini   /* Apply interface preconditioner
145750efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1458dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
145950efa1b5SStefano Zampini 
146050efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
146150efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
146250efa1b5SStefano Zampini 
146350efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
146450efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
146550efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
146650efa1b5SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
146750efa1b5SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1468b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1469b0147a47SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1470b0147a47SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1471b0147a47SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
147250efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
147350efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
147450efa1b5SStefano Zampini   PetscFunctionReturn(0);
147550efa1b5SStefano Zampini }
1476da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1477674ae819SStefano Zampini 
1478da1bb401SStefano Zampini #undef __FUNCT__
1479da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1480da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1481da1bb401SStefano Zampini {
1482da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1483da1bb401SStefano Zampini   PetscErrorCode ierr;
1484da1bb401SStefano Zampini 
1485da1bb401SStefano Zampini   PetscFunctionBegin;
1486da1bb401SStefano Zampini   /* free data created by PCIS */
1487da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1488674ae819SStefano Zampini   /* free BDDC custom data  */
1489674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1490674ae819SStefano Zampini   /* destroy objects related to topography */
1491674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1492674ae819SStefano Zampini   /* free allocated graph structure */
1493da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1494b96c3477SStefano Zampini   /* free allocated sub schurs structure */
1495b96c3477SStefano Zampini   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
149634a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1497674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
149834a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1499674ae819SStefano Zampini   /* free solvers stuff */
1500674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
150162a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
150262a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
150362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1504906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1505906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1506906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1507906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1508906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1509906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1510906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1511906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1512906d46d4SStefano Zampini   }
1513906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
15143425bc38SStefano Zampini   /* remove functions */
1515906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1516674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1517bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
15182b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1519b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
15202b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1521bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1522bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
152382ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1524bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
152582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1526bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
152782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1528bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1529785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1530bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
153163602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1532bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1533bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1534bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1535bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1536674ae819SStefano Zampini   /* Free the private data structure */
1537674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1538da1bb401SStefano Zampini   PetscFunctionReturn(0);
1539da1bb401SStefano Zampini }
15403425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
15411e6b0712SBarry Smith 
15423425bc38SStefano Zampini #undef __FUNCT__
15433425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
15443425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
15453425bc38SStefano Zampini {
1546674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
1547c08af4c6SStefano Zampini   Vec            copy_standard_rhs;
15483425bc38SStefano Zampini   PC_IS*         pcis;
15493425bc38SStefano Zampini   PC_BDDC*       pcbddc;
15503425bc38SStefano Zampini   PetscErrorCode ierr;
15510c7d97c5SJed Brown 
15523425bc38SStefano Zampini   PetscFunctionBegin;
15533425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15543425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
15553425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
15563425bc38SStefano Zampini 
1557c08af4c6SStefano Zampini   /*
1558c08af4c6SStefano Zampini      change of basis for physical rhs if needed
1559c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
1560c08af4c6SStefano Zampini      TODO: better management when FETIDP will have its own class
1561c08af4c6SStefano Zampini   */
1562c08af4c6SStefano Zampini   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1563c08af4c6SStefano Zampini   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1564c08af4c6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
15653425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
1566c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1567c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1568fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1569fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
1570c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1571c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1572674ae819SStefano Zampini   /* Apply partition of unity */
15733425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1574c08af4c6SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
15758eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
15763425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
15773425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15783425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
15793425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1580c08af4c6SStefano Zampini     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1581c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1582c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1583c08af4c6SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1584c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1585c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
15863425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
15873425bc38SStefano Zampini   }
1588c08af4c6SStefano Zampini   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
15893425bc38SStefano Zampini   /* BDDC rhs */
15903425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
15918eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
15923425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15933425bc38SStefano Zampini   }
15943425bc38SStefano Zampini   /* apply BDDC */
1595dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
15963425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
15973425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
15983425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
15993425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16003425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16013425bc38SStefano Zampini   PetscFunctionReturn(0);
16023425bc38SStefano Zampini }
16031e6b0712SBarry Smith 
16043425bc38SStefano Zampini #undef __FUNCT__
16053425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
16063425bc38SStefano Zampini /*@
160728509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
16083425bc38SStefano Zampini 
16093425bc38SStefano Zampini    Collective
16103425bc38SStefano Zampini 
16113425bc38SStefano Zampini    Input Parameters:
161228509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
161328509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
16143425bc38SStefano Zampini 
16153425bc38SStefano Zampini    Output Parameters:
161628509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
16173425bc38SStefano Zampini 
16183425bc38SStefano Zampini    Level: developer
16193425bc38SStefano Zampini 
16203425bc38SStefano Zampini    Notes:
16213425bc38SStefano Zampini 
162228509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
16233425bc38SStefano Zampini @*/
16243425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
16253425bc38SStefano Zampini {
1626674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
16273425bc38SStefano Zampini   PetscErrorCode ierr;
16283425bc38SStefano Zampini 
16293425bc38SStefano Zampini   PetscFunctionBegin;
16303425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16313425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
16323425bc38SStefano Zampini   PetscFunctionReturn(0);
16333425bc38SStefano Zampini }
16343425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
16351e6b0712SBarry Smith 
16363425bc38SStefano Zampini #undef __FUNCT__
16373425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
16383425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
16393425bc38SStefano Zampini {
1640674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
16413425bc38SStefano Zampini   PC_IS*         pcis;
16423425bc38SStefano Zampini   PC_BDDC*       pcbddc;
16433425bc38SStefano Zampini   PetscErrorCode ierr;
16443425bc38SStefano Zampini 
16453425bc38SStefano Zampini   PetscFunctionBegin;
16463425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16473425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
16483425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
16493425bc38SStefano Zampini 
16503425bc38SStefano Zampini   /* apply B_delta^T */
16513425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16523425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16533425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
16543425bc38SStefano Zampini   /* compute rhs for BDDC application */
16553425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
16568eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
16573425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16583425bc38SStefano Zampini   }
16593425bc38SStefano Zampini   /* apply BDDC */
1660dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
16613425bc38SStefano Zampini   /* put values into standard global vector */
16623425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16633425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16648eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
16653425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
16663425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
16673425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
16683425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16693425bc38SStefano Zampini   }
16703425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16713425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
16723425bc38SStefano Zampini   /* final change of basis if needed
16733425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
16743308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
16753425bc38SStefano Zampini   PetscFunctionReturn(0);
16763425bc38SStefano Zampini }
16771e6b0712SBarry Smith 
16783425bc38SStefano Zampini #undef __FUNCT__
16793425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
16803425bc38SStefano Zampini /*@
168128509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
16823425bc38SStefano Zampini 
16833425bc38SStefano Zampini    Collective
16843425bc38SStefano Zampini 
16853425bc38SStefano Zampini    Input Parameters:
168628509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
168728509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
16883425bc38SStefano Zampini 
16893425bc38SStefano Zampini    Output Parameters:
169028509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
16913425bc38SStefano Zampini 
16923425bc38SStefano Zampini    Level: developer
16933425bc38SStefano Zampini 
16943425bc38SStefano Zampini    Notes:
16953425bc38SStefano Zampini 
169628509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
16973425bc38SStefano Zampini @*/
16983425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
16993425bc38SStefano Zampini {
1700674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17013425bc38SStefano Zampini   PetscErrorCode ierr;
17023425bc38SStefano Zampini 
17033425bc38SStefano Zampini   PetscFunctionBegin;
17043425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
17053425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
17063425bc38SStefano Zampini   PetscFunctionReturn(0);
17073425bc38SStefano Zampini }
17083425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
17091e6b0712SBarry Smith 
1710f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1711edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1712f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1713f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1714edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1715f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1716674ae819SStefano Zampini 
17173425bc38SStefano Zampini #undef __FUNCT__
17183425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
17193425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
17203425bc38SStefano Zampini {
1721674ae819SStefano Zampini 
1722674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
17233425bc38SStefano Zampini   Mat            newmat;
1724674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
17253425bc38SStefano Zampini   PC             newpc;
1726ce94432eSBarry Smith   MPI_Comm       comm;
17273425bc38SStefano Zampini   PetscErrorCode ierr;
17283425bc38SStefano Zampini 
17293425bc38SStefano Zampini   PetscFunctionBegin;
1730ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
17313425bc38SStefano Zampini   /* FETIDP linear matrix */
17323425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
17333425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
17343425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
17353425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1736edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
17373425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
17383425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
17393425bc38SStefano Zampini   /* FETIDP preconditioner */
17403425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
17413425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
17423425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
17433425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
17443425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
17453425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1746edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
17473425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
174823ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
17493425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
17503425bc38SStefano Zampini   /* return pointers for objects created */
17513425bc38SStefano Zampini   *fetidp_mat=newmat;
17523425bc38SStefano Zampini   *fetidp_pc=newpc;
17533425bc38SStefano Zampini   PetscFunctionReturn(0);
17543425bc38SStefano Zampini }
17551e6b0712SBarry Smith 
17563425bc38SStefano Zampini #undef __FUNCT__
17573425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
17583425bc38SStefano Zampini /*@
175928509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
17603425bc38SStefano Zampini 
17613425bc38SStefano Zampini    Collective
17623425bc38SStefano Zampini 
17633425bc38SStefano Zampini    Input Parameters:
176428509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
176528509bceSStefano Zampini 
176628509bceSStefano Zampini    Output Parameters:
176728509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
176828509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
176928509bceSStefano Zampini 
177028509bceSStefano Zampini    Options Database Keys:
177128509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
17723425bc38SStefano Zampini 
17733425bc38SStefano Zampini    Level: developer
17743425bc38SStefano Zampini 
17753425bc38SStefano Zampini    Notes:
177628509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
17773425bc38SStefano Zampini 
17783425bc38SStefano Zampini .seealso: PCBDDC
17793425bc38SStefano Zampini @*/
17803425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
17813425bc38SStefano Zampini {
17823425bc38SStefano Zampini   PetscErrorCode ierr;
17833425bc38SStefano Zampini 
17843425bc38SStefano Zampini   PetscFunctionBegin;
17853425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17863425bc38SStefano Zampini   if (pc->setupcalled) {
1787516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1788f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
17893425bc38SStefano Zampini   PetscFunctionReturn(0);
17903425bc38SStefano Zampini }
17910c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1792da1bb401SStefano Zampini /*MC
1793da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
17940c7d97c5SJed Brown 
179528509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
179628509bceSStefano Zampini 
179728509bceSStefano Zampini .vb
179828509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
179928509bceSStefano 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
180028509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
180128509bceSStefano Zampini .ve
180228509bceSStefano Zampini 
180328509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
180428509bceSStefano Zampini 
1805b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
180628509bceSStefano Zampini 
180728509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
180828509bceSStefano Zampini 
1809b6fdb6dfSStefano 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.
1810b6fdb6dfSStefano Zampini 
181128509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
181228509bceSStefano Zampini 
18136a818285SBarry Smith    Boundary nodes are split in vertices, edges and faces 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()
181428509bceSStefano Zampini 
18156a818285SBarry Smith    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().
181628509bceSStefano Zampini 
1817b6fdb6dfSStefano Zampini    Change of basis is performed similarly to [2] when requested. When more the 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.
181828509bceSStefano Zampini 
181928509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
182028509bceSStefano Zampini 
1821da1bb401SStefano Zampini    Options Database Keys:
182228509bceSStefano Zampini 
182328509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
182428509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1825b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
182628509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
182728509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
182828509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
182928509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
183028509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
183128509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
183228509bceSStefano Zampini 
183328509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
183428509bceSStefano Zampini .vb
183528509bceSStefano Zampini       -pc_bddc_dirichlet_
183628509bceSStefano Zampini       -pc_bddc_neumann_
183728509bceSStefano Zampini       -pc_bddc_coarse_
183828509bceSStefano Zampini .ve
183928509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
184028509bceSStefano Zampini 
184128509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
184228509bceSStefano Zampini .vb
1843312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1844312be037SStefano Zampini       -pc_bddc_neumann_lN_
1845312be037SStefano Zampini       -pc_bddc_coarse_lN_
184628509bceSStefano Zampini .ve
1847312be037SStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N.
1848da1bb401SStefano Zampini 
1849da1bb401SStefano Zampini    Level: intermediate
1850da1bb401SStefano Zampini 
1851b6fdb6dfSStefano Zampini    Developer notes:
1852da1bb401SStefano Zampini 
185328509bceSStefano Zampini    New deluxe scaling operator will be available soon.
1854da1bb401SStefano Zampini 
1855da1bb401SStefano Zampini    Contributed by Stefano Zampini
1856da1bb401SStefano Zampini 
1857da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1858da1bb401SStefano Zampini M*/
1859b2573a8aSBarry Smith 
1860da1bb401SStefano Zampini #undef __FUNCT__
1861da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
18628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1863da1bb401SStefano Zampini {
1864da1bb401SStefano Zampini   PetscErrorCode      ierr;
1865da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1866da1bb401SStefano Zampini 
1867da1bb401SStefano Zampini   PetscFunctionBegin;
1868da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1869b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1870da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1871da1bb401SStefano Zampini 
1872da1bb401SStefano Zampini   /* create PCIS data structure */
1873da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1874da1bb401SStefano Zampini 
187547d04d0dSStefano Zampini   /* BDDC customization */
187608a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
187747d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
187847d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
187947d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
188047d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
188147d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
188247d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
1883fa434743SStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE;
1884fa434743SStefano Zampini   pcbddc->use_qr_single       = PETSC_FALSE;
188547d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1886b9d89cd5SStefano Zampini   /* private */
1887b9d89cd5SStefano Zampini   pcbddc->issym                      = PETSC_FALSE;
1888727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1889e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1890e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1891e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1892727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1893fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
189468457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
1895f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1896727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1897f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1898f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1899f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1900674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
19010bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
19023972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1903534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1904534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1905534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1906b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
1907906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
1908da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1909da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1910da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1911da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1912da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
191315aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
191415aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1915da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1916da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1917da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1918da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1919da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1920da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1921da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1922da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1923da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1924da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1925785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1926785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1927785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
192860d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
192960d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
193063602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
1931da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
193263602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
1933da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
193485c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
193547d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
193647d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
19374fad6a16SStefano Zampini   pcbddc->current_level              = 0;
19382b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1939323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
194074e2c79eSStefano Zampini   pcbddc->redistribute_coarse        = 0;
1941f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
1942323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
1943da1bb401SStefano Zampini 
1944674ae819SStefano Zampini   /* create local graph structure */
1945674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1946674ae819SStefano Zampini 
1947674ae819SStefano Zampini   /* scaling */
1948674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
194934a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
1950b96c3477SStefano Zampini 
1951b96c3477SStefano Zampini   /* create sub schurs structure */
1952b96c3477SStefano Zampini   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
1953b96c3477SStefano Zampini   pcbddc->sub_schurs_threshold   = -1;
1954b96c3477SStefano Zampini   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
1955b96c3477SStefano Zampini   pcbddc->sub_schurs_layers      = -1;
1956b96c3477SStefano Zampini   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
1957b96c3477SStefano Zampini 
1958b96c3477SStefano Zampini   pcbddc->computed_rowadj = PETSC_FALSE;
1959da1bb401SStefano Zampini 
1960*b7eb3628SStefano Zampini   /* adaptivity */
1961*b7eb3628SStefano Zampini   pcbddc->adaptive_threshold = -1.0;
1962*b7eb3628SStefano Zampini 
1963da1bb401SStefano Zampini   /* function pointers */
1964da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
196593bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1966da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1967da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1968da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1969da1bb401SStefano Zampini   pc->ops->view                = 0;
1970da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1971da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1972da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1973534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1974534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1975da1bb401SStefano Zampini 
1976da1bb401SStefano Zampini   /* composing function */
1977906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
1978674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1979bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
19802b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1981b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
19822b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1983bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1984bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
198582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1986bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
198782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1988bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
198982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1990bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
199182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1992bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
199363602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1994bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1995bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1996bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1997bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1998da1bb401SStefano Zampini   PetscFunctionReturn(0);
1999da1bb401SStefano Zampini }
20003425bc38SStefano Zampini 
2001