xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 8efcfb2393bd971d0cec8bd3715dcd19d570e690)
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    - ** GetRid of PCBDDCApplySchur, use MatSchur instead
22b9b85e73SStefano Zampini    - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
23eb97c9d2SStefano Zampini    - add support for computing h,H and related using coordinates?
24c0b83709SStefano Zampini    - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap)
25eb97c9d2SStefano Zampini    - Better management in PCIS code
26eb97c9d2SStefano Zampini    - BDDC with MG framework?
27eb97c9d2SStefano Zampini 
28eb97c9d2SStefano Zampini    FETIDP
29eb97c9d2SStefano Zampini    - Move FETIDP code to its own classes
30eb97c9d2SStefano Zampini 
31eb97c9d2SStefano Zampini    MATIS related operations contained in BDDC code
32eb97c9d2SStefano Zampini    - Provide general case for subassembling
33b9b85e73SStefano Zampini    - *** Preallocation routines in MatISGetMPIAXAIJ
34eb97c9d2SStefano Zampini 
3553cdbc3dSStefano Zampini */
360c7d97c5SJed Brown 
3753cdbc3dSStefano Zampini /* ----------------------------------------------------------------------------------------------------------------------------------------------
380c7d97c5SJed Brown    Implementation of BDDC preconditioner based on:
390c7d97c5SJed Brown    C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
4053cdbc3dSStefano Zampini    ---------------------------------------------------------------------------------------------------------------------------------------------- */
4153cdbc3dSStefano Zampini 
42ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
43ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
443b03a366Sstefano_zampini #include <petscblaslapack.h>
45674ae819SStefano Zampini 
460c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
470c7d97c5SJed Brown #undef __FUNCT__
480c7d97c5SJed Brown #define __FUNCT__ "PCSetFromOptions_BDDC"
490c7d97c5SJed Brown PetscErrorCode PCSetFromOptions_BDDC(PC pc)
500c7d97c5SJed Brown {
510c7d97c5SJed Brown   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
520c7d97c5SJed Brown   PetscErrorCode ierr;
530c7d97c5SJed Brown 
540c7d97c5SJed Brown   PetscFunctionBegin;
550c7d97c5SJed Brown   ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr);
568eeda7d8SStefano Zampini   /* Verbose debugging */
578eeda7d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr);
588eeda7d8SStefano Zampini   /* Primal space cumstomization */
5908a5cf49SStefano 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);
608eeda7d8SStefano 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);
618eeda7d8SStefano 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);
628eeda7d8SStefano 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);
638eeda7d8SStefano Zampini   /* Change of basis */
64b9b85e73SStefano 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);
65b9b85e73SStefano 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);
66674ae819SStefano Zampini   if (!pcbddc->use_change_of_basis) {
67674ae819SStefano Zampini     pcbddc->use_change_on_faces = PETSC_FALSE;
68674ae819SStefano Zampini   }
698eeda7d8SStefano Zampini   /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
708eeda7d8SStefano 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);
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);
751cef56d8SStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_deluxe_threshold","Deluxe subproblems (Schur principal minors) smaller than this value are explicilty computed (-1 computes all)","none",pcbddc->deluxe_threshold,&pcbddc->deluxe_threshold,NULL);CHKERRQ(ierr);
7634a97f8cSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_deluxe_rebuild","Whether or not the interface graph for deluxe has to be rebuilt (i.e. use a standard definition of the interface)","none",pcbddc->deluxe_rebuild,&pcbddc->deluxe_rebuild,NULL);CHKERRQ(ierr);
7734a97f8cSStefano Zampini   ierr = PetscOptionsInt("-pc_bddc_deluxe_layers","Number of dofs' layers for the application of deluxe cheap version (i.e. -1 uses all dofs)","none",pcbddc->deluxe_layers,&pcbddc->deluxe_layers,NULL);CHKERRQ(ierr);
7834a97f8cSStefano Zampini   ierr = PetscOptionsBool("-pc_bddc_deluxe_use_useradj","Whether or not the CSR graph specified by the user should be used for computing layers (default is to use adj of local mat)","none",pcbddc->deluxe_use_useradj,&pcbddc->deluxe_use_useradj,NULL);CHKERRQ(ierr);
790c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
800c7d97c5SJed Brown   PetscFunctionReturn(0);
810c7d97c5SJed Brown }
820c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
83674ae819SStefano Zampini #undef __FUNCT__
84906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
85906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
86b9b85e73SStefano Zampini {
87b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
88b9b85e73SStefano Zampini   PetscErrorCode ierr;
89b9b85e73SStefano Zampini 
90b9b85e73SStefano Zampini   PetscFunctionBegin;
91b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
92b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
93b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
94b9b85e73SStefano Zampini   PetscFunctionReturn(0);
95b9b85e73SStefano Zampini }
96b9b85e73SStefano Zampini #undef __FUNCT__
97906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
98b9b85e73SStefano Zampini /*@
99906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
100b9b85e73SStefano Zampini 
101b9b85e73SStefano Zampini    Collective on PC
102b9b85e73SStefano Zampini 
103b9b85e73SStefano Zampini    Input Parameters:
104b9b85e73SStefano Zampini +  pc - the preconditioning context
105906d46d4SStefano Zampini -  change - the change of basis matrix
106b9b85e73SStefano Zampini 
107b9b85e73SStefano Zampini    Level: intermediate
108b9b85e73SStefano Zampini 
109b9b85e73SStefano Zampini    Notes:
110b9b85e73SStefano Zampini 
111b9b85e73SStefano Zampini .seealso: PCBDDC
112b9b85e73SStefano Zampini @*/
113906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
114b9b85e73SStefano Zampini {
115b9b85e73SStefano Zampini   PetscErrorCode ierr;
116b9b85e73SStefano Zampini 
117b9b85e73SStefano Zampini   PetscFunctionBegin;
118b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
120906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
121906d46d4SStefano Zampini   if (pc->mat) {
122906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
123906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
124906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
125906d46d4SStefano Zampini     if (rows_c != rows) {
126906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
127906d46d4SStefano Zampini     }
128906d46d4SStefano Zampini     if (cols_c != cols) {
129906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
130906d46d4SStefano Zampini     }
131906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
132906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
133906d46d4SStefano Zampini     if (rows_c != rows) {
134906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
135906d46d4SStefano Zampini     }
136906d46d4SStefano Zampini     if (cols_c != cols) {
137906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
138906d46d4SStefano Zampini     }
139906d46d4SStefano Zampini   }
140906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
141b9b85e73SStefano Zampini   PetscFunctionReturn(0);
142b9b85e73SStefano Zampini }
143b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
144b9b85e73SStefano Zampini #undef __FUNCT__
145674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
146674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
147674ae819SStefano Zampini {
148674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
149674ae819SStefano Zampini   PetscErrorCode ierr;
1501e6b0712SBarry Smith 
151674ae819SStefano Zampini   PetscFunctionBegin;
152674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
153674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
154674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
155674ae819SStefano Zampini   PetscFunctionReturn(0);
156674ae819SStefano Zampini }
157674ae819SStefano Zampini #undef __FUNCT__
158674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
159674ae819SStefano Zampini /*@
16028509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
161674ae819SStefano Zampini 
162674ae819SStefano Zampini    Not collective
163674ae819SStefano Zampini 
164674ae819SStefano Zampini    Input Parameters:
165674ae819SStefano Zampini +  pc - the preconditioning context
16628509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
167674ae819SStefano Zampini 
168674ae819SStefano Zampini    Level: intermediate
169674ae819SStefano Zampini 
170674ae819SStefano Zampini    Notes:
171674ae819SStefano Zampini 
172674ae819SStefano Zampini .seealso: PCBDDC
173674ae819SStefano Zampini @*/
174674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
175674ae819SStefano Zampini {
176674ae819SStefano Zampini   PetscErrorCode ierr;
177674ae819SStefano Zampini 
178674ae819SStefano Zampini   PetscFunctionBegin;
179674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
180674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
181674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
182674ae819SStefano Zampini   PetscFunctionReturn(0);
183674ae819SStefano Zampini }
184674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1850c7d97c5SJed Brown #undef __FUNCT__
1864fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1874fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1884fad6a16SStefano Zampini {
1894fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1904fad6a16SStefano Zampini 
1914fad6a16SStefano Zampini   PetscFunctionBegin;
1924fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1934fad6a16SStefano Zampini   PetscFunctionReturn(0);
1944fad6a16SStefano Zampini }
1951e6b0712SBarry Smith 
1964fad6a16SStefano Zampini #undef __FUNCT__
1974fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1984fad6a16SStefano Zampini /*@
19928509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
2004fad6a16SStefano Zampini 
2014fad6a16SStefano Zampini    Logically collective on PC
2024fad6a16SStefano Zampini 
2034fad6a16SStefano Zampini    Input Parameters:
2044fad6a16SStefano Zampini +  pc - the preconditioning context
20528509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2064fad6a16SStefano Zampini 
20728509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2084fad6a16SStefano Zampini 
2094fad6a16SStefano Zampini    Level: intermediate
2104fad6a16SStefano Zampini 
2114fad6a16SStefano Zampini    Notes:
2124fad6a16SStefano Zampini 
2134fad6a16SStefano Zampini .seealso: PCBDDC
2144fad6a16SStefano Zampini @*/
2154fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2164fad6a16SStefano Zampini {
2174fad6a16SStefano Zampini   PetscErrorCode ierr;
2184fad6a16SStefano Zampini 
2194fad6a16SStefano Zampini   PetscFunctionBegin;
2204fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2212b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2224fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2234fad6a16SStefano Zampini   PetscFunctionReturn(0);
2244fad6a16SStefano Zampini }
2252b510759SStefano Zampini 
226b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2272b510759SStefano Zampini #undef __FUNCT__
228b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
229b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
230b8ffe317SStefano Zampini {
231b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
232b8ffe317SStefano Zampini 
233b8ffe317SStefano Zampini   PetscFunctionBegin;
23485c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
235b8ffe317SStefano Zampini   PetscFunctionReturn(0);
236b8ffe317SStefano Zampini }
237b8ffe317SStefano Zampini 
238b8ffe317SStefano Zampini #undef __FUNCT__
239b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
240b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2412b510759SStefano Zampini {
2422b510759SStefano Zampini   PetscErrorCode ierr;
2432b510759SStefano Zampini 
2442b510759SStefano Zampini   PetscFunctionBegin;
2452b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
246b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
247b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2482b510759SStefano Zampini   PetscFunctionReturn(0);
2492b510759SStefano Zampini }
2501e6b0712SBarry Smith 
2514fad6a16SStefano Zampini #undef __FUNCT__
2522b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2532b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2544fad6a16SStefano Zampini {
2554fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2564fad6a16SStefano Zampini 
2574fad6a16SStefano Zampini   PetscFunctionBegin;
2582b510759SStefano Zampini   pcbddc->current_level = level;
2594fad6a16SStefano Zampini   PetscFunctionReturn(0);
2604fad6a16SStefano Zampini }
2611e6b0712SBarry Smith 
2624fad6a16SStefano Zampini #undef __FUNCT__
263b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
264b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
265b8ffe317SStefano Zampini {
266b8ffe317SStefano Zampini   PetscErrorCode ierr;
267b8ffe317SStefano Zampini 
268b8ffe317SStefano Zampini   PetscFunctionBegin;
269b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
270b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
271b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
272b8ffe317SStefano Zampini   PetscFunctionReturn(0);
273b8ffe317SStefano Zampini }
274b8ffe317SStefano Zampini 
275b8ffe317SStefano Zampini #undef __FUNCT__
2762b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2772b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2782b510759SStefano Zampini {
2792b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2802b510759SStefano Zampini 
2812b510759SStefano Zampini   PetscFunctionBegin;
2822b510759SStefano Zampini   pcbddc->max_levels = levels;
2832b510759SStefano Zampini   PetscFunctionReturn(0);
2842b510759SStefano Zampini }
2852b510759SStefano Zampini 
2862b510759SStefano Zampini #undef __FUNCT__
2872b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2884fad6a16SStefano Zampini /*@
28928509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2904fad6a16SStefano Zampini 
2914fad6a16SStefano Zampini    Logically collective on PC
2924fad6a16SStefano Zampini 
2934fad6a16SStefano Zampini    Input Parameters:
2944fad6a16SStefano Zampini +  pc - the preconditioning context
29528509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2964fad6a16SStefano Zampini 
29728509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2984fad6a16SStefano Zampini 
2994fad6a16SStefano Zampini    Level: intermediate
3004fad6a16SStefano Zampini 
3014fad6a16SStefano Zampini    Notes:
3024fad6a16SStefano Zampini 
3034fad6a16SStefano Zampini .seealso: PCBDDC
3044fad6a16SStefano Zampini @*/
3052b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3064fad6a16SStefano Zampini {
3074fad6a16SStefano Zampini   PetscErrorCode ierr;
3084fad6a16SStefano Zampini 
3094fad6a16SStefano Zampini   PetscFunctionBegin;
3104fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3112b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
312312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3132b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3144fad6a16SStefano Zampini   PetscFunctionReturn(0);
3154fad6a16SStefano Zampini }
3164fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3171e6b0712SBarry Smith 
3184fad6a16SStefano Zampini #undef __FUNCT__
3190bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3200bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3210bdf917eSStefano Zampini {
3220bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3230bdf917eSStefano Zampini   PetscErrorCode ierr;
3240bdf917eSStefano Zampini 
3250bdf917eSStefano Zampini   PetscFunctionBegin;
3260bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3270bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3280bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3290bdf917eSStefano Zampini   PetscFunctionReturn(0);
3300bdf917eSStefano Zampini }
3311e6b0712SBarry Smith 
3320bdf917eSStefano Zampini #undef __FUNCT__
3330bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3340bdf917eSStefano Zampini /*@
33528509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3360bdf917eSStefano Zampini 
3370bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3380bdf917eSStefano Zampini 
3390bdf917eSStefano Zampini    Input Parameters:
3400bdf917eSStefano Zampini +  pc - the preconditioning context
34128509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3420bdf917eSStefano Zampini 
3430bdf917eSStefano Zampini    Level: intermediate
3440bdf917eSStefano Zampini 
3450bdf917eSStefano Zampini    Notes:
3460bdf917eSStefano Zampini 
3470bdf917eSStefano Zampini .seealso: PCBDDC
3480bdf917eSStefano Zampini @*/
3490bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3500bdf917eSStefano Zampini {
3510bdf917eSStefano Zampini   PetscErrorCode ierr;
3520bdf917eSStefano Zampini 
3530bdf917eSStefano Zampini   PetscFunctionBegin;
3540bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
355674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
3562b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
3570bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
3580bdf917eSStefano Zampini   PetscFunctionReturn(0);
3590bdf917eSStefano Zampini }
3600bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3611e6b0712SBarry Smith 
3620bdf917eSStefano Zampini #undef __FUNCT__
3633b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3643b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3653b03a366Sstefano_zampini {
3663b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3673b03a366Sstefano_zampini   PetscErrorCode ierr;
3683b03a366Sstefano_zampini 
3693b03a366Sstefano_zampini   PetscFunctionBegin;
370785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
371785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3723b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
37336e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
37436e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
375fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3763b03a366Sstefano_zampini   PetscFunctionReturn(0);
3773b03a366Sstefano_zampini }
3781e6b0712SBarry Smith 
3793b03a366Sstefano_zampini #undef __FUNCT__
3803b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3813b03a366Sstefano_zampini /*@
38228509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3833b03a366Sstefano_zampini 
384785d1243SStefano Zampini    Collective
3853b03a366Sstefano_zampini 
3863b03a366Sstefano_zampini    Input Parameters:
3873b03a366Sstefano_zampini +  pc - the preconditioning context
388785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3893b03a366Sstefano_zampini 
3903b03a366Sstefano_zampini    Level: intermediate
3913b03a366Sstefano_zampini 
392785d1243SStefano Zampini    Notes: Any process can list any global node
3933b03a366Sstefano_zampini 
3943b03a366Sstefano_zampini .seealso: PCBDDC
3953b03a366Sstefano_zampini @*/
3963b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3973b03a366Sstefano_zampini {
3983b03a366Sstefano_zampini   PetscErrorCode ierr;
3993b03a366Sstefano_zampini 
4003b03a366Sstefano_zampini   PetscFunctionBegin;
4013b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
402674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
403785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4043b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4053b03a366Sstefano_zampini   PetscFunctionReturn(0);
4063b03a366Sstefano_zampini }
4073b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4081e6b0712SBarry Smith 
4093b03a366Sstefano_zampini #undef __FUNCT__
41082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
41182ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4120c7d97c5SJed Brown {
4130c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4140c7d97c5SJed Brown   PetscErrorCode ierr;
4150c7d97c5SJed Brown 
4160c7d97c5SJed Brown   PetscFunctionBegin;
417785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
418785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4190c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4200c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
421785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4227d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4230c7d97c5SJed Brown   PetscFunctionReturn(0);
4240c7d97c5SJed Brown }
4250c7d97c5SJed Brown 
4260c7d97c5SJed Brown #undef __FUNCT__
42782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4289c0446d6SStefano Zampini /*@
42982ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4309c0446d6SStefano Zampini 
431785d1243SStefano Zampini    Collective
43253cdbc3dSStefano Zampini 
43353cdbc3dSStefano Zampini    Input Parameters:
43453cdbc3dSStefano Zampini +  pc - the preconditioning context
43582ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
43653cdbc3dSStefano Zampini 
43753cdbc3dSStefano Zampini    Level: intermediate
43853cdbc3dSStefano Zampini 
4399c0446d6SStefano Zampini    Notes:
44053cdbc3dSStefano Zampini 
44153cdbc3dSStefano Zampini .seealso: PCBDDC
44253cdbc3dSStefano Zampini @*/
44382ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4440c7d97c5SJed Brown {
4450c7d97c5SJed Brown   PetscErrorCode ierr;
4460c7d97c5SJed Brown 
4470c7d97c5SJed Brown   PetscFunctionBegin;
4480c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4490c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
45082ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
45182ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4520c7d97c5SJed Brown   PetscFunctionReturn(0);
4530c7d97c5SJed Brown }
4540c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4550c7d97c5SJed Brown 
4560c7d97c5SJed Brown #undef __FUNCT__
4570c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
45853cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4590c7d97c5SJed Brown {
4600c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
46153cdbc3dSStefano Zampini   PetscErrorCode ierr;
4620c7d97c5SJed Brown 
4630c7d97c5SJed Brown   PetscFunctionBegin;
464785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
465785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
46653cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
46736e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
46836e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
469fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4700c7d97c5SJed Brown   PetscFunctionReturn(0);
4710c7d97c5SJed Brown }
4721e6b0712SBarry Smith 
4730c7d97c5SJed Brown #undef __FUNCT__
4740c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
47557527edcSJed Brown /*@
47628509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
47757527edcSJed Brown 
478785d1243SStefano Zampini    Collective
47957527edcSJed Brown 
48057527edcSJed Brown    Input Parameters:
48157527edcSJed Brown +  pc - the preconditioning context
482785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
48357527edcSJed Brown 
48457527edcSJed Brown    Level: intermediate
48557527edcSJed Brown 
486785d1243SStefano Zampini    Notes: Any process can list any global node
48757527edcSJed Brown 
48857527edcSJed Brown .seealso: PCBDDC
48957527edcSJed Brown @*/
49053cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4910c7d97c5SJed Brown {
4920c7d97c5SJed Brown   PetscErrorCode ierr;
4930c7d97c5SJed Brown 
4940c7d97c5SJed Brown   PetscFunctionBegin;
4950c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
496674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
497785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
49853cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
49953cdbc3dSStefano Zampini   PetscFunctionReturn(0);
50053cdbc3dSStefano Zampini }
50153cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
5021e6b0712SBarry Smith 
50353cdbc3dSStefano Zampini #undef __FUNCT__
50482ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
50582ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5060c7d97c5SJed Brown {
5070c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5080c7d97c5SJed Brown   PetscErrorCode ierr;
5090c7d97c5SJed Brown 
5100c7d97c5SJed Brown   PetscFunctionBegin;
511785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
512785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5130c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5140c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
515785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5167d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5170c7d97c5SJed Brown   PetscFunctionReturn(0);
5180c7d97c5SJed Brown }
5190c7d97c5SJed Brown 
5200c7d97c5SJed Brown #undef __FUNCT__
52182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5220c7d97c5SJed Brown /*@
52382ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5240c7d97c5SJed Brown 
525785d1243SStefano Zampini    Collective
5260c7d97c5SJed Brown 
5270c7d97c5SJed Brown    Input Parameters:
5280c7d97c5SJed Brown +  pc - the preconditioning context
52982ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5300c7d97c5SJed Brown 
5310c7d97c5SJed Brown    Level: intermediate
5320c7d97c5SJed Brown 
5330c7d97c5SJed Brown    Notes:
5340c7d97c5SJed Brown 
5350c7d97c5SJed Brown .seealso: PCBDDC
5360c7d97c5SJed Brown @*/
53782ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5380c7d97c5SJed Brown {
5390c7d97c5SJed Brown   PetscErrorCode ierr;
5400c7d97c5SJed Brown 
5410c7d97c5SJed Brown   PetscFunctionBegin;
5420c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5430c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
54482ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54582ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
54653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
54753cdbc3dSStefano Zampini }
54853cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
54953cdbc3dSStefano Zampini 
55053cdbc3dSStefano Zampini #undef __FUNCT__
551da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
552da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
553da1bb401SStefano Zampini {
554da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
555da1bb401SStefano Zampini 
556da1bb401SStefano Zampini   PetscFunctionBegin;
557da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
558da1bb401SStefano Zampini   PetscFunctionReturn(0);
559da1bb401SStefano Zampini }
5601e6b0712SBarry Smith 
561da1bb401SStefano Zampini #undef __FUNCT__
562da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
563da1bb401SStefano Zampini /*@
564785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
565da1bb401SStefano Zampini 
566785d1243SStefano Zampini    Collective
567785d1243SStefano Zampini 
568785d1243SStefano Zampini    Input Parameters:
569785d1243SStefano Zampini .  pc - the preconditioning context
570785d1243SStefano Zampini 
571785d1243SStefano Zampini    Output Parameters:
572785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
573785d1243SStefano Zampini 
574785d1243SStefano Zampini    Level: intermediate
575785d1243SStefano Zampini 
576785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
577785d1243SStefano Zampini 
578785d1243SStefano Zampini .seealso: PCBDDC
579785d1243SStefano Zampini @*/
580785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
581785d1243SStefano Zampini {
582785d1243SStefano Zampini   PetscErrorCode ierr;
583785d1243SStefano Zampini 
584785d1243SStefano Zampini   PetscFunctionBegin;
585785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
586785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
587785d1243SStefano Zampini   PetscFunctionReturn(0);
588785d1243SStefano Zampini }
589785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
590785d1243SStefano Zampini 
591785d1243SStefano Zampini #undef __FUNCT__
592785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
593785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
594785d1243SStefano Zampini {
595785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
596785d1243SStefano Zampini 
597785d1243SStefano Zampini   PetscFunctionBegin;
598785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
599785d1243SStefano Zampini   PetscFunctionReturn(0);
600785d1243SStefano Zampini }
601785d1243SStefano Zampini 
602785d1243SStefano Zampini #undef __FUNCT__
60382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
604da1bb401SStefano Zampini /*@
60582ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
606da1bb401SStefano Zampini 
607785d1243SStefano Zampini    Collective
608da1bb401SStefano Zampini 
609da1bb401SStefano Zampini    Input Parameters:
61028509bceSStefano Zampini .  pc - the preconditioning context
611da1bb401SStefano Zampini 
612da1bb401SStefano Zampini    Output Parameters:
61328509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
614da1bb401SStefano Zampini 
615da1bb401SStefano Zampini    Level: intermediate
616da1bb401SStefano Zampini 
617da1bb401SStefano Zampini    Notes:
618da1bb401SStefano Zampini 
619da1bb401SStefano Zampini .seealso: PCBDDC
620da1bb401SStefano Zampini @*/
62182ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
622da1bb401SStefano Zampini {
623da1bb401SStefano Zampini   PetscErrorCode ierr;
624da1bb401SStefano Zampini 
625da1bb401SStefano Zampini   PetscFunctionBegin;
626da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
62782ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
628da1bb401SStefano Zampini   PetscFunctionReturn(0);
629da1bb401SStefano Zampini }
630da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6311e6b0712SBarry Smith 
632da1bb401SStefano Zampini #undef __FUNCT__
63353cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
63453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
63553cdbc3dSStefano Zampini {
63653cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
63753cdbc3dSStefano Zampini 
63853cdbc3dSStefano Zampini   PetscFunctionBegin;
63953cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
64053cdbc3dSStefano Zampini   PetscFunctionReturn(0);
64153cdbc3dSStefano Zampini }
6421e6b0712SBarry Smith 
64353cdbc3dSStefano Zampini #undef __FUNCT__
64453cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
64553cdbc3dSStefano Zampini /*@
646785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
64753cdbc3dSStefano Zampini 
648785d1243SStefano Zampini    Collective
649785d1243SStefano Zampini 
650785d1243SStefano Zampini    Input Parameters:
651785d1243SStefano Zampini .  pc - the preconditioning context
652785d1243SStefano Zampini 
653785d1243SStefano Zampini    Output Parameters:
654785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
655785d1243SStefano Zampini 
656785d1243SStefano Zampini    Level: intermediate
657785d1243SStefano Zampini 
658785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
659785d1243SStefano Zampini 
660785d1243SStefano Zampini .seealso: PCBDDC
661785d1243SStefano Zampini @*/
662785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
663785d1243SStefano Zampini {
664785d1243SStefano Zampini   PetscErrorCode ierr;
665785d1243SStefano Zampini 
666785d1243SStefano Zampini   PetscFunctionBegin;
667785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
668785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
669785d1243SStefano Zampini   PetscFunctionReturn(0);
670785d1243SStefano Zampini }
671785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
672785d1243SStefano Zampini 
673785d1243SStefano Zampini #undef __FUNCT__
674785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
675785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
676785d1243SStefano Zampini {
677785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
678785d1243SStefano Zampini 
679785d1243SStefano Zampini   PetscFunctionBegin;
680785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
681785d1243SStefano Zampini   PetscFunctionReturn(0);
682785d1243SStefano Zampini }
683785d1243SStefano Zampini 
684785d1243SStefano Zampini #undef __FUNCT__
68582ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
68653cdbc3dSStefano Zampini /*@
68782ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
68853cdbc3dSStefano Zampini 
689785d1243SStefano Zampini    Collective
69053cdbc3dSStefano Zampini 
69153cdbc3dSStefano Zampini    Input Parameters:
69228509bceSStefano Zampini .  pc - the preconditioning context
69353cdbc3dSStefano Zampini 
69453cdbc3dSStefano Zampini    Output Parameters:
69528509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
69653cdbc3dSStefano Zampini 
69753cdbc3dSStefano Zampini    Level: intermediate
69853cdbc3dSStefano Zampini 
69953cdbc3dSStefano Zampini    Notes:
70053cdbc3dSStefano Zampini 
70153cdbc3dSStefano Zampini .seealso: PCBDDC
70253cdbc3dSStefano Zampini @*/
70382ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
70453cdbc3dSStefano Zampini {
70553cdbc3dSStefano Zampini   PetscErrorCode ierr;
70653cdbc3dSStefano Zampini 
70753cdbc3dSStefano Zampini   PetscFunctionBegin;
70853cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
70982ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7100c7d97c5SJed Brown   PetscFunctionReturn(0);
7110c7d97c5SJed Brown }
71236e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7131e6b0712SBarry Smith 
71436e030ebSStefano Zampini #undef __FUNCT__
715da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7161a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
71736e030ebSStefano Zampini {
71836e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
719da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
720da1bb401SStefano Zampini   PetscErrorCode ierr;
72136e030ebSStefano Zampini 
72236e030ebSStefano Zampini   PetscFunctionBegin;
723674ae819SStefano Zampini   /* free old CSR */
724674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
725fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
726674ae819SStefano Zampini   /* get CSR into graph structure */
727da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
728785e854fSJed Brown     ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr);
729785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
730674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
731674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
732da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7331a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7341a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
735674ae819SStefano Zampini   }
736575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
73736e030ebSStefano Zampini   PetscFunctionReturn(0);
73836e030ebSStefano Zampini }
7391e6b0712SBarry Smith 
74036e030ebSStefano Zampini #undef __FUNCT__
741da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
74236e030ebSStefano Zampini /*@
74328509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
74436e030ebSStefano Zampini 
74536e030ebSStefano Zampini    Not collective
74636e030ebSStefano Zampini 
74736e030ebSStefano Zampini    Input Parameters:
74836e030ebSStefano Zampini +  pc - the preconditioning context
74928509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
75028509bceSStefano Zampini .  xadj, adjncy - the CSR graph
75128509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
75236e030ebSStefano Zampini 
75336e030ebSStefano Zampini    Level: intermediate
75436e030ebSStefano Zampini 
75536e030ebSStefano Zampini    Notes:
75636e030ebSStefano Zampini 
75728509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
75836e030ebSStefano Zampini @*/
7591a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
76036e030ebSStefano Zampini {
761575ad6abSStefano Zampini   void (*f)(void) = 0;
76236e030ebSStefano Zampini   PetscErrorCode ierr;
76336e030ebSStefano Zampini 
76436e030ebSStefano Zampini   PetscFunctionBegin;
76536e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
766674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
767d7de1dd9SStefano Zampini   PetscValidIntPointer(adjncy,4);
768674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
769674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
770da1bb401SStefano Zampini   }
77136e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
772575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
773575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
774575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
775575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
776575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
77736e030ebSStefano Zampini   }
77836e030ebSStefano Zampini   PetscFunctionReturn(0);
77936e030ebSStefano Zampini }
7809c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7811e6b0712SBarry Smith 
7829c0446d6SStefano Zampini #undef __FUNCT__
78363602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
78463602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
78563602bcaSStefano Zampini {
78663602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
78763602bcaSStefano Zampini   PetscInt i;
78863602bcaSStefano Zampini   PetscErrorCode ierr;
78963602bcaSStefano Zampini 
79063602bcaSStefano Zampini   PetscFunctionBegin;
79163602bcaSStefano Zampini   /* Destroy ISes if they were already set */
79263602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
79363602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
79463602bcaSStefano Zampini   }
79563602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
79663602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
79763602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
79863602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
79963602bcaSStefano Zampini   }
80063602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
80163602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
80263602bcaSStefano Zampini   /* allocate space then set */
80363602bcaSStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
80463602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
80563602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
80663602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
80763602bcaSStefano Zampini   }
80863602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
80963602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8101a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
81163602bcaSStefano Zampini   PetscFunctionReturn(0);
81263602bcaSStefano Zampini }
81363602bcaSStefano Zampini 
81463602bcaSStefano Zampini #undef __FUNCT__
81563602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
81663602bcaSStefano Zampini /*@
81763602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
81863602bcaSStefano Zampini 
81963602bcaSStefano Zampini    Collective
82063602bcaSStefano Zampini 
82163602bcaSStefano Zampini    Input Parameters:
82263602bcaSStefano Zampini +  pc - the preconditioning context
82363602bcaSStefano Zampini -  n_is - number of index sets defining the fields
82463602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in local ordering
82563602bcaSStefano Zampini 
82663602bcaSStefano Zampini    Level: intermediate
82763602bcaSStefano Zampini 
82863602bcaSStefano 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.
82963602bcaSStefano Zampini 
83063602bcaSStefano Zampini .seealso: PCBDDC
83163602bcaSStefano Zampini @*/
83263602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
83363602bcaSStefano Zampini {
83463602bcaSStefano Zampini   PetscInt       i;
83563602bcaSStefano Zampini   PetscErrorCode ierr;
83663602bcaSStefano Zampini 
83763602bcaSStefano Zampini   PetscFunctionBegin;
83863602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
83963602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
84063602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
84163602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
84263602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
84363602bcaSStefano Zampini   }
84463602bcaSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
84563602bcaSStefano Zampini   PetscFunctionReturn(0);
84663602bcaSStefano Zampini }
84763602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
84863602bcaSStefano Zampini 
84963602bcaSStefano Zampini #undef __FUNCT__
8509c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8519c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8529c0446d6SStefano Zampini {
8539c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8549c0446d6SStefano Zampini   PetscInt i;
8559c0446d6SStefano Zampini   PetscErrorCode ierr;
8569c0446d6SStefano Zampini 
8579c0446d6SStefano Zampini   PetscFunctionBegin;
858da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8599c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8609c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8619c0446d6SStefano Zampini   }
862d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
86363602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
86463602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
86563602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
86663602bcaSStefano Zampini   }
86763602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
86863602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
869da1bb401SStefano Zampini   /* allocate space then set */
870785e854fSJed Brown   ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
8719c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
872da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
873da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8749c0446d6SStefano Zampini   }
8759c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
87663602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8771a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8789c0446d6SStefano Zampini   PetscFunctionReturn(0);
8799c0446d6SStefano Zampini }
8801e6b0712SBarry Smith 
8819c0446d6SStefano Zampini #undef __FUNCT__
8829c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8839c0446d6SStefano Zampini /*@
88463602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8859c0446d6SStefano Zampini 
88663602bcaSStefano Zampini    Collective
8879c0446d6SStefano Zampini 
8889c0446d6SStefano Zampini    Input Parameters:
8899c0446d6SStefano Zampini +  pc - the preconditioning context
89028509bceSStefano Zampini -  n_is - number of index sets defining the fields
89163602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in global ordering
8929c0446d6SStefano Zampini 
8939c0446d6SStefano Zampini    Level: intermediate
8949c0446d6SStefano Zampini 
89563602bcaSStefano Zampini    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
8969c0446d6SStefano Zampini 
8979c0446d6SStefano Zampini .seealso: PCBDDC
8989c0446d6SStefano Zampini @*/
8999c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
9009c0446d6SStefano Zampini {
9012b510759SStefano Zampini   PetscInt       i;
9029c0446d6SStefano Zampini   PetscErrorCode ierr;
9039c0446d6SStefano Zampini 
9049c0446d6SStefano Zampini   PetscFunctionBegin;
9059c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
90663602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9072b510759SStefano Zampini   for (i=0;i<n_is;i++) {
90863602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
90963602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9102b510759SStefano Zampini   }
9119c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9129c0446d6SStefano Zampini   PetscFunctionReturn(0);
9139c0446d6SStefano Zampini }
914906d46d4SStefano Zampini 
915da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
916534831adSStefano Zampini #undef __FUNCT__
917534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
918534831adSStefano Zampini /* -------------------------------------------------------------------------- */
919534831adSStefano Zampini /*
920534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
921534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9229c0446d6SStefano Zampini 
923534831adSStefano Zampini    Input Parameter:
924534831adSStefano Zampini +  pc - the preconditioner contex
925534831adSStefano Zampini 
926534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
927534831adSStefano Zampini 
928534831adSStefano Zampini    Notes:
929534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
930534831adSStefano Zampini    the user, but instead is called by KSPSolve().
931534831adSStefano Zampini */
932534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
933534831adSStefano Zampini {
934534831adSStefano Zampini   PetscErrorCode ierr;
935534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
936534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9373972b0daSStefano Zampini   IS             dirIS;
9383972b0daSStefano Zampini   Vec            used_vec;
939*8efcfb23SStefano Zampini   PetscBool      remove_rhs = PETSC_FALSE;
940534831adSStefano Zampini 
941534831adSStefano Zampini   PetscFunctionBegin;
94285c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
94385c4d303SStefano Zampini   if (ksp) {
94485c4d303SStefano Zampini     PetscBool iscg;
94585c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
94685c4d303SStefano Zampini     if (!iscg) {
94785c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
94885c4d303SStefano Zampini     }
94985c4d303SStefano Zampini   }
95085c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
95162a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
95262a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
95362a6ff1dSStefano Zampini   }
95462a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
95562a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
95662a6ff1dSStefano Zampini   }
9573972b0daSStefano Zampini   if (x) {
9583972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9593972b0daSStefano Zampini     used_vec = x;
9603972b0daSStefano Zampini   } else {
9613972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9623972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9633972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9643972b0daSStefano Zampini   }
965*8efcfb23SStefano Zampini 
966*8efcfb23SStefano Zampini   /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
9673972b0daSStefano Zampini   if (ksp) {
968*8efcfb23SStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
969*8efcfb23SStefano Zampini     if (!pcbddc->ksp_guess_nonzero) {
9703972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9713972b0daSStefano Zampini     }
9723972b0daSStefano Zampini   }
9733308cffdSStefano Zampini 
9743972b0daSStefano Zampini   /* store the original rhs */
9753972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
9763972b0daSStefano Zampini 
9773972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
978785d1243SStefano Zampini   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
97982ba6b80SStefano Zampini   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
98082ba6b80SStefano Zampini   if (rhs && dirIS) {
981906d46d4SStefano Zampini     Mat_IS      *matis = (Mat_IS*)pc->pmat->data;
982785d1243SStefano Zampini     PetscInt    dirsize,i,*is_indices;
983785d1243SStefano Zampini     PetscScalar *array_x,*array_diagonal;
984785d1243SStefano Zampini 
9853972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
9863972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
9873972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9883972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9893972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9903972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
99182ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
9923972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9933972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9943972b0daSStefano Zampini     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9952fa5cd67SKarl Rupp     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
9963972b0daSStefano Zampini     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9973972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9983972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9993972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10003972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1001*8efcfb23SStefano Zampini     remove_rhs = PETSC_TRUE;
1002*8efcfb23SStefano Zampini   }
1003b76ba322SStefano Zampini 
1004*8efcfb23SStefano Zampini   /* remove the computed solution or the initial guess from the rhs */
1005*8efcfb23SStefano Zampini   if (remove_rhs || (ksp && pcbddc->ksp_guess_nonzero) ) {
10063972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10073972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10083972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
1009*8efcfb23SStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10103308cffdSStefano Zampini   }
1011*8efcfb23SStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1012b76ba322SStefano Zampini 
1013b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
10143972b0daSStefano Zampini   if (x) {
1015*8efcfb23SStefano Zampini     ierr = VecSet(x,0.0);CHKERRQ(ierr);
101685c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
1017b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1018b76ba322SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1019b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1020*8efcfb23SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1021*8efcfb23SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1022*8efcfb23SStefano Zampini       if (ksp && !pcbddc->ksp_guess_nonzero) {
1023b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1024b76ba322SStefano Zampini       }
1025b76ba322SStefano Zampini     }
10263972b0daSStefano Zampini   }
1027b76ba322SStefano Zampini 
1028b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1029906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1030906d46d4SStefano Zampini 
1031906d46d4SStefano Zampini     /* get change ctx */
1032906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1033906d46d4SStefano Zampini 
1034906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1035906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1036906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1037906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1038906d46d4SStefano Zampini 
1039906d46d4SStefano Zampini     /* change iteration matrix */
1040906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1041906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1042906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1043906d46d4SStefano Zampini 
1044906d46d4SStefano Zampini     /* change rhs */
1045906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1046906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
104792e3dcfbSStefano Zampini   }
104892e3dcfbSStefano Zampini 
104992e3dcfbSStefano Zampini   /* remove nullspace if present */
1050*8efcfb23SStefano Zampini   if (ksp && x && pcbddc->NullSpace) {
1051*8efcfb23SStefano Zampini     ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr);
1052d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1053b76ba322SStefano Zampini   }
1054534831adSStefano Zampini   PetscFunctionReturn(0);
1055534831adSStefano Zampini }
1056906d46d4SStefano Zampini 
1057534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1058534831adSStefano Zampini #undef __FUNCT__
1059534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1060534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1061534831adSStefano Zampini /*
1062534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1063534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1064534831adSStefano Zampini 
1065534831adSStefano Zampini    Input Parameter:
1066534831adSStefano Zampini +  pc - the preconditioner contex
1067534831adSStefano Zampini 
1068534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1069534831adSStefano Zampini 
1070534831adSStefano Zampini    Notes:
1071534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
1072534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1073534831adSStefano Zampini */
1074534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1075534831adSStefano Zampini {
1076534831adSStefano Zampini   PetscErrorCode ierr;
1077534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1078534831adSStefano Zampini 
1079534831adSStefano Zampini   PetscFunctionBegin;
1080b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1081906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1082906d46d4SStefano Zampini 
1083906d46d4SStefano Zampini     /* get change ctx */
1084906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1085906d46d4SStefano Zampini 
1086906d46d4SStefano Zampini     /* restore iteration matrix */
1087906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1088906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1089906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1090906d46d4SStefano Zampini 
1091906d46d4SStefano Zampini     /* get solution in original basis */
1092906d46d4SStefano Zampini     if (x) {
1093906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1094906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1095906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
10963425bc38SStefano Zampini     }
1097534831adSStefano Zampini   }
1098906d46d4SStefano Zampini 
10993972b0daSStefano Zampini   /* add solution removed in presolve */
11003425bc38SStefano Zampini   if (x) {
11013425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
11023425bc38SStefano Zampini   }
1103906d46d4SStefano Zampini 
1104fb223d50SStefano Zampini   /* restore rhs to its original state */
1105fb223d50SStefano Zampini   if (rhs) {
1106fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1107fb223d50SStefano Zampini   }
1108*8efcfb23SStefano Zampini 
1109*8efcfb23SStefano Zampini   /* restore ksp guess state */
1110*8efcfb23SStefano Zampini   if (ksp) {
1111*8efcfb23SStefano Zampini     ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr);
1112*8efcfb23SStefano Zampini   }
1113534831adSStefano Zampini   PetscFunctionReturn(0);
1114534831adSStefano Zampini }
1115534831adSStefano Zampini /* -------------------------------------------------------------------------- */
111653cdbc3dSStefano Zampini #undef __FUNCT__
111753cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
11180c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
11190c7d97c5SJed Brown /*
11200c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
11210c7d97c5SJed Brown                   by setting data structures and options.
11220c7d97c5SJed Brown 
11230c7d97c5SJed Brown    Input Parameter:
112453cdbc3dSStefano Zampini +  pc - the preconditioner context
11250c7d97c5SJed Brown 
11260c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
11270c7d97c5SJed Brown 
11280c7d97c5SJed Brown    Notes:
11290c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
11300c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
11310c7d97c5SJed Brown */
113253cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
11330c7d97c5SJed Brown {
11340c7d97c5SJed Brown   PetscErrorCode   ierr;
11350c7d97c5SJed Brown   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1136f4ddd8eeSStefano Zampini   MatNullSpace     nearnullspace;
1137674ae819SStefano Zampini   PetscBool        computeis,computetopography,computesolvers;
1138165b64e2SStefano Zampini   PetscBool        new_nearnullspace_provided;
11390c7d97c5SJed Brown 
11400c7d97c5SJed Brown   PetscFunctionBegin;
1141f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
11423b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1143fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1144f4ddd8eeSStefano Zampini   /* split work */
1145674ae819SStefano Zampini   if (pc->setupcalled) {
1146674ae819SStefano Zampini     computeis = PETSC_FALSE;
1147d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1148674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1149674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1150674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1151674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1152674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1153674ae819SStefano Zampini     }
1154674ae819SStefano Zampini   } else {
1155674ae819SStefano Zampini     computeis = PETSC_TRUE;
1156674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1157674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1158674ae819SStefano Zampini   }
1159fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1160fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1161fb180af4SStefano Zampini   }
1162f4ddd8eeSStefano Zampini 
1163f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
116470cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
116570cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
116658a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1167f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1168f4ddd8eeSStefano Zampini     }
116958a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1170f4ddd8eeSStefano Zampini   }
1171f4ddd8eeSStefano Zampini 
1172fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
1173674ae819SStefano Zampini   if (computeis) {
1174fcd409f5SStefano Zampini     /* HACK INTO PCIS */
1175fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
1176fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
1177674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
117839e2fb2aSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1179674ae819SStefano Zampini   }
1180f4ddd8eeSStefano Zampini 
1181b9b85e73SStefano Zampini   /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1182906d46d4SStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
1183b9b85e73SStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
1184b9b85e73SStefano Zampini   }
1185906d46d4SStefano Zampini 
1186f4ddd8eeSStefano Zampini   /* Analyze interface */
1187674ae819SStefano Zampini   if (computetopography) {
1188674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
118934a97f8cSStefano Zampini     /* Schurs on subsets should be reset */
119034a97f8cSStefano Zampini     if (pcbddc->deluxe_ctx) {
119134a97f8cSStefano Zampini       ierr = PCBDDCSubSchursReset(pcbddc->deluxe_ctx->sub_schurs);CHKERRQ(ierr);
119234a97f8cSStefano Zampini     }
1193674ae819SStefano Zampini   }
1194fb8d54d4SStefano Zampini 
1195f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1196fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1197f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1198f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1199f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1200f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1201f4ddd8eeSStefano Zampini     } else {
1202f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1203f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1204f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1205165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1206f4ddd8eeSStefano Zampini         PetscInt         i;
1207165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1208165b64e2SStefano Zampini         PetscObjectState state;
1209165b64e2SStefano Zampini         PetscInt         nnsp_size;
1210165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1211f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1212f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1213165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1214f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1215f4ddd8eeSStefano Zampini             break;
1216f4ddd8eeSStefano Zampini           }
1217f4ddd8eeSStefano Zampini         }
1218f4ddd8eeSStefano Zampini       }
1219f4ddd8eeSStefano Zampini     }
1220f4ddd8eeSStefano Zampini   } else {
1221f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1222f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1223f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1224f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1225f4ddd8eeSStefano Zampini     }
1226f4ddd8eeSStefano Zampini   }
1227f4ddd8eeSStefano Zampini 
1228f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1229727cdba6SStefano Zampini   /* reset primal space flags */
1230f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1231727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
1232fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
1233727cdba6SStefano Zampini     /* It also sets the primal space flags */
1234674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1235e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1236f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
123734a97f8cSStefano Zampini     /* Schurs on subsets should be reset */
123834a97f8cSStefano Zampini     if (pcbddc->deluxe_ctx) {
123934a97f8cSStefano Zampini       ierr = PCBDDCSubSchursReset(pcbddc->deluxe_ctx->sub_schurs);CHKERRQ(ierr);
124034a97f8cSStefano Zampini     }
1241674ae819SStefano Zampini   }
1242fb8d54d4SStefano Zampini 
1243f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
124499cc7994SStefano Zampini     /* Create coarse and local stuffs */
124599cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
124634a97f8cSStefano Zampini     /* Create scaling operators */
1247674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
12480c7d97c5SJed Brown   }
124970cf5478SStefano Zampini 
125058a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
125158a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
12522b510759SStefano Zampini   }
12530c7d97c5SJed Brown   PetscFunctionReturn(0);
12540c7d97c5SJed Brown }
12550c7d97c5SJed Brown 
12560c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
12570c7d97c5SJed Brown /*
125850efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
12590c7d97c5SJed Brown 
12600c7d97c5SJed Brown    Input Parameters:
12610c7d97c5SJed Brown .  pc - the preconditioner context
12620c7d97c5SJed Brown .  r - input vector (global)
12630c7d97c5SJed Brown 
12640c7d97c5SJed Brown    Output Parameter:
12650c7d97c5SJed Brown .  z - output vector (global)
12660c7d97c5SJed Brown 
12670c7d97c5SJed Brown    Application Interface Routine: PCApply()
12680c7d97c5SJed Brown  */
12690c7d97c5SJed Brown #undef __FUNCT__
12700c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
127153cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
12720c7d97c5SJed Brown {
12730c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
12740c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
12750c7d97c5SJed Brown   PetscErrorCode    ierr;
12763b03a366Sstefano_zampini   const PetscScalar one = 1.0;
12773b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
12782617d88aSStefano Zampini   const PetscScalar zero = 0.0;
12790c7d97c5SJed Brown 
12800c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
12810c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
12828eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
12830c7d97c5SJed Brown 
12840c7d97c5SJed Brown   PetscFunctionBegin;
128585c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
12860c7d97c5SJed Brown     /* First Dirichlet solve */
12870c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12880c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
128953cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
12900c7d97c5SJed Brown     /*
12910c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1292674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1293674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
12940c7d97c5SJed Brown     */
12950c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12960c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
12978eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
12980c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12990c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
13000c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13010c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1302674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1303b76ba322SStefano Zampini   } else {
13040bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1305b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1306674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1307b76ba322SStefano Zampini   }
1308b76ba322SStefano Zampini 
13092617d88aSStefano Zampini   /* Apply interface preconditioner
13102617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1311dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
13122617d88aSStefano Zampini 
1313674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1314674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
13150c7d97c5SJed Brown 
13163b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
13170c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13180c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
13190c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
13208eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1321df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1322df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1323df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1324df187020SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
13250c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13260c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13270c7d97c5SJed Brown   PetscFunctionReturn(0);
13280c7d97c5SJed Brown }
132950efa1b5SStefano Zampini 
133050efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
133150efa1b5SStefano Zampini /*
133250efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
133350efa1b5SStefano Zampini 
133450efa1b5SStefano Zampini    Input Parameters:
133550efa1b5SStefano Zampini .  pc - the preconditioner context
133650efa1b5SStefano Zampini .  r - input vector (global)
133750efa1b5SStefano Zampini 
133850efa1b5SStefano Zampini    Output Parameter:
133950efa1b5SStefano Zampini .  z - output vector (global)
134050efa1b5SStefano Zampini 
134150efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
134250efa1b5SStefano Zampini  */
134350efa1b5SStefano Zampini #undef __FUNCT__
134450efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
134550efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
134650efa1b5SStefano Zampini {
134750efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
134850efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
134950efa1b5SStefano Zampini   PetscErrorCode    ierr;
135050efa1b5SStefano Zampini   const PetscScalar one = 1.0;
135150efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
135250efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
135350efa1b5SStefano Zampini 
135450efa1b5SStefano Zampini   PetscFunctionBegin;
135550efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
135650efa1b5SStefano Zampini     /* First Dirichlet solve */
135750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
135850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
135950efa1b5SStefano Zampini     ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
136050efa1b5SStefano Zampini     /*
136150efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
136250efa1b5SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
136350efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
136450efa1b5SStefano Zampini     */
136550efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
136650efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
136750efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
136850efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
136950efa1b5SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
137050efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
137150efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
137250efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
137350efa1b5SStefano Zampini   } else {
137450efa1b5SStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
137550efa1b5SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
137650efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
137750efa1b5SStefano Zampini   }
137850efa1b5SStefano Zampini 
137950efa1b5SStefano Zampini   /* Apply interface preconditioner
138050efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1381dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
138250efa1b5SStefano Zampini 
138350efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
138450efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
138550efa1b5SStefano Zampini 
138650efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
138750efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
138850efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
138950efa1b5SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
139050efa1b5SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1391b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1392b0147a47SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1393b0147a47SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1394b0147a47SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
139550efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
139650efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
139750efa1b5SStefano Zampini   PetscFunctionReturn(0);
139850efa1b5SStefano Zampini }
1399da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1400674ae819SStefano Zampini 
1401da1bb401SStefano Zampini #undef __FUNCT__
1402da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1403da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1404da1bb401SStefano Zampini {
1405da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1406da1bb401SStefano Zampini   PetscErrorCode ierr;
1407da1bb401SStefano Zampini 
1408da1bb401SStefano Zampini   PetscFunctionBegin;
1409da1bb401SStefano Zampini   /* free data created by PCIS */
1410da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1411674ae819SStefano Zampini   /* free BDDC custom data  */
1412674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1413674ae819SStefano Zampini   /* destroy objects related to topography */
1414674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1415674ae819SStefano Zampini   /* free allocated graph structure */
1416da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
141734a97f8cSStefano Zampini   /* destroy objects for scaling operator */
1418674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
141934a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
1420674ae819SStefano Zampini   /* free solvers stuff */
1421674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
142262a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
142362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
142462a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1425906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1426906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1427906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1428906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1429906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1430906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1431906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1432906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1433906d46d4SStefano Zampini   }
1434906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
1435906d46d4SStefano Zampini   /* remove map from local boundary to local numbering */
143639e2fb2aSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
14373425bc38SStefano Zampini   /* remove functions */
1438906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1439674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1440bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
14412b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1442b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
14432b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1444bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1445bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
144682ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1447bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
144882ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1449bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
145082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1451bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1452785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1453bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
145463602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1455bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1456bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1457bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1458bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1459674ae819SStefano Zampini   /* Free the private data structure */
1460674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1461da1bb401SStefano Zampini   PetscFunctionReturn(0);
1462da1bb401SStefano Zampini }
14633425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
14641e6b0712SBarry Smith 
14653425bc38SStefano Zampini #undef __FUNCT__
14663425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
14673425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
14683425bc38SStefano Zampini {
1469674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
14703425bc38SStefano Zampini   PC_IS*         pcis;
14713425bc38SStefano Zampini   PC_BDDC*       pcbddc;
14723425bc38SStefano Zampini   PetscErrorCode ierr;
14730c7d97c5SJed Brown 
14743425bc38SStefano Zampini   PetscFunctionBegin;
14753425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
14763425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
14773425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
14783425bc38SStefano Zampini 
14793425bc38SStefano Zampini   /* change of basis for physical rhs if needed
14803425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
14813308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
14823425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
14833425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14843425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1485fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1486fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
14873425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14883425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1489674ae819SStefano Zampini   /* Apply partition of unity */
14903425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1491674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
14928eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
14933425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
14943425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14953425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
14963425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
14973425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
14983425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14993425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1500674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
15013425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
15023425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
15033425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
15043425bc38SStefano Zampini   }
15053425bc38SStefano Zampini   /* BDDC rhs */
15063425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
15078eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
15083425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15093425bc38SStefano Zampini   }
15103425bc38SStefano Zampini   /* apply BDDC */
1511dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
15123425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
15133425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
15143425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
15153425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
15163425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
15173425bc38SStefano Zampini   /* restore original rhs */
15183425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
15193425bc38SStefano Zampini   PetscFunctionReturn(0);
15203425bc38SStefano Zampini }
15211e6b0712SBarry Smith 
15223425bc38SStefano Zampini #undef __FUNCT__
15233425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
15243425bc38SStefano Zampini /*@
152528509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
15263425bc38SStefano Zampini 
15273425bc38SStefano Zampini    Collective
15283425bc38SStefano Zampini 
15293425bc38SStefano Zampini    Input Parameters:
153028509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
153128509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
15323425bc38SStefano Zampini 
15333425bc38SStefano Zampini    Output Parameters:
153428509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
15353425bc38SStefano Zampini 
15363425bc38SStefano Zampini    Level: developer
15373425bc38SStefano Zampini 
15383425bc38SStefano Zampini    Notes:
15393425bc38SStefano Zampini 
154028509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
15413425bc38SStefano Zampini @*/
15423425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
15433425bc38SStefano Zampini {
1544674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15453425bc38SStefano Zampini   PetscErrorCode ierr;
15463425bc38SStefano Zampini 
15473425bc38SStefano Zampini   PetscFunctionBegin;
15483425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15493425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
15503425bc38SStefano Zampini   PetscFunctionReturn(0);
15513425bc38SStefano Zampini }
15523425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
15531e6b0712SBarry Smith 
15543425bc38SStefano Zampini #undef __FUNCT__
15553425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
15563425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
15573425bc38SStefano Zampini {
1558674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15593425bc38SStefano Zampini   PC_IS*         pcis;
15603425bc38SStefano Zampini   PC_BDDC*       pcbddc;
15613425bc38SStefano Zampini   PetscErrorCode ierr;
15623425bc38SStefano Zampini 
15633425bc38SStefano Zampini   PetscFunctionBegin;
15643425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15653425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
15663425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
15673425bc38SStefano Zampini 
15683425bc38SStefano Zampini   /* apply B_delta^T */
15693425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15703425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15713425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
15723425bc38SStefano Zampini   /* compute rhs for BDDC application */
15733425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
15748eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
15753425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15763425bc38SStefano Zampini   }
15773425bc38SStefano Zampini   /* apply BDDC */
1578dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
15793425bc38SStefano Zampini   /* put values into standard global vector */
15803425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15813425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15828eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
15833425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
15843425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
15853425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
15863425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15873425bc38SStefano Zampini   }
15883425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15893425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15903425bc38SStefano Zampini   /* final change of basis if needed
15913425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
15923308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
15933425bc38SStefano Zampini   PetscFunctionReturn(0);
15943425bc38SStefano Zampini }
15951e6b0712SBarry Smith 
15963425bc38SStefano Zampini #undef __FUNCT__
15973425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
15983425bc38SStefano Zampini /*@
159928509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
16003425bc38SStefano Zampini 
16013425bc38SStefano Zampini    Collective
16023425bc38SStefano Zampini 
16033425bc38SStefano Zampini    Input Parameters:
160428509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
160528509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
16063425bc38SStefano Zampini 
16073425bc38SStefano Zampini    Output Parameters:
160828509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
16093425bc38SStefano Zampini 
16103425bc38SStefano Zampini    Level: developer
16113425bc38SStefano Zampini 
16123425bc38SStefano Zampini    Notes:
16133425bc38SStefano Zampini 
161428509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
16153425bc38SStefano Zampini @*/
16163425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
16173425bc38SStefano Zampini {
1618674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
16193425bc38SStefano Zampini   PetscErrorCode ierr;
16203425bc38SStefano Zampini 
16213425bc38SStefano Zampini   PetscFunctionBegin;
16223425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16233425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
16243425bc38SStefano Zampini   PetscFunctionReturn(0);
16253425bc38SStefano Zampini }
16263425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
16271e6b0712SBarry Smith 
1628f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1629edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1630f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1631f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1632edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1633f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1634674ae819SStefano Zampini 
16353425bc38SStefano Zampini #undef __FUNCT__
16363425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
16373425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
16383425bc38SStefano Zampini {
1639674ae819SStefano Zampini 
1640674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
16413425bc38SStefano Zampini   Mat            newmat;
1642674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
16433425bc38SStefano Zampini   PC             newpc;
1644ce94432eSBarry Smith   MPI_Comm       comm;
16453425bc38SStefano Zampini   PetscErrorCode ierr;
16463425bc38SStefano Zampini 
16473425bc38SStefano Zampini   PetscFunctionBegin;
1648ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
16493425bc38SStefano Zampini   /* FETIDP linear matrix */
16503425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
16513425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
16523425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
16533425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1654edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
16553425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
16563425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
16573425bc38SStefano Zampini   /* FETIDP preconditioner */
16583425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
16593425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
16603425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
16613425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
16623425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
16633425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1664edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
16653425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
166623ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
16673425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
16683425bc38SStefano Zampini   /* return pointers for objects created */
16693425bc38SStefano Zampini   *fetidp_mat=newmat;
16703425bc38SStefano Zampini   *fetidp_pc=newpc;
16713425bc38SStefano Zampini   PetscFunctionReturn(0);
16723425bc38SStefano Zampini }
16731e6b0712SBarry Smith 
16743425bc38SStefano Zampini #undef __FUNCT__
16753425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
16763425bc38SStefano Zampini /*@
167728509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
16783425bc38SStefano Zampini 
16793425bc38SStefano Zampini    Collective
16803425bc38SStefano Zampini 
16813425bc38SStefano Zampini    Input Parameters:
168228509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
168328509bceSStefano Zampini 
168428509bceSStefano Zampini    Output Parameters:
168528509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
168628509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
168728509bceSStefano Zampini 
168828509bceSStefano Zampini    Options Database Keys:
168928509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
16903425bc38SStefano Zampini 
16913425bc38SStefano Zampini    Level: developer
16923425bc38SStefano Zampini 
16933425bc38SStefano Zampini    Notes:
169428509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
16953425bc38SStefano Zampini 
16963425bc38SStefano Zampini .seealso: PCBDDC
16973425bc38SStefano Zampini @*/
16983425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
16993425bc38SStefano Zampini {
17003425bc38SStefano Zampini   PetscErrorCode ierr;
17013425bc38SStefano Zampini 
17023425bc38SStefano Zampini   PetscFunctionBegin;
17033425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17043425bc38SStefano Zampini   if (pc->setupcalled) {
1705516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1706f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
17073425bc38SStefano Zampini   PetscFunctionReturn(0);
17083425bc38SStefano Zampini }
17090c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1710da1bb401SStefano Zampini /*MC
1711da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
17120c7d97c5SJed Brown 
171328509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
171428509bceSStefano Zampini 
171528509bceSStefano Zampini .vb
171628509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
171728509bceSStefano 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
171828509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
171928509bceSStefano Zampini .ve
172028509bceSStefano Zampini 
172128509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
172228509bceSStefano Zampini 
1723b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
172428509bceSStefano Zampini 
172528509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
172628509bceSStefano Zampini 
1727b6fdb6dfSStefano 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.
1728b6fdb6dfSStefano Zampini 
172928509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
173028509bceSStefano Zampini 
17316a818285SBarry 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()
173228509bceSStefano Zampini 
17336a818285SBarry Smith    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().
173428509bceSStefano Zampini 
1735b6fdb6dfSStefano 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.
173628509bceSStefano Zampini 
173728509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
173828509bceSStefano Zampini 
1739da1bb401SStefano Zampini    Options Database Keys:
174028509bceSStefano Zampini 
174128509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
174228509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1743b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
174428509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
174528509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
174628509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
174728509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
174828509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
174928509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
175028509bceSStefano Zampini 
175128509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
175228509bceSStefano Zampini .vb
175328509bceSStefano Zampini       -pc_bddc_dirichlet_
175428509bceSStefano Zampini       -pc_bddc_neumann_
175528509bceSStefano Zampini       -pc_bddc_coarse_
175628509bceSStefano Zampini .ve
175728509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
175828509bceSStefano Zampini 
175928509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
176028509bceSStefano Zampini .vb
1761312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1762312be037SStefano Zampini       -pc_bddc_neumann_lN_
1763312be037SStefano Zampini       -pc_bddc_coarse_lN_
176428509bceSStefano Zampini .ve
1765312be037SStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N.
1766da1bb401SStefano Zampini 
1767da1bb401SStefano Zampini    Level: intermediate
1768da1bb401SStefano Zampini 
1769b6fdb6dfSStefano Zampini    Developer notes:
1770da1bb401SStefano Zampini 
177128509bceSStefano Zampini    New deluxe scaling operator will be available soon.
1772da1bb401SStefano Zampini 
1773da1bb401SStefano Zampini    Contributed by Stefano Zampini
1774da1bb401SStefano Zampini 
1775da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1776da1bb401SStefano Zampini M*/
1777b2573a8aSBarry Smith 
1778da1bb401SStefano Zampini #undef __FUNCT__
1779da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
17808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1781da1bb401SStefano Zampini {
1782da1bb401SStefano Zampini   PetscErrorCode      ierr;
1783da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1784da1bb401SStefano Zampini 
1785da1bb401SStefano Zampini   PetscFunctionBegin;
1786da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1787b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1788da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1789da1bb401SStefano Zampini 
1790da1bb401SStefano Zampini   /* create PCIS data structure */
1791da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1792da1bb401SStefano Zampini 
179347d04d0dSStefano Zampini   /* BDDC customization */
179408a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
179547d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
179647d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
179747d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
179847d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
179947d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
180047d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
180147d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
180247d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1803b9d89cd5SStefano Zampini   /* private */
1804b9d89cd5SStefano Zampini   pcbddc->issym                      = PETSC_FALSE;
180539e2fb2aSStefano Zampini   pcbddc->BtoNmap                    = 0;
1806727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1807e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1808e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1809e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1810727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1811fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
181268457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
1813f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1814727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1815f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1816f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1817f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1818674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
18190bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
18203972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1821534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1822534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1823534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1824b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
1825906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
1826da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1827da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1828da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1829da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1830da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
183115aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
183215aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1833da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1834da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1835da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1836da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1837da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1838da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1839da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1840da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1841da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1842da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1843785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1844785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1845785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
184660d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
184760d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
184863602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
1849da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
185063602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
1851da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
185285c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
185347d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
185447d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
18554fad6a16SStefano Zampini   pcbddc->current_level              = 0;
18562b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1857323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
1858f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
1859323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
1860da1bb401SStefano Zampini 
1861674ae819SStefano Zampini   /* create local graph structure */
1862674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1863674ae819SStefano Zampini 
1864674ae819SStefano Zampini   /* scaling */
1865674ae819SStefano Zampini   pcbddc->work_scaling          = 0;
186634a97f8cSStefano Zampini   pcbddc->use_deluxe_scaling    = PETSC_FALSE;
186734a97f8cSStefano Zampini   pcbddc->deluxe_threshold      = 1;
186834a97f8cSStefano Zampini   pcbddc->deluxe_rebuild        = PETSC_FALSE;
186934a97f8cSStefano Zampini   pcbddc->deluxe_layers         = -1;
187034a97f8cSStefano Zampini   pcbddc->deluxe_use_useradj    = PETSC_FALSE;
187134a97f8cSStefano Zampini   pcbddc->deluxe_compute_rowadj = PETSC_TRUE;
1872da1bb401SStefano Zampini 
1873da1bb401SStefano Zampini   /* function pointers */
1874da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
187593bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1876da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1877da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1878da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1879da1bb401SStefano Zampini   pc->ops->view                = 0;
1880da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1881da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1882da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1883534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1884534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1885da1bb401SStefano Zampini 
1886da1bb401SStefano Zampini   /* composing function */
1887906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
1888674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1889bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
18902b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1891b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
18922b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1893bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1894bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
189582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1896bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
189782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1898bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
189982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1900bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
190182ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1902bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
190363602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1904bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1905bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1906bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1907bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1908da1bb401SStefano Zampini   PetscFunctionReturn(0);
1909da1bb401SStefano Zampini }
19103425bc38SStefano Zampini 
1911