xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision d5574798c8464fdb92aede3bf97f4457bb71bbdc)
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"
468c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc)
470c7d97c5SJed Brown {
480c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
490c7d97c5SJed Brown   PetscErrorCode ierr;
500c7d97c5SJed Brown 
510c7d97c5SJed Brown   PetscFunctionBegin;
52e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"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);
61ac632422SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr);
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 = 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);
76b96c3477SStefano 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);
77b96c3477SStefano 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);
78ac632422SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);CHKERRQ(ierr);
794c6709b3SStefano Zampini   ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr);
8008122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr);
8108122e43SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr);
829552c7c7SStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_adaptive_invert_stildas","Whether or not to invert the Stildas matrices","none",pcbddc->adaptive_invert_Stildas,&pcbddc->adaptive_invert_Stildas,NULL);CHKERRQ(ierr);
830c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
840c7d97c5SJed Brown   PetscFunctionReturn(0);
850c7d97c5SJed Brown }
860c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
87674ae819SStefano Zampini #undef __FUNCT__
88906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
89906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
90b9b85e73SStefano Zampini {
91b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
92b9b85e73SStefano Zampini   PetscErrorCode ierr;
93b9b85e73SStefano Zampini 
94b9b85e73SStefano Zampini   PetscFunctionBegin;
95b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
96b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
97b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
98b9b85e73SStefano Zampini   PetscFunctionReturn(0);
99b9b85e73SStefano Zampini }
100b9b85e73SStefano Zampini #undef __FUNCT__
101906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
102b9b85e73SStefano Zampini /*@
103906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
104b9b85e73SStefano Zampini 
105b9b85e73SStefano Zampini    Collective on PC
106b9b85e73SStefano Zampini 
107b9b85e73SStefano Zampini    Input Parameters:
108b9b85e73SStefano Zampini +  pc - the preconditioning context
109906d46d4SStefano Zampini -  change - the change of basis matrix
110b9b85e73SStefano Zampini 
111b9b85e73SStefano Zampini    Level: intermediate
112b9b85e73SStefano Zampini 
113b9b85e73SStefano Zampini    Notes:
114b9b85e73SStefano Zampini 
115b9b85e73SStefano Zampini .seealso: PCBDDC
116b9b85e73SStefano Zampini @*/
117906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
118b9b85e73SStefano Zampini {
119b9b85e73SStefano Zampini   PetscErrorCode ierr;
120b9b85e73SStefano Zampini 
121b9b85e73SStefano Zampini   PetscFunctionBegin;
122b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
123b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
124906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
125906d46d4SStefano Zampini   if (pc->mat) {
126906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
127906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
128906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
129906d46d4SStefano Zampini     if (rows_c != rows) {
130906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
131906d46d4SStefano Zampini     }
132906d46d4SStefano Zampini     if (cols_c != cols) {
133906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
134906d46d4SStefano Zampini     }
135906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
136906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
137906d46d4SStefano Zampini     if (rows_c != rows) {
138906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
139906d46d4SStefano Zampini     }
140906d46d4SStefano Zampini     if (cols_c != cols) {
141906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
142906d46d4SStefano Zampini     }
143906d46d4SStefano Zampini   }
144906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
145b9b85e73SStefano Zampini   PetscFunctionReturn(0);
146b9b85e73SStefano Zampini }
147b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
148b9b85e73SStefano Zampini #undef __FUNCT__
149674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
150674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
151674ae819SStefano Zampini {
152674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
153674ae819SStefano Zampini   PetscErrorCode ierr;
1541e6b0712SBarry Smith 
155674ae819SStefano Zampini   PetscFunctionBegin;
156674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
157674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
158674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
159674ae819SStefano Zampini   PetscFunctionReturn(0);
160674ae819SStefano Zampini }
161674ae819SStefano Zampini #undef __FUNCT__
162674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
163674ae819SStefano Zampini /*@
16428509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
165674ae819SStefano Zampini 
166674ae819SStefano Zampini    Not collective
167674ae819SStefano Zampini 
168674ae819SStefano Zampini    Input Parameters:
169674ae819SStefano Zampini +  pc - the preconditioning context
17028509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
171674ae819SStefano Zampini 
172674ae819SStefano Zampini    Level: intermediate
173674ae819SStefano Zampini 
174674ae819SStefano Zampini    Notes:
175674ae819SStefano Zampini 
176674ae819SStefano Zampini .seealso: PCBDDC
177674ae819SStefano Zampini @*/
178674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
179674ae819SStefano Zampini {
180674ae819SStefano Zampini   PetscErrorCode ierr;
181674ae819SStefano Zampini 
182674ae819SStefano Zampini   PetscFunctionBegin;
183674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
184674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
185674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
186674ae819SStefano Zampini   PetscFunctionReturn(0);
187674ae819SStefano Zampini }
188674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1890c7d97c5SJed Brown #undef __FUNCT__
1904fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1914fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1924fad6a16SStefano Zampini {
1934fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1944fad6a16SStefano Zampini 
1954fad6a16SStefano Zampini   PetscFunctionBegin;
1964fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1974fad6a16SStefano Zampini   PetscFunctionReturn(0);
1984fad6a16SStefano Zampini }
1991e6b0712SBarry Smith 
2004fad6a16SStefano Zampini #undef __FUNCT__
2014fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
2024fad6a16SStefano Zampini /*@
20328509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
2044fad6a16SStefano Zampini 
2054fad6a16SStefano Zampini    Logically collective on PC
2064fad6a16SStefano Zampini 
2074fad6a16SStefano Zampini    Input Parameters:
2084fad6a16SStefano Zampini +  pc - the preconditioning context
20928509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2104fad6a16SStefano Zampini 
21128509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2124fad6a16SStefano Zampini 
2134fad6a16SStefano Zampini    Level: intermediate
2144fad6a16SStefano Zampini 
2154fad6a16SStefano Zampini    Notes:
2164fad6a16SStefano Zampini 
2174fad6a16SStefano Zampini .seealso: PCBDDC
2184fad6a16SStefano Zampini @*/
2194fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2204fad6a16SStefano Zampini {
2214fad6a16SStefano Zampini   PetscErrorCode ierr;
2224fad6a16SStefano Zampini 
2234fad6a16SStefano Zampini   PetscFunctionBegin;
2244fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2252b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2264fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2274fad6a16SStefano Zampini   PetscFunctionReturn(0);
2284fad6a16SStefano Zampini }
2292b510759SStefano Zampini 
230b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2312b510759SStefano Zampini #undef __FUNCT__
232b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
233b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
234b8ffe317SStefano Zampini {
235b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
236b8ffe317SStefano Zampini 
237b8ffe317SStefano Zampini   PetscFunctionBegin;
23885c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
239b8ffe317SStefano Zampini   PetscFunctionReturn(0);
240b8ffe317SStefano Zampini }
241b8ffe317SStefano Zampini 
242b8ffe317SStefano Zampini #undef __FUNCT__
243b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
244b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2452b510759SStefano Zampini {
2462b510759SStefano Zampini   PetscErrorCode ierr;
2472b510759SStefano Zampini 
2482b510759SStefano Zampini   PetscFunctionBegin;
2492b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
250b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
251b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2522b510759SStefano Zampini   PetscFunctionReturn(0);
2532b510759SStefano Zampini }
2541e6b0712SBarry Smith 
2554fad6a16SStefano Zampini #undef __FUNCT__
2562b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2572b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2584fad6a16SStefano Zampini {
2594fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2604fad6a16SStefano Zampini 
2614fad6a16SStefano Zampini   PetscFunctionBegin;
2622b510759SStefano Zampini   pcbddc->current_level = level;
2634fad6a16SStefano Zampini   PetscFunctionReturn(0);
2644fad6a16SStefano Zampini }
2651e6b0712SBarry Smith 
2664fad6a16SStefano Zampini #undef __FUNCT__
267b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
268b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
269b8ffe317SStefano Zampini {
270b8ffe317SStefano Zampini   PetscErrorCode ierr;
271b8ffe317SStefano Zampini 
272b8ffe317SStefano Zampini   PetscFunctionBegin;
273b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
274b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
275b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
276b8ffe317SStefano Zampini   PetscFunctionReturn(0);
277b8ffe317SStefano Zampini }
278b8ffe317SStefano Zampini 
279b8ffe317SStefano Zampini #undef __FUNCT__
2802b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2812b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2822b510759SStefano Zampini {
2832b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2842b510759SStefano Zampini 
2852b510759SStefano Zampini   PetscFunctionBegin;
2862b510759SStefano Zampini   pcbddc->max_levels = levels;
2872b510759SStefano Zampini   PetscFunctionReturn(0);
2882b510759SStefano Zampini }
2892b510759SStefano Zampini 
2902b510759SStefano Zampini #undef __FUNCT__
2912b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2924fad6a16SStefano Zampini /*@
29328509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2944fad6a16SStefano Zampini 
2954fad6a16SStefano Zampini    Logically collective on PC
2964fad6a16SStefano Zampini 
2974fad6a16SStefano Zampini    Input Parameters:
2984fad6a16SStefano Zampini +  pc - the preconditioning context
29928509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
3004fad6a16SStefano Zampini 
30128509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
3024fad6a16SStefano Zampini 
3034fad6a16SStefano Zampini    Level: intermediate
3044fad6a16SStefano Zampini 
3054fad6a16SStefano Zampini    Notes:
3064fad6a16SStefano Zampini 
3074fad6a16SStefano Zampini .seealso: PCBDDC
3084fad6a16SStefano Zampini @*/
3092b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3104fad6a16SStefano Zampini {
3114fad6a16SStefano Zampini   PetscErrorCode ierr;
3124fad6a16SStefano Zampini 
3134fad6a16SStefano Zampini   PetscFunctionBegin;
3144fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3152b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
316312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3172b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3184fad6a16SStefano Zampini   PetscFunctionReturn(0);
3194fad6a16SStefano Zampini }
3204fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3211e6b0712SBarry Smith 
3224fad6a16SStefano Zampini #undef __FUNCT__
3230bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3240bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3250bdf917eSStefano Zampini {
3260bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3270bdf917eSStefano Zampini   PetscErrorCode ierr;
3280bdf917eSStefano Zampini 
3290bdf917eSStefano Zampini   PetscFunctionBegin;
3300bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3310bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3320bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3330bdf917eSStefano Zampini   PetscFunctionReturn(0);
3340bdf917eSStefano Zampini }
3351e6b0712SBarry Smith 
3360bdf917eSStefano Zampini #undef __FUNCT__
3370bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3380bdf917eSStefano Zampini /*@
33928509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3400bdf917eSStefano Zampini 
3410bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3420bdf917eSStefano Zampini 
3430bdf917eSStefano Zampini    Input Parameters:
3440bdf917eSStefano Zampini +  pc - the preconditioning context
34528509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3460bdf917eSStefano Zampini 
3470bdf917eSStefano Zampini    Level: intermediate
3480bdf917eSStefano Zampini 
3490bdf917eSStefano Zampini    Notes:
3500bdf917eSStefano Zampini 
3510bdf917eSStefano Zampini .seealso: PCBDDC
3520bdf917eSStefano Zampini @*/
3530bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3540bdf917eSStefano Zampini {
3550bdf917eSStefano Zampini   PetscErrorCode ierr;
3560bdf917eSStefano Zampini 
3570bdf917eSStefano Zampini   PetscFunctionBegin;
3580bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
359674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
3602b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
3610bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
3620bdf917eSStefano Zampini   PetscFunctionReturn(0);
3630bdf917eSStefano Zampini }
3640bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3651e6b0712SBarry Smith 
3660bdf917eSStefano Zampini #undef __FUNCT__
3673b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3683b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3693b03a366Sstefano_zampini {
3703b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3713b03a366Sstefano_zampini   PetscErrorCode ierr;
3723b03a366Sstefano_zampini 
3733b03a366Sstefano_zampini   PetscFunctionBegin;
374785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
375785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3763b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
37736e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
37836e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
379fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3803b03a366Sstefano_zampini   PetscFunctionReturn(0);
3813b03a366Sstefano_zampini }
3821e6b0712SBarry Smith 
3833b03a366Sstefano_zampini #undef __FUNCT__
3843b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3853b03a366Sstefano_zampini /*@
38628509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3873b03a366Sstefano_zampini 
388785d1243SStefano Zampini    Collective
3893b03a366Sstefano_zampini 
3903b03a366Sstefano_zampini    Input Parameters:
3913b03a366Sstefano_zampini +  pc - the preconditioning context
392785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3933b03a366Sstefano_zampini 
3943b03a366Sstefano_zampini    Level: intermediate
3953b03a366Sstefano_zampini 
396785d1243SStefano Zampini    Notes: Any process can list any global node
3973b03a366Sstefano_zampini 
3983b03a366Sstefano_zampini .seealso: PCBDDC
3993b03a366Sstefano_zampini @*/
4003b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
4013b03a366Sstefano_zampini {
4023b03a366Sstefano_zampini   PetscErrorCode ierr;
4033b03a366Sstefano_zampini 
4043b03a366Sstefano_zampini   PetscFunctionBegin;
4053b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
406674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
407785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4083b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4093b03a366Sstefano_zampini   PetscFunctionReturn(0);
4103b03a366Sstefano_zampini }
4113b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4121e6b0712SBarry Smith 
4133b03a366Sstefano_zampini #undef __FUNCT__
41482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
41582ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4160c7d97c5SJed Brown {
4170c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4180c7d97c5SJed Brown   PetscErrorCode ierr;
4190c7d97c5SJed Brown 
4200c7d97c5SJed Brown   PetscFunctionBegin;
421785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
422785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4230c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4240c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
425785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4267d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4270c7d97c5SJed Brown   PetscFunctionReturn(0);
4280c7d97c5SJed Brown }
4290c7d97c5SJed Brown 
4300c7d97c5SJed Brown #undef __FUNCT__
43182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4329c0446d6SStefano Zampini /*@
43382ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4349c0446d6SStefano Zampini 
435785d1243SStefano Zampini    Collective
43653cdbc3dSStefano Zampini 
43753cdbc3dSStefano Zampini    Input Parameters:
43853cdbc3dSStefano Zampini +  pc - the preconditioning context
43982ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
44053cdbc3dSStefano Zampini 
44153cdbc3dSStefano Zampini    Level: intermediate
44253cdbc3dSStefano Zampini 
4439c0446d6SStefano Zampini    Notes:
44453cdbc3dSStefano Zampini 
44553cdbc3dSStefano Zampini .seealso: PCBDDC
44653cdbc3dSStefano Zampini @*/
44782ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4480c7d97c5SJed Brown {
4490c7d97c5SJed Brown   PetscErrorCode ierr;
4500c7d97c5SJed Brown 
4510c7d97c5SJed Brown   PetscFunctionBegin;
4520c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4530c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
45482ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
45582ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4560c7d97c5SJed Brown   PetscFunctionReturn(0);
4570c7d97c5SJed Brown }
4580c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4590c7d97c5SJed Brown 
4600c7d97c5SJed Brown #undef __FUNCT__
4610c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
46253cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4630c7d97c5SJed Brown {
4640c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
46553cdbc3dSStefano Zampini   PetscErrorCode ierr;
4660c7d97c5SJed Brown 
4670c7d97c5SJed Brown   PetscFunctionBegin;
468785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
469785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
47053cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
47136e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
47236e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
473fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4740c7d97c5SJed Brown   PetscFunctionReturn(0);
4750c7d97c5SJed Brown }
4761e6b0712SBarry Smith 
4770c7d97c5SJed Brown #undef __FUNCT__
4780c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
47957527edcSJed Brown /*@
48028509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
48157527edcSJed Brown 
482785d1243SStefano Zampini    Collective
48357527edcSJed Brown 
48457527edcSJed Brown    Input Parameters:
48557527edcSJed Brown +  pc - the preconditioning context
486785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
48757527edcSJed Brown 
48857527edcSJed Brown    Level: intermediate
48957527edcSJed Brown 
490785d1243SStefano Zampini    Notes: Any process can list any global node
49157527edcSJed Brown 
49257527edcSJed Brown .seealso: PCBDDC
49357527edcSJed Brown @*/
49453cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4950c7d97c5SJed Brown {
4960c7d97c5SJed Brown   PetscErrorCode ierr;
4970c7d97c5SJed Brown 
4980c7d97c5SJed Brown   PetscFunctionBegin;
4990c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
500674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
501785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
50253cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
50353cdbc3dSStefano Zampini   PetscFunctionReturn(0);
50453cdbc3dSStefano Zampini }
50553cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5061e6b0712SBarry Smith 
50753cdbc3dSStefano Zampini #undef __FUNCT__
50882ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
50982ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5100c7d97c5SJed Brown {
5110c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5120c7d97c5SJed Brown   PetscErrorCode ierr;
5130c7d97c5SJed Brown 
5140c7d97c5SJed Brown   PetscFunctionBegin;
515785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
516785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5170c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5180c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
519785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5207d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5210c7d97c5SJed Brown   PetscFunctionReturn(0);
5220c7d97c5SJed Brown }
5230c7d97c5SJed Brown 
5240c7d97c5SJed Brown #undef __FUNCT__
52582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5260c7d97c5SJed Brown /*@
52782ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5280c7d97c5SJed Brown 
529785d1243SStefano Zampini    Collective
5300c7d97c5SJed Brown 
5310c7d97c5SJed Brown    Input Parameters:
5320c7d97c5SJed Brown +  pc - the preconditioning context
53382ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5340c7d97c5SJed Brown 
5350c7d97c5SJed Brown    Level: intermediate
5360c7d97c5SJed Brown 
5370c7d97c5SJed Brown    Notes:
5380c7d97c5SJed Brown 
5390c7d97c5SJed Brown .seealso: PCBDDC
5400c7d97c5SJed Brown @*/
54182ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5420c7d97c5SJed Brown {
5430c7d97c5SJed Brown   PetscErrorCode ierr;
5440c7d97c5SJed Brown 
5450c7d97c5SJed Brown   PetscFunctionBegin;
5460c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5470c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
54882ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54982ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
55053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
55153cdbc3dSStefano Zampini }
55253cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
55353cdbc3dSStefano Zampini 
55453cdbc3dSStefano Zampini #undef __FUNCT__
555da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
556da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
557da1bb401SStefano Zampini {
558da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
559da1bb401SStefano Zampini 
560da1bb401SStefano Zampini   PetscFunctionBegin;
561da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
562da1bb401SStefano Zampini   PetscFunctionReturn(0);
563da1bb401SStefano Zampini }
5641e6b0712SBarry Smith 
565da1bb401SStefano Zampini #undef __FUNCT__
566da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
567da1bb401SStefano Zampini /*@
568785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
569da1bb401SStefano Zampini 
570785d1243SStefano Zampini    Collective
571785d1243SStefano Zampini 
572785d1243SStefano Zampini    Input Parameters:
573785d1243SStefano Zampini .  pc - the preconditioning context
574785d1243SStefano Zampini 
575785d1243SStefano Zampini    Output Parameters:
576785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
577785d1243SStefano Zampini 
578785d1243SStefano Zampini    Level: intermediate
579785d1243SStefano Zampini 
580785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
581785d1243SStefano Zampini 
582785d1243SStefano Zampini .seealso: PCBDDC
583785d1243SStefano Zampini @*/
584785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
585785d1243SStefano Zampini {
586785d1243SStefano Zampini   PetscErrorCode ierr;
587785d1243SStefano Zampini 
588785d1243SStefano Zampini   PetscFunctionBegin;
589785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
590785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
591785d1243SStefano Zampini   PetscFunctionReturn(0);
592785d1243SStefano Zampini }
593785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
594785d1243SStefano Zampini 
595785d1243SStefano Zampini #undef __FUNCT__
596785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
597785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
598785d1243SStefano Zampini {
599785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
600785d1243SStefano Zampini 
601785d1243SStefano Zampini   PetscFunctionBegin;
602785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
603785d1243SStefano Zampini   PetscFunctionReturn(0);
604785d1243SStefano Zampini }
605785d1243SStefano Zampini 
606785d1243SStefano Zampini #undef __FUNCT__
60782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
608da1bb401SStefano Zampini /*@
60982ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
610da1bb401SStefano Zampini 
611785d1243SStefano Zampini    Collective
612da1bb401SStefano Zampini 
613da1bb401SStefano Zampini    Input Parameters:
61428509bceSStefano Zampini .  pc - the preconditioning context
615da1bb401SStefano Zampini 
616da1bb401SStefano Zampini    Output Parameters:
61728509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
618da1bb401SStefano Zampini 
619da1bb401SStefano Zampini    Level: intermediate
620da1bb401SStefano Zampini 
621da1bb401SStefano Zampini    Notes:
622da1bb401SStefano Zampini 
623da1bb401SStefano Zampini .seealso: PCBDDC
624da1bb401SStefano Zampini @*/
62582ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
626da1bb401SStefano Zampini {
627da1bb401SStefano Zampini   PetscErrorCode ierr;
628da1bb401SStefano Zampini 
629da1bb401SStefano Zampini   PetscFunctionBegin;
630da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
63182ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
632da1bb401SStefano Zampini   PetscFunctionReturn(0);
633da1bb401SStefano Zampini }
634da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6351e6b0712SBarry Smith 
636da1bb401SStefano Zampini #undef __FUNCT__
63753cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
63853cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
63953cdbc3dSStefano Zampini {
64053cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
64153cdbc3dSStefano Zampini 
64253cdbc3dSStefano Zampini   PetscFunctionBegin;
64353cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
64453cdbc3dSStefano Zampini   PetscFunctionReturn(0);
64553cdbc3dSStefano Zampini }
6461e6b0712SBarry Smith 
64753cdbc3dSStefano Zampini #undef __FUNCT__
64853cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
64953cdbc3dSStefano Zampini /*@
650785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
65153cdbc3dSStefano Zampini 
652785d1243SStefano Zampini    Collective
653785d1243SStefano Zampini 
654785d1243SStefano Zampini    Input Parameters:
655785d1243SStefano Zampini .  pc - the preconditioning context
656785d1243SStefano Zampini 
657785d1243SStefano Zampini    Output Parameters:
658785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
659785d1243SStefano Zampini 
660785d1243SStefano Zampini    Level: intermediate
661785d1243SStefano Zampini 
662785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
663785d1243SStefano Zampini 
664785d1243SStefano Zampini .seealso: PCBDDC
665785d1243SStefano Zampini @*/
666785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
667785d1243SStefano Zampini {
668785d1243SStefano Zampini   PetscErrorCode ierr;
669785d1243SStefano Zampini 
670785d1243SStefano Zampini   PetscFunctionBegin;
671785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
672785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
673785d1243SStefano Zampini   PetscFunctionReturn(0);
674785d1243SStefano Zampini }
675785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
676785d1243SStefano Zampini 
677785d1243SStefano Zampini #undef __FUNCT__
678785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
679785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
680785d1243SStefano Zampini {
681785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
682785d1243SStefano Zampini 
683785d1243SStefano Zampini   PetscFunctionBegin;
684785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
685785d1243SStefano Zampini   PetscFunctionReturn(0);
686785d1243SStefano Zampini }
687785d1243SStefano Zampini 
688785d1243SStefano Zampini #undef __FUNCT__
68982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
69053cdbc3dSStefano Zampini /*@
69182ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
69253cdbc3dSStefano Zampini 
693785d1243SStefano Zampini    Collective
69453cdbc3dSStefano Zampini 
69553cdbc3dSStefano Zampini    Input Parameters:
69628509bceSStefano Zampini .  pc - the preconditioning context
69753cdbc3dSStefano Zampini 
69853cdbc3dSStefano Zampini    Output Parameters:
69928509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
70053cdbc3dSStefano Zampini 
70153cdbc3dSStefano Zampini    Level: intermediate
70253cdbc3dSStefano Zampini 
70353cdbc3dSStefano Zampini    Notes:
70453cdbc3dSStefano Zampini 
70553cdbc3dSStefano Zampini .seealso: PCBDDC
70653cdbc3dSStefano Zampini @*/
70782ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
70853cdbc3dSStefano Zampini {
70953cdbc3dSStefano Zampini   PetscErrorCode ierr;
71053cdbc3dSStefano Zampini 
71153cdbc3dSStefano Zampini   PetscFunctionBegin;
71253cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
71382ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7140c7d97c5SJed Brown   PetscFunctionReturn(0);
7150c7d97c5SJed Brown }
71636e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7171e6b0712SBarry Smith 
71836e030ebSStefano Zampini #undef __FUNCT__
719da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7201a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
72136e030ebSStefano Zampini {
72236e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
723da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
724da1bb401SStefano Zampini   PetscErrorCode ierr;
72536e030ebSStefano Zampini 
72636e030ebSStefano Zampini   PetscFunctionBegin;
727674ae819SStefano Zampini   /* free old CSR */
728674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
729fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
730674ae819SStefano Zampini   /* get CSR into graph structure */
731da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
732854ce69bSBarry Smith     ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr);
733785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
734674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
735674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
736da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7371a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7381a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
739674ae819SStefano Zampini   }
740575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
74136e030ebSStefano Zampini   PetscFunctionReturn(0);
74236e030ebSStefano Zampini }
7431e6b0712SBarry Smith 
74436e030ebSStefano Zampini #undef __FUNCT__
745da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
74636e030ebSStefano Zampini /*@
74728509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
74836e030ebSStefano Zampini 
74936e030ebSStefano Zampini    Not collective
75036e030ebSStefano Zampini 
75136e030ebSStefano Zampini    Input Parameters:
75236e030ebSStefano Zampini +  pc - the preconditioning context
75328509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
75428509bceSStefano Zampini .  xadj, adjncy - the CSR graph
75528509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
75636e030ebSStefano Zampini 
75736e030ebSStefano Zampini    Level: intermediate
75836e030ebSStefano Zampini 
75936e030ebSStefano Zampini    Notes:
76036e030ebSStefano Zampini 
76128509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
76236e030ebSStefano Zampini @*/
7631a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
76436e030ebSStefano Zampini {
765575ad6abSStefano Zampini   void (*f)(void) = 0;
76636e030ebSStefano Zampini   PetscErrorCode ierr;
76736e030ebSStefano Zampini 
76836e030ebSStefano Zampini   PetscFunctionBegin;
76936e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
770674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
771d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
772674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
773674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
774da1bb401SStefano Zampini   }
77536e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
776575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
777575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
778575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
779575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
780575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
78136e030ebSStefano Zampini   }
78236e030ebSStefano Zampini   PetscFunctionReturn(0);
78336e030ebSStefano Zampini }
7849c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7851e6b0712SBarry Smith 
7869c0446d6SStefano Zampini #undef __FUNCT__
78763602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
78863602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
78963602bcaSStefano Zampini {
79063602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
79163602bcaSStefano Zampini   PetscInt i;
79263602bcaSStefano Zampini   PetscErrorCode ierr;
79363602bcaSStefano Zampini 
79463602bcaSStefano Zampini   PetscFunctionBegin;
79563602bcaSStefano Zampini   /* Destroy ISes if they were already set */
79663602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
79763602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
79863602bcaSStefano Zampini   }
79963602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
80063602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
80163602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
80263602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
80363602bcaSStefano Zampini   }
80463602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
80563602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
80663602bcaSStefano Zampini   /* allocate space then set */
807d02579f5SStefano Zampini   if (n_is) {
808d02579f5SStefano Zampini     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
809d02579f5SStefano Zampini   }
81063602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
81163602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
81263602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
81363602bcaSStefano Zampini   }
81463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
81563602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8161a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
81763602bcaSStefano Zampini   PetscFunctionReturn(0);
81863602bcaSStefano Zampini }
81963602bcaSStefano Zampini 
82063602bcaSStefano Zampini #undef __FUNCT__
82163602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
82263602bcaSStefano Zampini /*@
82363602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
82463602bcaSStefano Zampini 
82563602bcaSStefano Zampini    Collective
82663602bcaSStefano Zampini 
82763602bcaSStefano Zampini    Input Parameters:
82863602bcaSStefano Zampini +  pc - the preconditioning context
82963602bcaSStefano Zampini -  n_is - number of index sets defining the fields
83063602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in local ordering
83163602bcaSStefano Zampini 
83263602bcaSStefano Zampini    Level: intermediate
83363602bcaSStefano Zampini 
83463602bcaSStefano 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.
83563602bcaSStefano Zampini 
83663602bcaSStefano Zampini .seealso: PCBDDC
83763602bcaSStefano Zampini @*/
83863602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
83963602bcaSStefano Zampini {
84063602bcaSStefano Zampini   PetscInt       i;
84163602bcaSStefano Zampini   PetscErrorCode ierr;
84263602bcaSStefano Zampini 
84363602bcaSStefano Zampini   PetscFunctionBegin;
84463602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
84563602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
84663602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
84763602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
84863602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
84963602bcaSStefano Zampini   }
850e71e7a71SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
85163602bcaSStefano Zampini   PetscFunctionReturn(0);
85263602bcaSStefano Zampini }
85363602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
85463602bcaSStefano Zampini 
85563602bcaSStefano Zampini #undef __FUNCT__
8569c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8579c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8589c0446d6SStefano Zampini {
8599c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8609c0446d6SStefano Zampini   PetscInt i;
8619c0446d6SStefano Zampini   PetscErrorCode ierr;
8629c0446d6SStefano Zampini 
8639c0446d6SStefano Zampini   PetscFunctionBegin;
864da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8659c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8669c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8679c0446d6SStefano Zampini   }
868d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
86963602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
87063602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
87163602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
87263602bcaSStefano Zampini   }
87363602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
87463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
875da1bb401SStefano Zampini   /* allocate space then set */
876d02579f5SStefano Zampini   if (n_is) {
877785e854fSJed Brown     ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
878d02579f5SStefano Zampini   }
8799c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
880da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
881da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8829c0446d6SStefano Zampini   }
8839c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
88463602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8851a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8869c0446d6SStefano Zampini   PetscFunctionReturn(0);
8879c0446d6SStefano Zampini }
8881e6b0712SBarry Smith 
8899c0446d6SStefano Zampini #undef __FUNCT__
8909c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8919c0446d6SStefano Zampini /*@
89263602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8939c0446d6SStefano Zampini 
89463602bcaSStefano Zampini    Collective
8959c0446d6SStefano Zampini 
8969c0446d6SStefano Zampini    Input Parameters:
8979c0446d6SStefano Zampini +  pc - the preconditioning context
89828509bceSStefano Zampini -  n_is - number of index sets defining the fields
89963602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in global ordering
9009c0446d6SStefano Zampini 
9019c0446d6SStefano Zampini    Level: intermediate
9029c0446d6SStefano Zampini 
90363602bcaSStefano Zampini    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
9049c0446d6SStefano Zampini 
9059c0446d6SStefano Zampini .seealso: PCBDDC
9069c0446d6SStefano Zampini @*/
9079c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9089c0446d6SStefano Zampini {
9092b510759SStefano Zampini   PetscInt       i;
9109c0446d6SStefano Zampini   PetscErrorCode ierr;
9119c0446d6SStefano Zampini 
9129c0446d6SStefano Zampini   PetscFunctionBegin;
9139c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
91463602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9152b510759SStefano Zampini   for (i=0;i<n_is;i++) {
91663602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
91763602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9182b510759SStefano Zampini   }
9199c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9209c0446d6SStefano Zampini   PetscFunctionReturn(0);
9219c0446d6SStefano Zampini }
922906d46d4SStefano Zampini 
923da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
924534831adSStefano Zampini #undef __FUNCT__
925534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
926534831adSStefano Zampini /* -------------------------------------------------------------------------- */
927534831adSStefano Zampini /*
928534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
929534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9309c0446d6SStefano Zampini 
931534831adSStefano Zampini    Input Parameter:
932534831adSStefano Zampini +  pc - the preconditioner contex
933534831adSStefano Zampini 
934534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
935534831adSStefano Zampini 
936534831adSStefano Zampini    Notes:
937534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
938534831adSStefano Zampini    the user, but instead is called by KSPSolve().
939534831adSStefano Zampini */
940534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
941534831adSStefano Zampini {
942534831adSStefano Zampini   PetscErrorCode ierr;
943534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
944534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9453972b0daSStefano Zampini   Vec            used_vec;
9468d00608fSStefano Zampini   PetscBool      copy_rhs = PETSC_TRUE;
947534831adSStefano Zampini 
948534831adSStefano Zampini   PetscFunctionBegin;
94985c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
95085c4d303SStefano Zampini   if (ksp) {
95185c4d303SStefano Zampini     PetscBool iscg;
95285c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
95385c4d303SStefano Zampini     if (!iscg) {
95485c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
95585c4d303SStefano Zampini     }
95685c4d303SStefano Zampini   }
95785c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
95862a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
95962a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
96062a6ff1dSStefano Zampini   }
96162a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
96262a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
96362a6ff1dSStefano Zampini   }
9648d00608fSStefano Zampini 
9653972b0daSStefano Zampini   if (x) {
9663972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9673972b0daSStefano Zampini     used_vec = x;
9688d00608fSStefano Zampini   } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
9693972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9703972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9713972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9723972b0daSStefano Zampini   }
9738efcfb23SStefano Zampini 
9748efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
9753972b0daSStefano Zampini   if (ksp) {
976a0cb1b98SStefano Zampini     /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
9778efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
9788efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
9793972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9803972b0daSStefano Zampini     }
9813972b0daSStefano Zampini   }
9823308cffdSStefano Zampini 
9838d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
9843972b0daSStefano Zampini 
9853972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
986b097fa66SStefano Zampini   if (rhs && pcbddc->DirichletBoundariesLocal) {
9873975b054SStefano Zampini     IS dirIS = NULL;
9883975b054SStefano Zampini 
9893975b054SStefano Zampini     /* DirichletBoundariesLocal may not be consistent among neighbours */
9903975b054SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
9913975b054SStefano Zampini     if (dirIS) {
992906d46d4SStefano Zampini       Mat_IS            *matis = (Mat_IS*)pc->pmat->data;
993785d1243SStefano Zampini       PetscInt          dirsize,i,*is_indices;
9942b095fd8SStefano Zampini       PetscScalar       *array_x;
9952b095fd8SStefano Zampini       const PetscScalar *array_diagonal;
996785d1243SStefano Zampini 
9973972b0daSStefano Zampini       ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
9983972b0daSStefano Zampini       ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
9993972b0daSStefano Zampini       ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10003972b0daSStefano Zampini       ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10013972b0daSStefano Zampini       ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10023972b0daSStefano Zampini       ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
100382ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
10043972b0daSStefano Zampini       ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10052b095fd8SStefano Zampini       ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10063972b0daSStefano Zampini       ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10072fa5cd67SKarl Rupp       for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
10083972b0daSStefano Zampini       ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
10092b095fd8SStefano Zampini       ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
10103972b0daSStefano Zampini       ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
10113972b0daSStefano Zampini       ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10123972b0daSStefano Zampini       ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10138d00608fSStefano Zampini       pcbddc->rhs_change = PETSC_TRUE;
10143975b054SStefano Zampini     }
10151b968477SStefano Zampini     ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
10168efcfb23SStefano Zampini   }
1017b76ba322SStefano Zampini 
10188efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
10198d00608fSStefano Zampini   if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
10208d00608fSStefano Zampini     /* store the original rhs */
10218d00608fSStefano Zampini     if (copy_rhs) {
10228d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10238d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10248d00608fSStefano Zampini     }
10258d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
10263972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10273972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10283972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10298efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10307acc28cbSStefano Zampini     if (ksp) {
10317acc28cbSStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr);
10327acc28cbSStefano Zampini     }
10333308cffdSStefano Zampini   }
10348efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1035b76ba322SStefano Zampini 
1036b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
1037b097fa66SStefano Zampini   if (x && pcbddc->use_exact_dirichlet_trick) {
10388efcfb23SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1039b76ba322SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1040b76ba322SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1041b76ba322SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
10428efcfb23SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10438efcfb23SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1044a0cb1b98SStefano Zampini     if (ksp) {
1045b76ba322SStefano Zampini       ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1046b76ba322SStefano Zampini     }
1047b76ba322SStefano Zampini   }
1048b76ba322SStefano Zampini 
1049b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1050906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1051906d46d4SStefano Zampini 
1052906d46d4SStefano Zampini     /* get change ctx */
1053906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1054906d46d4SStefano Zampini 
1055906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1056906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1057906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1058906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1059906d46d4SStefano Zampini 
1060906d46d4SStefano Zampini     /* change iteration matrix */
1061906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1062906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1063906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1064906d46d4SStefano Zampini 
10658d00608fSStefano Zampini     /* store the original rhs */
10668d00608fSStefano Zampini     if (copy_rhs) {
10678d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10688d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10698d00608fSStefano Zampini     }
10708d00608fSStefano Zampini 
1071906d46d4SStefano Zampini     /* change rhs */
1072906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1073906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
10748d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
107592e3dcfbSStefano Zampini   }
107692e3dcfbSStefano Zampini 
107792e3dcfbSStefano Zampini   /* remove nullspace if present */
10788efcfb23SStefano Zampini   if (ksp && x && pcbddc->NullSpace) {
10798efcfb23SStefano Zampini     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
10808d00608fSStefano Zampini     /* store the original rhs */
10818d00608fSStefano Zampini     if (copy_rhs) {
10828d00608fSStefano Zampini       ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
10838d00608fSStefano Zampini       copy_rhs = PETSC_FALSE;
10848d00608fSStefano Zampini     }
10858d00608fSStefano Zampini     pcbddc->rhs_change = PETSC_TRUE;
1086d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1087b76ba322SStefano Zampini   }
1088534831adSStefano Zampini   PetscFunctionReturn(0);
1089534831adSStefano Zampini }
1090906d46d4SStefano Zampini 
1091534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1092534831adSStefano Zampini #undef __FUNCT__
1093534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1094534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1095534831adSStefano Zampini /*
1096534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1097534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1098534831adSStefano Zampini 
1099534831adSStefano Zampini    Input Parameter:
1100534831adSStefano Zampini +  pc - the preconditioner contex
1101534831adSStefano Zampini 
1102534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1103534831adSStefano Zampini 
1104534831adSStefano Zampini    Notes:
1105534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
1106534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1107534831adSStefano Zampini */
1108534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1109534831adSStefano Zampini {
1110534831adSStefano Zampini   PetscErrorCode ierr;
1111534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1112534831adSStefano Zampini 
1113534831adSStefano Zampini   PetscFunctionBegin;
1114b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1115906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1116906d46d4SStefano Zampini 
1117906d46d4SStefano Zampini     /* get change ctx */
1118906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1119906d46d4SStefano Zampini 
1120906d46d4SStefano Zampini     /* restore iteration matrix */
1121906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1122906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1123906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1124906d46d4SStefano Zampini 
1125906d46d4SStefano Zampini     /* get solution in original basis */
1126906d46d4SStefano Zampini     if (x) {
1127906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1128906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1129906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
11303425bc38SStefano Zampini     }
1131534831adSStefano Zampini   }
1132906d46d4SStefano Zampini 
11333972b0daSStefano Zampini   /* add solution removed in presolve */
11346bcfc461SStefano Zampini   if (x && pcbddc->rhs_change) {
11353425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
11363425bc38SStefano Zampini   }
1137906d46d4SStefano Zampini 
1138fb223d50SStefano Zampini   /* restore rhs to its original state */
11398d00608fSStefano Zampini   if (rhs && pcbddc->rhs_change) {
1140fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1141fb223d50SStefano Zampini   }
11428d00608fSStefano Zampini   pcbddc->rhs_change = PETSC_FALSE;
11438efcfb23SStefano Zampini 
11448efcfb23SStefano Zampini   /* restore ksp guess state */
11458efcfb23SStefano Zampini   if (ksp) {
11468efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
11478efcfb23SStefano Zampini   }
1148534831adSStefano Zampini   PetscFunctionReturn(0);
1149534831adSStefano Zampini }
1150534831adSStefano Zampini /* -------------------------------------------------------------------------- */
115153cdbc3dSStefano Zampini #undef __FUNCT__
115253cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
11530c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
11540c7d97c5SJed Brown /*
11550c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
11560c7d97c5SJed Brown                   by setting data structures and options.
11570c7d97c5SJed Brown 
11580c7d97c5SJed Brown    Input Parameter:
115953cdbc3dSStefano Zampini +  pc - the preconditioner context
11600c7d97c5SJed Brown 
11610c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
11620c7d97c5SJed Brown 
11630c7d97c5SJed Brown    Notes:
11640c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
11650c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
11660c7d97c5SJed Brown */
116753cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
11680c7d97c5SJed Brown {
11690c7d97c5SJed Brown   PetscErrorCode ierr;
11700c7d97c5SJed Brown   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
11715e8657edSStefano Zampini   Mat_IS*        matis;
117208122e43SStefano Zampini   MatNullSpace   nearnullspace;
117308122e43SStefano Zampini   PetscBool      computetopography,computesolvers,computesubschurs;
11748de1fae6SStefano Zampini   PetscBool      computeconstraintsmatrix;
11755e8657edSStefano Zampini   PetscBool      new_nearnullspace_provided,ismatis;
11760c7d97c5SJed Brown 
11770c7d97c5SJed Brown   PetscFunctionBegin;
11785e8657edSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr);
11795e8657edSStefano Zampini   if (!ismatis) {
11805e8657edSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
11815e8657edSStefano Zampini   }
11825e8657edSStefano Zampini   matis = (Mat_IS*)pc->pmat->data;
1183f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
11843b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1185fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1186f4ddd8eeSStefano Zampini   /* split work */
1187674ae819SStefano Zampini   if (pc->setupcalled) {
1188d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1189674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1190674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1191674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1192674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1193674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1194674ae819SStefano Zampini     }
1195674ae819SStefano Zampini   } else {
1196674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1197674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1198674ae819SStefano Zampini   }
1199fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1200fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1201fb180af4SStefano Zampini   }
12028de1fae6SStefano Zampini   computeconstraintsmatrix = PETSC_FALSE;
12035a95e1ceSStefano Zampini   if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) {
12045a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling");
12055a95e1ceSStefano Zampini   }
12065a95e1ceSStefano Zampini   pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling);
12075a95e1ceSStefano Zampini   computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1208f4ddd8eeSStefano Zampini 
1209f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
121070cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
121170cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
121258a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1213f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1214f4ddd8eeSStefano Zampini     }
121558a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1216f4ddd8eeSStefano Zampini   }
1217f4ddd8eeSStefano Zampini 
1218d65f70fdSStefano Zampini   /*ierr = MatIsSymmetric(pc->pmat,0.,&pcbddc->issym);CHKERRQ(ierr);*/
1219d65f70fdSStefano Zampini   { /* this is a temporary workaround since seqbaij matrices does not have support for symmetry checking */
1220d65f70fdSStefano Zampini     PetscBool setsym;
1221d65f70fdSStefano Zampini     ierr = MatIsSymmetricKnown(pc->pmat,&setsym,&pcbddc->issym);CHKERRQ(ierr);
1222d65f70fdSStefano Zampini     if (!setsym) pcbddc->issym = PETSC_FALSE;
1223d65f70fdSStefano Zampini   }
1224d65f70fdSStefano Zampini 
12255e8657edSStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
12265e8657edSStefano Zampini     /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
12275e8657edSStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
12285e8657edSStefano Zampini     ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
12295e8657edSStefano Zampini   } else {
1230b96c3477SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
12315e8657edSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
12325e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
1233674ae819SStefano Zampini   }
1234f4ddd8eeSStefano Zampini 
1235e496cd5dSStefano Zampini   /* workaround for reals */
1236e496cd5dSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
1237e496cd5dSStefano Zampini   ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,pcbddc->issym);CHKERRQ(ierr);
1238e496cd5dSStefano Zampini #endif
1239e496cd5dSStefano Zampini 
12405e8657edSStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
12415e8657edSStefano Zampini   {
12425e8657edSStefano Zampini     Mat temp_mat;
12435e8657edSStefano Zampini 
12445e8657edSStefano Zampini     temp_mat = matis->A;
12455e8657edSStefano Zampini     matis->A = pcbddc->local_mat;
12465e8657edSStefano Zampini     ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr);
12475e8657edSStefano Zampini     pcbddc->local_mat = matis->A;
12485e8657edSStefano Zampini     matis->A = temp_mat;
12495e8657edSStefano Zampini   }
1250684f6988SStefano Zampini 
1251b96c3477SStefano Zampini   /* Analyze interface and setup sub_schurs data */
1252674ae819SStefano Zampini   if (computetopography) {
1253674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
12548de1fae6SStefano Zampini     computeconstraintsmatrix = PETSC_TRUE;
1255674ae819SStefano Zampini   }
1256fb8d54d4SStefano Zampini 
1257b96c3477SStefano Zampini   /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1258684f6988SStefano Zampini   if (computesolvers) {
1259*d5574798SStefano Zampini     PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs;
1260*d5574798SStefano Zampini 
1261*d5574798SStefano Zampini     if (computesubschurs && computetopography) {
126208122e43SStefano Zampini       ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr);
1263b1b3d7a2SStefano Zampini     }
1264*d5574798SStefano Zampini     if (sub_schurs->use_mumps) {
126508122e43SStefano Zampini       ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1266*d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1267*d5574798SStefano Zampini     } else {
1268*d5574798SStefano Zampini       ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1269*d5574798SStefano Zampini       ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr);
1270*d5574798SStefano Zampini     }
127108122e43SStefano Zampini     if (pcbddc->adaptive_selection) {
127208122e43SStefano Zampini       ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr);
12738de1fae6SStefano Zampini       computeconstraintsmatrix = PETSC_TRUE;
1274b7eb3628SStefano Zampini     }
1275b1b3d7a2SStefano Zampini   }
1276684f6988SStefano Zampini 
1277f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1278fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1279f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1280f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1281f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1282f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1283f4ddd8eeSStefano Zampini     } else {
1284f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1285f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1286f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1287165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1288f4ddd8eeSStefano Zampini         PetscInt         i;
1289165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1290165b64e2SStefano Zampini         PetscObjectState state;
1291165b64e2SStefano Zampini         PetscInt         nnsp_size;
1292165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1293f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1294f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1295165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1296f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1297f4ddd8eeSStefano Zampini             break;
1298f4ddd8eeSStefano Zampini           }
1299f4ddd8eeSStefano Zampini         }
1300f4ddd8eeSStefano Zampini       }
1301f4ddd8eeSStefano Zampini     }
1302f4ddd8eeSStefano Zampini   } else {
1303f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1304f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1305f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1306f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1307f4ddd8eeSStefano Zampini     }
1308f4ddd8eeSStefano Zampini   }
1309f4ddd8eeSStefano Zampini 
1310f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1311727cdba6SStefano Zampini   /* reset primal space flags */
1312f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1313727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
13148de1fae6SStefano Zampini   if (computeconstraintsmatrix || new_nearnullspace_provided) {
1315727cdba6SStefano Zampini     /* It also sets the primal space flags */
1316674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1317e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1318f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
13193975b054SStefano Zampini   }
13205e8657edSStefano Zampini 
13213975b054SStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
13225e8657edSStefano Zampini     if (pcbddc->use_change_of_basis) {
13235e8657edSStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
13245e8657edSStefano Zampini 
13255e8657edSStefano Zampini       ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
13265e8657edSStefano Zampini       /* get submatrices */
13275e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
13285e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
13295e8657edSStefano Zampini       ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
13305e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr);
13315e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr);
13325e8657edSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr);
13333975b054SStefano Zampini       /* set flag in pcis to not reuse submatrices during PCISCreate */
13343975b054SStefano Zampini       pcis->reusesubmatrices = PETSC_FALSE;
13355e8657edSStefano Zampini     } else if (!pcbddc->user_ChangeOfBasisMatrix) {
1336b96c3477SStefano Zampini       ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
13375e8657edSStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
13385e8657edSStefano Zampini       pcbddc->local_mat = matis->A;
13395e8657edSStefano Zampini     }
1340b96c3477SStefano Zampini     /* SetUp coarse and local Neumann solvers */
134199cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1342b96c3477SStefano Zampini     /* SetUp Scaling operator */
1343674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
13440c7d97c5SJed Brown   }
134570cf5478SStefano Zampini 
134658a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
134758a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
13482b510759SStefano Zampini   }
13490c7d97c5SJed Brown   PetscFunctionReturn(0);
13500c7d97c5SJed Brown }
13510c7d97c5SJed Brown 
13520c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
13530c7d97c5SJed Brown /*
135450efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
13550c7d97c5SJed Brown 
13560c7d97c5SJed Brown    Input Parameters:
13570c7d97c5SJed Brown .  pc - the preconditioner context
13580c7d97c5SJed Brown .  r - input vector (global)
13590c7d97c5SJed Brown 
13600c7d97c5SJed Brown    Output Parameter:
13610c7d97c5SJed Brown .  z - output vector (global)
13620c7d97c5SJed Brown 
13630c7d97c5SJed Brown    Application Interface Routine: PCApply()
13640c7d97c5SJed Brown  */
13650c7d97c5SJed Brown #undef __FUNCT__
13660c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
136753cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
13680c7d97c5SJed Brown {
13690c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
13700c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1371b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
13720c7d97c5SJed Brown   PetscErrorCode    ierr;
13733b03a366Sstefano_zampini   const PetscScalar one = 1.0;
13743b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
13752617d88aSStefano Zampini   const PetscScalar zero = 0.0;
13760c7d97c5SJed Brown 
13770c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
13780c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
1379b097fa66SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
13800c7d97c5SJed Brown 
13810c7d97c5SJed Brown   PetscFunctionBegin;
138285c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1383b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
13840c7d97c5SJed Brown     /* First Dirichlet solve */
13850c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13860c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13870c7d97c5SJed Brown     /*
13880c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1389b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1390674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
13910c7d97c5SJed Brown     */
1392b097fa66SStefano Zampini     if (n_D) {
1393b097fa66SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
13940c7d97c5SJed Brown       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
13958eeda7d8SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1396b097fa66SStefano Zampini       ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1397b097fa66SStefano Zampini     } else {
1398b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1399b097fa66SStefano Zampini     }
14000c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14010c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1402674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1403b76ba322SStefano Zampini   } else {
1404b097fa66SStefano Zampini     if (pcbddc->switch_static) {
14050bdf917eSStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1406b097fa66SStefano Zampini     }
1407674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1408b76ba322SStefano Zampini   }
1409b76ba322SStefano Zampini 
14102617d88aSStefano Zampini   /* Apply interface preconditioner
14112617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1412dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
14132617d88aSStefano Zampini 
1414674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1415674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
14160c7d97c5SJed Brown 
14173b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
14180c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14190c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1420b097fa66SStefano Zampini   if (n_B) {
14210c7d97c5SJed Brown     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
14228eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1423b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1424b097fa66SStefano Zampini     ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1425b097fa66SStefano Zampini   }
1426df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1427b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1428b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1429b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1430b097fa66SStefano Zampini     } else {
1431b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1432b097fa66SStefano Zampini     }
14330c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14340c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1435b097fa66SStefano Zampini   } else {
1436b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1437b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1438b097fa66SStefano Zampini     } else {
1439b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1440b097fa66SStefano Zampini     }
1441b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1442b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1443b097fa66SStefano Zampini   }
14440c7d97c5SJed Brown   PetscFunctionReturn(0);
14450c7d97c5SJed Brown }
144650efa1b5SStefano Zampini 
144750efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
144850efa1b5SStefano Zampini /*
144950efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
145050efa1b5SStefano Zampini 
145150efa1b5SStefano Zampini    Input Parameters:
145250efa1b5SStefano Zampini .  pc - the preconditioner context
145350efa1b5SStefano Zampini .  r - input vector (global)
145450efa1b5SStefano Zampini 
145550efa1b5SStefano Zampini    Output Parameter:
145650efa1b5SStefano Zampini .  z - output vector (global)
145750efa1b5SStefano Zampini 
145850efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
145950efa1b5SStefano Zampini  */
146050efa1b5SStefano Zampini #undef __FUNCT__
146150efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
146250efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
146350efa1b5SStefano Zampini {
146450efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
146550efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
1466b097fa66SStefano Zampini   PetscInt          n_B = pcis->n_B, n_D = pcis->n - n_B;
146750efa1b5SStefano Zampini   PetscErrorCode    ierr;
146850efa1b5SStefano Zampini   const PetscScalar one = 1.0;
146950efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
147050efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
147150efa1b5SStefano Zampini 
147250efa1b5SStefano Zampini   PetscFunctionBegin;
147350efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1474b097fa66SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
147550efa1b5SStefano Zampini     /* First Dirichlet solve */
147650efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
147750efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
147850efa1b5SStefano Zampini     /*
147950efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
1480b097fa66SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
148150efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
148250efa1b5SStefano Zampini     */
1483b097fa66SStefano Zampini     if (n_D) {
1484b097fa66SStefano Zampini       ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
148550efa1b5SStefano Zampini       ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
148650efa1b5SStefano Zampini       if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
1487b097fa66SStefano Zampini       ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
1488b097fa66SStefano Zampini     } else {
1489b097fa66SStefano Zampini       ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
1490b097fa66SStefano Zampini     }
149150efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
149250efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
149350efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
149450efa1b5SStefano Zampini   } else {
1495b097fa66SStefano Zampini     if (pcbddc->switch_static) {
149650efa1b5SStefano Zampini       ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1497b097fa66SStefano Zampini     }
149850efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
149950efa1b5SStefano Zampini   }
150050efa1b5SStefano Zampini 
150150efa1b5SStefano Zampini   /* Apply interface preconditioner
150250efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1503dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
150450efa1b5SStefano Zampini 
150550efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
150650efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
150750efa1b5SStefano Zampini 
150850efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
150950efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
151050efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1511b097fa66SStefano Zampini   if (n_B) {
151250efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
151350efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1514b097fa66SStefano Zampini   } else if (pcbddc->switch_static) {
1515b097fa66SStefano Zampini     ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr);
1516b097fa66SStefano Zampini   }
1517b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1518b097fa66SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
1519b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1520b097fa66SStefano Zampini       ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr);
1521b097fa66SStefano Zampini     } else {
1522b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr);
1523b097fa66SStefano Zampini     }
152450efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
152550efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1526b097fa66SStefano Zampini   } else {
1527b097fa66SStefano Zampini     if (pcbddc->switch_static) {
1528b097fa66SStefano Zampini       ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr);
1529b097fa66SStefano Zampini     } else {
1530b097fa66SStefano Zampini       ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1531b097fa66SStefano Zampini     }
1532b097fa66SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1533b097fa66SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1534b097fa66SStefano Zampini   }
153550efa1b5SStefano Zampini   PetscFunctionReturn(0);
153650efa1b5SStefano Zampini }
1537da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1538674ae819SStefano Zampini 
1539da1bb401SStefano Zampini #undef __FUNCT__
1540da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1541da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1542da1bb401SStefano Zampini {
1543da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1544da1bb401SStefano Zampini   PetscErrorCode ierr;
1545da1bb401SStefano Zampini 
1546da1bb401SStefano Zampini   PetscFunctionBegin;
1547da1bb401SStefano Zampini   /* free data created by PCIS */
1548da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1549674ae819SStefano Zampini   /* free BDDC custom data  */
1550674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1551674ae819SStefano Zampini   /* destroy objects related to topography */
1552674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1553674ae819SStefano Zampini   /* free allocated graph structure */
1554da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1555b96c3477SStefano Zampini   /* free allocated sub schurs structure */
1556b96c3477SStefano Zampini   ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr);
155734a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1558674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
155934a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1560674ae819SStefano Zampini   /* free solvers stuff */
1561674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
156262a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
156362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
156462a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1565906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1566906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1567906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1568906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1569906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1570906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1571906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1572906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1573906d46d4SStefano Zampini   }
1574906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
15753425bc38SStefano Zampini   /* remove functions */
1576906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1577674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1578bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
15792b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1580b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
15812b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1582bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1583bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
158482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1585bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
158682ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1587bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
158882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1589bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1590785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1591bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
159263602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1593bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1594bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1595bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1596bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1597674ae819SStefano Zampini   /* Free the private data structure */
1598674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1599da1bb401SStefano Zampini   PetscFunctionReturn(0);
1600da1bb401SStefano Zampini }
16013425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
16021e6b0712SBarry Smith 
16033425bc38SStefano Zampini #undef __FUNCT__
16043425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
16053425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
16063425bc38SStefano Zampini {
1607674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
1608c08af4c6SStefano Zampini   Vec            copy_standard_rhs;
16093425bc38SStefano Zampini   PC_IS*         pcis;
16103425bc38SStefano Zampini   PC_BDDC*       pcbddc;
16113425bc38SStefano Zampini   PetscErrorCode ierr;
16120c7d97c5SJed Brown 
16133425bc38SStefano Zampini   PetscFunctionBegin;
16143425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16153425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
16163425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
16173425bc38SStefano Zampini 
1618c08af4c6SStefano Zampini   /*
1619c08af4c6SStefano Zampini      change of basis for physical rhs if needed
1620c08af4c6SStefano Zampini      It also changes the rhs in case of dirichlet boundaries
1621c08af4c6SStefano Zampini      TODO: better management when FETIDP will have its own class
1622c08af4c6SStefano Zampini   */
1623c08af4c6SStefano Zampini   ierr = VecDuplicate(standard_rhs,&copy_standard_rhs);CHKERRQ(ierr);
1624c08af4c6SStefano Zampini   ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr);
1625c08af4c6SStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr);
16263425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
1627c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1628c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1629fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1630fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
1631c08af4c6SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1632c08af4c6SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1633674ae819SStefano Zampini   /* Apply partition of unity */
16343425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1635c08af4c6SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
16368eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
16373425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
16383425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16393425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
16403425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
1641c08af4c6SStefano Zampini     ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr);
1642c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1643c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1644c08af4c6SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
1645c08af4c6SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1646c08af4c6SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16473425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
16483425bc38SStefano Zampini   }
1649c08af4c6SStefano Zampini   ierr = VecDestroy(&copy_standard_rhs);CHKERRQ(ierr);
16503425bc38SStefano Zampini   /* BDDC rhs */
16513425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
16528eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
16533425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
16543425bc38SStefano Zampini   }
16553425bc38SStefano Zampini   /* apply BDDC */
1656dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
16573425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
16583425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
16593425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
16603425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16613425bc38SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16623425bc38SStefano Zampini   PetscFunctionReturn(0);
16633425bc38SStefano Zampini }
16641e6b0712SBarry Smith 
16653425bc38SStefano Zampini #undef __FUNCT__
16663425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
16673425bc38SStefano Zampini /*@
166828509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
16693425bc38SStefano Zampini 
16703425bc38SStefano Zampini    Collective
16713425bc38SStefano Zampini 
16723425bc38SStefano Zampini    Input Parameters:
167328509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
167428509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
16753425bc38SStefano Zampini 
16763425bc38SStefano Zampini    Output Parameters:
167728509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
16783425bc38SStefano Zampini 
16793425bc38SStefano Zampini    Level: developer
16803425bc38SStefano Zampini 
16813425bc38SStefano Zampini    Notes:
16823425bc38SStefano Zampini 
168328509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
16843425bc38SStefano Zampini @*/
16853425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
16863425bc38SStefano Zampini {
1687674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
16883425bc38SStefano Zampini   PetscErrorCode ierr;
16893425bc38SStefano Zampini 
16903425bc38SStefano Zampini   PetscFunctionBegin;
16913425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16923425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
16933425bc38SStefano Zampini   PetscFunctionReturn(0);
16943425bc38SStefano Zampini }
16953425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
16961e6b0712SBarry Smith 
16973425bc38SStefano Zampini #undef __FUNCT__
16983425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
16993425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
17003425bc38SStefano Zampini {
1701674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17023425bc38SStefano Zampini   PC_IS*         pcis;
17033425bc38SStefano Zampini   PC_BDDC*       pcbddc;
17043425bc38SStefano Zampini   PetscErrorCode ierr;
17053425bc38SStefano Zampini 
17063425bc38SStefano Zampini   PetscFunctionBegin;
17073425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
17083425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
17093425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
17103425bc38SStefano Zampini 
17113425bc38SStefano Zampini   /* apply B_delta^T */
17123425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17133425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17143425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
17153425bc38SStefano Zampini   /* compute rhs for BDDC application */
17163425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
17178eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
17183425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17193425bc38SStefano Zampini   }
17203425bc38SStefano Zampini   /* apply BDDC */
1721dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
17223425bc38SStefano Zampini   /* put values into standard global vector */
17233425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17243425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17258eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
17263425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
17273425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
17283425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
17293425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
17303425bc38SStefano Zampini   }
17313425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17323425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17333425bc38SStefano Zampini   /* final change of basis if needed
17343425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
17353308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
17363425bc38SStefano Zampini   PetscFunctionReturn(0);
17373425bc38SStefano Zampini }
17381e6b0712SBarry Smith 
17393425bc38SStefano Zampini #undef __FUNCT__
17403425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
17413425bc38SStefano Zampini /*@
174228509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
17433425bc38SStefano Zampini 
17443425bc38SStefano Zampini    Collective
17453425bc38SStefano Zampini 
17463425bc38SStefano Zampini    Input Parameters:
174728509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
174828509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
17493425bc38SStefano Zampini 
17503425bc38SStefano Zampini    Output Parameters:
175128509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
17523425bc38SStefano Zampini 
17533425bc38SStefano Zampini    Level: developer
17543425bc38SStefano Zampini 
17553425bc38SStefano Zampini    Notes:
17563425bc38SStefano Zampini 
175728509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
17583425bc38SStefano Zampini @*/
17593425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
17603425bc38SStefano Zampini {
1761674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
17623425bc38SStefano Zampini   PetscErrorCode ierr;
17633425bc38SStefano Zampini 
17643425bc38SStefano Zampini   PetscFunctionBegin;
17653425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
17663425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
17673425bc38SStefano Zampini   PetscFunctionReturn(0);
17683425bc38SStefano Zampini }
17693425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
17701e6b0712SBarry Smith 
1771f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1772edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1773f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1774f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1775edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1776f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1777674ae819SStefano Zampini 
17783425bc38SStefano Zampini #undef __FUNCT__
17793425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
17803425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
17813425bc38SStefano Zampini {
1782674ae819SStefano Zampini 
1783674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
17843425bc38SStefano Zampini   Mat            newmat;
1785674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
17863425bc38SStefano Zampini   PC             newpc;
1787ce94432eSBarry Smith   MPI_Comm       comm;
17883425bc38SStefano Zampini   PetscErrorCode ierr;
17893425bc38SStefano Zampini 
17903425bc38SStefano Zampini   PetscFunctionBegin;
1791ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
17923425bc38SStefano Zampini   /* FETIDP linear matrix */
17933425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
17943425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
17953425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
17963425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1797edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
17983425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
17993425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
18003425bc38SStefano Zampini   /* FETIDP preconditioner */
18013425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
18023425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
18033425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
18043425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
18053425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
18063425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1807edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
18083425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
180923ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
18103425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
18113425bc38SStefano Zampini   /* return pointers for objects created */
18123425bc38SStefano Zampini   *fetidp_mat=newmat;
18133425bc38SStefano Zampini   *fetidp_pc=newpc;
18143425bc38SStefano Zampini   PetscFunctionReturn(0);
18153425bc38SStefano Zampini }
18161e6b0712SBarry Smith 
18173425bc38SStefano Zampini #undef __FUNCT__
18183425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
18193425bc38SStefano Zampini /*@
182028509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
18213425bc38SStefano Zampini 
18223425bc38SStefano Zampini    Collective
18233425bc38SStefano Zampini 
18243425bc38SStefano Zampini    Input Parameters:
182528509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
182628509bceSStefano Zampini 
182728509bceSStefano Zampini    Output Parameters:
182828509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
182928509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
183028509bceSStefano Zampini 
183128509bceSStefano Zampini    Options Database Keys:
183228509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
18333425bc38SStefano Zampini 
18343425bc38SStefano Zampini    Level: developer
18353425bc38SStefano Zampini 
18363425bc38SStefano Zampini    Notes:
183728509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
18383425bc38SStefano Zampini 
18393425bc38SStefano Zampini .seealso: PCBDDC
18403425bc38SStefano Zampini @*/
18413425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
18423425bc38SStefano Zampini {
18433425bc38SStefano Zampini   PetscErrorCode ierr;
18443425bc38SStefano Zampini 
18453425bc38SStefano Zampini   PetscFunctionBegin;
18463425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18473425bc38SStefano Zampini   if (pc->setupcalled) {
1848516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1849f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
18503425bc38SStefano Zampini   PetscFunctionReturn(0);
18513425bc38SStefano Zampini }
18520c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1853da1bb401SStefano Zampini /*MC
1854da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
18550c7d97c5SJed Brown 
185628509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
185728509bceSStefano Zampini 
185828509bceSStefano Zampini .vb
185928509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
186028509bceSStefano 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
186128509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
186228509bceSStefano Zampini .ve
186328509bceSStefano Zampini 
186428509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
186528509bceSStefano Zampini 
1866b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
186728509bceSStefano Zampini 
186828509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
186928509bceSStefano Zampini 
1870b6fdb6dfSStefano 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.
1871b6fdb6dfSStefano Zampini 
187228509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
187328509bceSStefano Zampini 
18746a818285SBarry 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()
187528509bceSStefano Zampini 
18766a818285SBarry Smith    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().
187728509bceSStefano Zampini 
1878b6fdb6dfSStefano 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.
187928509bceSStefano Zampini 
188028509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
188128509bceSStefano Zampini 
1882da1bb401SStefano Zampini    Options Database Keys:
188328509bceSStefano Zampini 
188428509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
188528509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1886b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
188728509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
188828509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
188928509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
189028509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
189128509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
189228509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
189328509bceSStefano Zampini 
189428509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
189528509bceSStefano Zampini .vb
189628509bceSStefano Zampini       -pc_bddc_dirichlet_
189728509bceSStefano Zampini       -pc_bddc_neumann_
189828509bceSStefano Zampini       -pc_bddc_coarse_
189928509bceSStefano Zampini .ve
190028509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
190128509bceSStefano Zampini 
190228509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
190328509bceSStefano Zampini .vb
1904312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1905312be037SStefano Zampini       -pc_bddc_neumann_lN_
1906312be037SStefano Zampini       -pc_bddc_coarse_lN_
190728509bceSStefano Zampini .ve
1908312be037SStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N.
1909da1bb401SStefano Zampini 
1910da1bb401SStefano Zampini    Level: intermediate
1911da1bb401SStefano Zampini 
1912b6fdb6dfSStefano Zampini    Developer notes:
1913da1bb401SStefano Zampini 
191428509bceSStefano Zampini    New deluxe scaling operator will be available soon.
1915da1bb401SStefano Zampini 
1916da1bb401SStefano Zampini    Contributed by Stefano Zampini
1917da1bb401SStefano Zampini 
1918da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1919da1bb401SStefano Zampini M*/
1920b2573a8aSBarry Smith 
1921da1bb401SStefano Zampini #undef __FUNCT__
1922da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
19238cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1924da1bb401SStefano Zampini {
1925da1bb401SStefano Zampini   PetscErrorCode      ierr;
1926da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1927da1bb401SStefano Zampini 
1928da1bb401SStefano Zampini   PetscFunctionBegin;
1929da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1930b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1931da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1932da1bb401SStefano Zampini 
1933da1bb401SStefano Zampini   /* create PCIS data structure */
1934da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1935da1bb401SStefano Zampini 
193647d04d0dSStefano Zampini   /* BDDC customization */
193708a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
193847d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
193947d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
194047d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
194147d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
194247d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
194347d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
1944fa434743SStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE;
1945fa434743SStefano Zampini   pcbddc->use_qr_single       = PETSC_FALSE;
194647d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1947b9d89cd5SStefano Zampini   /* private */
1948b9d89cd5SStefano Zampini   pcbddc->issym                      = PETSC_FALSE;
1949727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1950e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1951e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1952e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1953727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1954fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
195568457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
1956f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1957727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1958f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1959f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1960f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1961674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
19620bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
19633972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1964534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1965534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1966534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1967b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
1968906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
1969da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1970da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1971da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1972da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
197315aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
197415aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1975da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1976da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1977da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1978da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1979da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1980da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1981da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1982da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1983da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1984da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1985785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1986785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1987785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
198860d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
198960d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
199063602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
1991da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
199263602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
1993da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
199485c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
199547d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
199647d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
19974fad6a16SStefano Zampini   pcbddc->current_level              = 0;
19982b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1999323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
200074e2c79eSStefano Zampini   pcbddc->redistribute_coarse        = 0;
2001f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
2002323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
2003da1bb401SStefano Zampini 
2004674ae819SStefano Zampini   /* create local graph structure */
2005674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
2006674ae819SStefano Zampini 
2007674ae819SStefano Zampini   /* scaling */
2008674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
200934a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
2010ac632422SStefano Zampini   pcbddc->faster_deluxe         = PETSC_FALSE;
2011b96c3477SStefano Zampini 
2012b96c3477SStefano Zampini   /* create sub schurs structure */
2013b96c3477SStefano Zampini   ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr);
2014b96c3477SStefano Zampini   pcbddc->sub_schurs_rebuild     = PETSC_FALSE;
2015b96c3477SStefano Zampini   pcbddc->sub_schurs_layers      = -1;
2016b96c3477SStefano Zampini   pcbddc->sub_schurs_use_useradj = PETSC_FALSE;
2017b96c3477SStefano Zampini 
2018b96c3477SStefano Zampini   pcbddc->computed_rowadj = PETSC_FALSE;
2019da1bb401SStefano Zampini 
2020b7eb3628SStefano Zampini   /* adaptivity */
2021b7eb3628SStefano Zampini   pcbddc->adaptive_threshold      = -1.0;
20229552c7c7SStefano Zampini   pcbddc->adaptive_invert_Stildas = PETSC_TRUE;
202308122e43SStefano Zampini   pcbddc->adaptive_nmax           = 0;
2024e0a4a720SStefano Zampini   pcbddc->adaptive_nmin           = -1;
2025b7eb3628SStefano Zampini 
2026da1bb401SStefano Zampini   /* function pointers */
2027da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
202893bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
2029da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
2030da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
2031da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
2032da1bb401SStefano Zampini   pc->ops->view                = 0;
2033da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
2034da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
2035da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
2036534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
2037534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
2038da1bb401SStefano Zampini 
2039da1bb401SStefano Zampini   /* composing function */
2040906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
2041674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
2042bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
20432b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
2044b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
20452b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
2046bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
2047bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
204882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2049bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
205082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2051bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
205282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
2053bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
205482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
2055bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
205663602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
2057bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
2058bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
2059bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
2060bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
2061da1bb401SStefano Zampini   PetscFunctionReturn(0);
2062da1bb401SStefano Zampini }
20633425bc38SStefano Zampini 
2064