xref: /petsc/src/ksp/pc/impls/bddc/bddc.c (revision 906d46d482912ff65fe95caf6add84f15a3f7e20)
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);
750c7d97c5SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
760c7d97c5SJed Brown   PetscFunctionReturn(0);
770c7d97c5SJed Brown }
780c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
79674ae819SStefano Zampini #undef __FUNCT__
80*906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC"
81*906d46d4SStefano Zampini static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change)
82b9b85e73SStefano Zampini {
83b9b85e73SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
84b9b85e73SStefano Zampini   PetscErrorCode ierr;
85b9b85e73SStefano Zampini 
86b9b85e73SStefano Zampini   PetscFunctionBegin;
87b9b85e73SStefano Zampini   ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr);
88b9b85e73SStefano Zampini   ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr);
89b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix = change;
90b9b85e73SStefano Zampini   PetscFunctionReturn(0);
91b9b85e73SStefano Zampini }
92b9b85e73SStefano Zampini #undef __FUNCT__
93*906d46d4SStefano Zampini #define __FUNCT__ "PCBDDCSetChangeOfBasisMat"
94b9b85e73SStefano Zampini /*@
95*906d46d4SStefano Zampini  PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
96b9b85e73SStefano Zampini 
97b9b85e73SStefano Zampini    Collective on PC
98b9b85e73SStefano Zampini 
99b9b85e73SStefano Zampini    Input Parameters:
100b9b85e73SStefano Zampini +  pc - the preconditioning context
101*906d46d4SStefano Zampini -  change - the change of basis matrix
102b9b85e73SStefano Zampini 
103b9b85e73SStefano Zampini    Level: intermediate
104b9b85e73SStefano Zampini 
105b9b85e73SStefano Zampini    Notes:
106b9b85e73SStefano Zampini 
107b9b85e73SStefano Zampini .seealso: PCBDDC
108b9b85e73SStefano Zampini @*/
109*906d46d4SStefano Zampini PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change)
110b9b85e73SStefano Zampini {
111b9b85e73SStefano Zampini   PetscErrorCode ierr;
112b9b85e73SStefano Zampini 
113b9b85e73SStefano Zampini   PetscFunctionBegin;
114b9b85e73SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115b9b85e73SStefano Zampini   PetscValidHeaderSpecific(change,MAT_CLASSID,2);
116*906d46d4SStefano Zampini   PetscCheckSameComm(pc,1,change,2);
117*906d46d4SStefano Zampini   if (pc->mat) {
118*906d46d4SStefano Zampini     PetscInt rows_c,cols_c,rows,cols;
119*906d46d4SStefano Zampini     ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
120*906d46d4SStefano Zampini     ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
121*906d46d4SStefano Zampini     if (rows_c != rows) {
122*906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows);
123*906d46d4SStefano Zampini     }
124*906d46d4SStefano Zampini     if (cols_c != cols) {
125*906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols);
126*906d46d4SStefano Zampini     }
127*906d46d4SStefano Zampini     ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr);
128*906d46d4SStefano Zampini     ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr);
129*906d46d4SStefano Zampini     if (rows_c != rows) {
130*906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows);
131*906d46d4SStefano Zampini     }
132*906d46d4SStefano Zampini     if (cols_c != cols) {
133*906d46d4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols);
134*906d46d4SStefano Zampini     }
135*906d46d4SStefano Zampini   }
136*906d46d4SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr);
137b9b85e73SStefano Zampini   PetscFunctionReturn(0);
138b9b85e73SStefano Zampini }
139b9b85e73SStefano Zampini /* -------------------------------------------------------------------------- */
140b9b85e73SStefano Zampini #undef __FUNCT__
141674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC"
142674ae819SStefano Zampini static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
143674ae819SStefano Zampini {
144674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
145674ae819SStefano Zampini   PetscErrorCode ierr;
1461e6b0712SBarry Smith 
147674ae819SStefano Zampini   PetscFunctionBegin;
148674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
149674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr);
150674ae819SStefano Zampini   pcbddc->user_primal_vertices = PrimalVertices;
151674ae819SStefano Zampini   PetscFunctionReturn(0);
152674ae819SStefano Zampini }
153674ae819SStefano Zampini #undef __FUNCT__
154674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS"
155674ae819SStefano Zampini /*@
15628509bceSStefano Zampini  PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
157674ae819SStefano Zampini 
158674ae819SStefano Zampini    Not collective
159674ae819SStefano Zampini 
160674ae819SStefano Zampini    Input Parameters:
161674ae819SStefano Zampini +  pc - the preconditioning context
16228509bceSStefano Zampini -  PrimalVertices - index set of primal vertices in local numbering
163674ae819SStefano Zampini 
164674ae819SStefano Zampini    Level: intermediate
165674ae819SStefano Zampini 
166674ae819SStefano Zampini    Notes:
167674ae819SStefano Zampini 
168674ae819SStefano Zampini .seealso: PCBDDC
169674ae819SStefano Zampini @*/
170674ae819SStefano Zampini PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
171674ae819SStefano Zampini {
172674ae819SStefano Zampini   PetscErrorCode ierr;
173674ae819SStefano Zampini 
174674ae819SStefano Zampini   PetscFunctionBegin;
175674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176674ae819SStefano Zampini   PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2);
177674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr);
178674ae819SStefano Zampini   PetscFunctionReturn(0);
179674ae819SStefano Zampini }
180674ae819SStefano Zampini /* -------------------------------------------------------------------------- */
1810c7d97c5SJed Brown #undef __FUNCT__
1824fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC"
1834fad6a16SStefano Zampini static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
1844fad6a16SStefano Zampini {
1854fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
1864fad6a16SStefano Zampini 
1874fad6a16SStefano Zampini   PetscFunctionBegin;
1884fad6a16SStefano Zampini   pcbddc->coarsening_ratio = k;
1894fad6a16SStefano Zampini   PetscFunctionReturn(0);
1904fad6a16SStefano Zampini }
1911e6b0712SBarry Smith 
1924fad6a16SStefano Zampini #undef __FUNCT__
1934fad6a16SStefano Zampini #define __FUNCT__ "PCBDDCSetCoarseningRatio"
1944fad6a16SStefano Zampini /*@
19528509bceSStefano Zampini  PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
1964fad6a16SStefano Zampini 
1974fad6a16SStefano Zampini    Logically collective on PC
1984fad6a16SStefano Zampini 
1994fad6a16SStefano Zampini    Input Parameters:
2004fad6a16SStefano Zampini +  pc - the preconditioning context
20128509bceSStefano Zampini -  k - coarsening ratio (H/h at the coarser level)
2024fad6a16SStefano Zampini 
20328509bceSStefano Zampini    Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
2044fad6a16SStefano Zampini 
2054fad6a16SStefano Zampini    Level: intermediate
2064fad6a16SStefano Zampini 
2074fad6a16SStefano Zampini    Notes:
2084fad6a16SStefano Zampini 
2094fad6a16SStefano Zampini .seealso: PCBDDC
2104fad6a16SStefano Zampini @*/
2114fad6a16SStefano Zampini PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
2124fad6a16SStefano Zampini {
2134fad6a16SStefano Zampini   PetscErrorCode ierr;
2144fad6a16SStefano Zampini 
2154fad6a16SStefano Zampini   PetscFunctionBegin;
2164fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2172b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,k,2);
2184fad6a16SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr);
2194fad6a16SStefano Zampini   PetscFunctionReturn(0);
2204fad6a16SStefano Zampini }
2212b510759SStefano Zampini 
222b8ffe317SStefano Zampini /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
2232b510759SStefano Zampini #undef __FUNCT__
224b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC"
225b8ffe317SStefano Zampini static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
226b8ffe317SStefano Zampini {
227b8ffe317SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
228b8ffe317SStefano Zampini 
229b8ffe317SStefano Zampini   PetscFunctionBegin;
23085c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick = flg;
231b8ffe317SStefano Zampini   PetscFunctionReturn(0);
232b8ffe317SStefano Zampini }
233b8ffe317SStefano Zampini 
234b8ffe317SStefano Zampini #undef __FUNCT__
235b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
236b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
2372b510759SStefano Zampini {
2382b510759SStefano Zampini   PetscErrorCode ierr;
2392b510759SStefano Zampini 
2402b510759SStefano Zampini   PetscFunctionBegin;
2412b510759SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
242b8ffe317SStefano Zampini   PetscValidLogicalCollectiveBool(pc,flg,2);
243b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
2442b510759SStefano Zampini   PetscFunctionReturn(0);
2452b510759SStefano Zampini }
2461e6b0712SBarry Smith 
2474fad6a16SStefano Zampini #undef __FUNCT__
2482b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel_BDDC"
2492b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
2504fad6a16SStefano Zampini {
2514fad6a16SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2524fad6a16SStefano Zampini 
2534fad6a16SStefano Zampini   PetscFunctionBegin;
2542b510759SStefano Zampini   pcbddc->current_level = level;
2554fad6a16SStefano Zampini   PetscFunctionReturn(0);
2564fad6a16SStefano Zampini }
2571e6b0712SBarry Smith 
2584fad6a16SStefano Zampini #undef __FUNCT__
259b8ffe317SStefano Zampini #define __FUNCT__ "PCBDDCSetLevel"
260b8ffe317SStefano Zampini PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
261b8ffe317SStefano Zampini {
262b8ffe317SStefano Zampini   PetscErrorCode ierr;
263b8ffe317SStefano Zampini 
264b8ffe317SStefano Zampini   PetscFunctionBegin;
265b8ffe317SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
266b8ffe317SStefano Zampini   PetscValidLogicalCollectiveInt(pc,level,2);
267b8ffe317SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr);
268b8ffe317SStefano Zampini   PetscFunctionReturn(0);
269b8ffe317SStefano Zampini }
270b8ffe317SStefano Zampini 
271b8ffe317SStefano Zampini #undef __FUNCT__
2722b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels_BDDC"
2732b510759SStefano Zampini static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
2742b510759SStefano Zampini {
2752b510759SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
2762b510759SStefano Zampini 
2772b510759SStefano Zampini   PetscFunctionBegin;
2782b510759SStefano Zampini   pcbddc->max_levels = levels;
2792b510759SStefano Zampini   PetscFunctionReturn(0);
2802b510759SStefano Zampini }
2812b510759SStefano Zampini 
2822b510759SStefano Zampini #undef __FUNCT__
2832b510759SStefano Zampini #define __FUNCT__ "PCBDDCSetLevels"
2844fad6a16SStefano Zampini /*@
28528509bceSStefano Zampini  PCBDDCSetLevels - Sets the maximum number of levels for multilevel
2864fad6a16SStefano Zampini 
2874fad6a16SStefano Zampini    Logically collective on PC
2884fad6a16SStefano Zampini 
2894fad6a16SStefano Zampini    Input Parameters:
2904fad6a16SStefano Zampini +  pc - the preconditioning context
29128509bceSStefano Zampini -  levels - the maximum number of levels (max 9)
2924fad6a16SStefano Zampini 
29328509bceSStefano Zampini    Default value is 0, i.e. traditional one-level BDDC
2944fad6a16SStefano Zampini 
2954fad6a16SStefano Zampini    Level: intermediate
2964fad6a16SStefano Zampini 
2974fad6a16SStefano Zampini    Notes:
2984fad6a16SStefano Zampini 
2994fad6a16SStefano Zampini .seealso: PCBDDC
3004fad6a16SStefano Zampini @*/
3012b510759SStefano Zampini PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
3024fad6a16SStefano Zampini {
3034fad6a16SStefano Zampini   PetscErrorCode ierr;
3044fad6a16SStefano Zampini 
3054fad6a16SStefano Zampini   PetscFunctionBegin;
3064fad6a16SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3072b510759SStefano Zampini   PetscValidLogicalCollectiveInt(pc,levels,2);
308312be037SStefano Zampini   if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n");
3092b510759SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
3104fad6a16SStefano Zampini   PetscFunctionReturn(0);
3114fad6a16SStefano Zampini }
3124fad6a16SStefano Zampini /* -------------------------------------------------------------------------- */
3131e6b0712SBarry Smith 
3144fad6a16SStefano Zampini #undef __FUNCT__
3150bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace_BDDC"
3160bdf917eSStefano Zampini static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace)
3170bdf917eSStefano Zampini {
3180bdf917eSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3190bdf917eSStefano Zampini   PetscErrorCode ierr;
3200bdf917eSStefano Zampini 
3210bdf917eSStefano Zampini   PetscFunctionBegin;
3220bdf917eSStefano Zampini   ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr);
3230bdf917eSStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
3240bdf917eSStefano Zampini   pcbddc->NullSpace = NullSpace;
3250bdf917eSStefano Zampini   PetscFunctionReturn(0);
3260bdf917eSStefano Zampini }
3271e6b0712SBarry Smith 
3280bdf917eSStefano Zampini #undef __FUNCT__
3290bdf917eSStefano Zampini #define __FUNCT__ "PCBDDCSetNullSpace"
3300bdf917eSStefano Zampini /*@
33128509bceSStefano Zampini  PCBDDCSetNullSpace - Set nullspace for BDDC operator
3320bdf917eSStefano Zampini 
3330bdf917eSStefano Zampini    Logically collective on PC and MatNullSpace
3340bdf917eSStefano Zampini 
3350bdf917eSStefano Zampini    Input Parameters:
3360bdf917eSStefano Zampini +  pc - the preconditioning context
33728509bceSStefano Zampini -  NullSpace - Null space of the linear operator to be preconditioned (Pmat)
3380bdf917eSStefano Zampini 
3390bdf917eSStefano Zampini    Level: intermediate
3400bdf917eSStefano Zampini 
3410bdf917eSStefano Zampini    Notes:
3420bdf917eSStefano Zampini 
3430bdf917eSStefano Zampini .seealso: PCBDDC
3440bdf917eSStefano Zampini @*/
3450bdf917eSStefano Zampini PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace)
3460bdf917eSStefano Zampini {
3470bdf917eSStefano Zampini   PetscErrorCode ierr;
3480bdf917eSStefano Zampini 
3490bdf917eSStefano Zampini   PetscFunctionBegin;
3500bdf917eSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
351674ae819SStefano Zampini   PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2);
3522b510759SStefano Zampini   PetscCheckSameComm(pc,1,NullSpace,2);
3530bdf917eSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr);
3540bdf917eSStefano Zampini   PetscFunctionReturn(0);
3550bdf917eSStefano Zampini }
3560bdf917eSStefano Zampini /* -------------------------------------------------------------------------- */
3571e6b0712SBarry Smith 
3580bdf917eSStefano Zampini #undef __FUNCT__
3593b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC"
3603b03a366Sstefano_zampini static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
3613b03a366Sstefano_zampini {
3623b03a366Sstefano_zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
3633b03a366Sstefano_zampini   PetscErrorCode ierr;
3643b03a366Sstefano_zampini 
3653b03a366Sstefano_zampini   PetscFunctionBegin;
366785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
367785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
3683b03a366Sstefano_zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
36936e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
37036e030ebSStefano Zampini   pcbddc->DirichletBoundaries = DirichletBoundaries;
371fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
3723b03a366Sstefano_zampini   PetscFunctionReturn(0);
3733b03a366Sstefano_zampini }
3741e6b0712SBarry Smith 
3753b03a366Sstefano_zampini #undef __FUNCT__
3763b03a366Sstefano_zampini #define __FUNCT__ "PCBDDCSetDirichletBoundaries"
3773b03a366Sstefano_zampini /*@
37828509bceSStefano Zampini  PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
3793b03a366Sstefano_zampini 
380785d1243SStefano Zampini    Collective
3813b03a366Sstefano_zampini 
3823b03a366Sstefano_zampini    Input Parameters:
3833b03a366Sstefano_zampini +  pc - the preconditioning context
384785d1243SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries
3853b03a366Sstefano_zampini 
3863b03a366Sstefano_zampini    Level: intermediate
3873b03a366Sstefano_zampini 
388785d1243SStefano Zampini    Notes: Any process can list any global node
3893b03a366Sstefano_zampini 
3903b03a366Sstefano_zampini .seealso: PCBDDC
3913b03a366Sstefano_zampini @*/
3923b03a366Sstefano_zampini PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
3933b03a366Sstefano_zampini {
3943b03a366Sstefano_zampini   PetscErrorCode ierr;
3953b03a366Sstefano_zampini 
3963b03a366Sstefano_zampini   PetscFunctionBegin;
3973b03a366Sstefano_zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
398674ae819SStefano Zampini   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
399785d1243SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
4003b03a366Sstefano_zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4013b03a366Sstefano_zampini   PetscFunctionReturn(0);
4023b03a366Sstefano_zampini }
4033b03a366Sstefano_zampini /* -------------------------------------------------------------------------- */
4041e6b0712SBarry Smith 
4053b03a366Sstefano_zampini #undef __FUNCT__
40682ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC"
40782ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
4080c7d97c5SJed Brown {
4090c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
4100c7d97c5SJed Brown   PetscErrorCode ierr;
4110c7d97c5SJed Brown 
4120c7d97c5SJed Brown   PetscFunctionBegin;
413785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
414785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr);
4150c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
4160c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr);
417785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
4187d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4190c7d97c5SJed Brown   PetscFunctionReturn(0);
4200c7d97c5SJed Brown }
4210c7d97c5SJed Brown 
4220c7d97c5SJed Brown #undef __FUNCT__
42382ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal"
4249c0446d6SStefano Zampini /*@
42582ba6b80SStefano Zampini  PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
4269c0446d6SStefano Zampini 
427785d1243SStefano Zampini    Collective
42853cdbc3dSStefano Zampini 
42953cdbc3dSStefano Zampini    Input Parameters:
43053cdbc3dSStefano Zampini +  pc - the preconditioning context
43182ba6b80SStefano Zampini -  DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
43253cdbc3dSStefano Zampini 
43353cdbc3dSStefano Zampini    Level: intermediate
43453cdbc3dSStefano Zampini 
4359c0446d6SStefano Zampini    Notes:
43653cdbc3dSStefano Zampini 
43753cdbc3dSStefano Zampini .seealso: PCBDDC
43853cdbc3dSStefano Zampini @*/
43982ba6b80SStefano Zampini PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
4400c7d97c5SJed Brown {
4410c7d97c5SJed Brown   PetscErrorCode ierr;
4420c7d97c5SJed Brown 
4430c7d97c5SJed Brown   PetscFunctionBegin;
4440c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4450c7d97c5SJed Brown   PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2);
44682ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,DirichletBoundaries,2);
44782ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr);
4480c7d97c5SJed Brown   PetscFunctionReturn(0);
4490c7d97c5SJed Brown }
4500c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
4510c7d97c5SJed Brown 
4520c7d97c5SJed Brown #undef __FUNCT__
4530c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC"
45453cdbc3dSStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
4550c7d97c5SJed Brown {
4560c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
45753cdbc3dSStefano Zampini   PetscErrorCode ierr;
4580c7d97c5SJed Brown 
4590c7d97c5SJed Brown   PetscFunctionBegin;
460785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
461785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
46253cdbc3dSStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
46336e030ebSStefano Zampini   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
46436e030ebSStefano Zampini   pcbddc->NeumannBoundaries = NeumannBoundaries;
465fb180af4SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
4660c7d97c5SJed Brown   PetscFunctionReturn(0);
4670c7d97c5SJed Brown }
4681e6b0712SBarry Smith 
4690c7d97c5SJed Brown #undef __FUNCT__
4700c7d97c5SJed Brown #define __FUNCT__ "PCBDDCSetNeumannBoundaries"
47157527edcSJed Brown /*@
47228509bceSStefano Zampini  PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
47357527edcSJed Brown 
474785d1243SStefano Zampini    Collective
47557527edcSJed Brown 
47657527edcSJed Brown    Input Parameters:
47757527edcSJed Brown +  pc - the preconditioning context
478785d1243SStefano Zampini -  NeumannBoundaries - parallel IS defining the Neumann boundaries
47957527edcSJed Brown 
48057527edcSJed Brown    Level: intermediate
48157527edcSJed Brown 
482785d1243SStefano Zampini    Notes: Any process can list any global node
48357527edcSJed Brown 
48457527edcSJed Brown .seealso: PCBDDC
48557527edcSJed Brown @*/
48653cdbc3dSStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
4870c7d97c5SJed Brown {
4880c7d97c5SJed Brown   PetscErrorCode ierr;
4890c7d97c5SJed Brown 
4900c7d97c5SJed Brown   PetscFunctionBegin;
4910c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
492674ae819SStefano Zampini   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
493785d1243SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
49453cdbc3dSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
49553cdbc3dSStefano Zampini   PetscFunctionReturn(0);
49653cdbc3dSStefano Zampini }
49753cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
4981e6b0712SBarry Smith 
49953cdbc3dSStefano Zampini #undef __FUNCT__
50082ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC"
50182ba6b80SStefano Zampini static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
5020c7d97c5SJed Brown {
5030c7d97c5SJed Brown   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
5040c7d97c5SJed Brown   PetscErrorCode ierr;
5050c7d97c5SJed Brown 
5060c7d97c5SJed Brown   PetscFunctionBegin;
507785d1243SStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
508785d1243SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr);
5090c7d97c5SJed Brown   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
5100c7d97c5SJed Brown   ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr);
511785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
5127d24d155SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
5130c7d97c5SJed Brown   PetscFunctionReturn(0);
5140c7d97c5SJed Brown }
5150c7d97c5SJed Brown 
5160c7d97c5SJed Brown #undef __FUNCT__
51782ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal"
5180c7d97c5SJed Brown /*@
51982ba6b80SStefano Zampini  PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
5200c7d97c5SJed Brown 
521785d1243SStefano Zampini    Collective
5220c7d97c5SJed Brown 
5230c7d97c5SJed Brown    Input Parameters:
5240c7d97c5SJed Brown +  pc - the preconditioning context
52582ba6b80SStefano Zampini -  NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
5260c7d97c5SJed Brown 
5270c7d97c5SJed Brown    Level: intermediate
5280c7d97c5SJed Brown 
5290c7d97c5SJed Brown    Notes:
5300c7d97c5SJed Brown 
5310c7d97c5SJed Brown .seealso: PCBDDC
5320c7d97c5SJed Brown @*/
53382ba6b80SStefano Zampini PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
5340c7d97c5SJed Brown {
5350c7d97c5SJed Brown   PetscErrorCode ierr;
5360c7d97c5SJed Brown 
5370c7d97c5SJed Brown   PetscFunctionBegin;
5380c7d97c5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5390c7d97c5SJed Brown   PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2);
54082ba6b80SStefano Zampini   PetscCheckSameComm(pc,1,NeumannBoundaries,2);
54182ba6b80SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr);
54253cdbc3dSStefano Zampini   PetscFunctionReturn(0);
54353cdbc3dSStefano Zampini }
54453cdbc3dSStefano Zampini /* -------------------------------------------------------------------------- */
54553cdbc3dSStefano Zampini 
54653cdbc3dSStefano Zampini #undef __FUNCT__
547da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC"
548da1bb401SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
549da1bb401SStefano Zampini {
550da1bb401SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
551da1bb401SStefano Zampini 
552da1bb401SStefano Zampini   PetscFunctionBegin;
553da1bb401SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundaries;
554da1bb401SStefano Zampini   PetscFunctionReturn(0);
555da1bb401SStefano Zampini }
5561e6b0712SBarry Smith 
557da1bb401SStefano Zampini #undef __FUNCT__
558da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundaries"
559da1bb401SStefano Zampini /*@
560785d1243SStefano Zampini  PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
561da1bb401SStefano Zampini 
562785d1243SStefano Zampini    Collective
563785d1243SStefano Zampini 
564785d1243SStefano Zampini    Input Parameters:
565785d1243SStefano Zampini .  pc - the preconditioning context
566785d1243SStefano Zampini 
567785d1243SStefano Zampini    Output Parameters:
568785d1243SStefano Zampini .  DirichletBoundaries - index set defining the Dirichlet boundaries
569785d1243SStefano Zampini 
570785d1243SStefano Zampini    Level: intermediate
571785d1243SStefano Zampini 
572785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
573785d1243SStefano Zampini 
574785d1243SStefano Zampini .seealso: PCBDDC
575785d1243SStefano Zampini @*/
576785d1243SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
577785d1243SStefano Zampini {
578785d1243SStefano Zampini   PetscErrorCode ierr;
579785d1243SStefano Zampini 
580785d1243SStefano Zampini   PetscFunctionBegin;
581785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
582785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
583785d1243SStefano Zampini   PetscFunctionReturn(0);
584785d1243SStefano Zampini }
585785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
586785d1243SStefano Zampini 
587785d1243SStefano Zampini #undef __FUNCT__
588785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC"
589785d1243SStefano Zampini static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
590785d1243SStefano Zampini {
591785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
592785d1243SStefano Zampini 
593785d1243SStefano Zampini   PetscFunctionBegin;
594785d1243SStefano Zampini   *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
595785d1243SStefano Zampini   PetscFunctionReturn(0);
596785d1243SStefano Zampini }
597785d1243SStefano Zampini 
598785d1243SStefano Zampini #undef __FUNCT__
59982ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal"
600da1bb401SStefano Zampini /*@
60182ba6b80SStefano Zampini  PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
602da1bb401SStefano Zampini 
603785d1243SStefano Zampini    Collective
604da1bb401SStefano Zampini 
605da1bb401SStefano Zampini    Input Parameters:
60628509bceSStefano Zampini .  pc - the preconditioning context
607da1bb401SStefano Zampini 
608da1bb401SStefano Zampini    Output Parameters:
60928509bceSStefano Zampini .  DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
610da1bb401SStefano Zampini 
611da1bb401SStefano Zampini    Level: intermediate
612da1bb401SStefano Zampini 
613da1bb401SStefano Zampini    Notes:
614da1bb401SStefano Zampini 
615da1bb401SStefano Zampini .seealso: PCBDDC
616da1bb401SStefano Zampini @*/
61782ba6b80SStefano Zampini PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
618da1bb401SStefano Zampini {
619da1bb401SStefano Zampini   PetscErrorCode ierr;
620da1bb401SStefano Zampini 
621da1bb401SStefano Zampini   PetscFunctionBegin;
622da1bb401SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
62382ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr);
624da1bb401SStefano Zampini   PetscFunctionReturn(0);
625da1bb401SStefano Zampini }
626da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
6271e6b0712SBarry Smith 
628da1bb401SStefano Zampini #undef __FUNCT__
62953cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC"
63053cdbc3dSStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
63153cdbc3dSStefano Zampini {
63253cdbc3dSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
63353cdbc3dSStefano Zampini 
63453cdbc3dSStefano Zampini   PetscFunctionBegin;
63553cdbc3dSStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundaries;
63653cdbc3dSStefano Zampini   PetscFunctionReturn(0);
63753cdbc3dSStefano Zampini }
6381e6b0712SBarry Smith 
63953cdbc3dSStefano Zampini #undef __FUNCT__
64053cdbc3dSStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundaries"
64153cdbc3dSStefano Zampini /*@
642785d1243SStefano Zampini  PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
64353cdbc3dSStefano Zampini 
644785d1243SStefano Zampini    Collective
645785d1243SStefano Zampini 
646785d1243SStefano Zampini    Input Parameters:
647785d1243SStefano Zampini .  pc - the preconditioning context
648785d1243SStefano Zampini 
649785d1243SStefano Zampini    Output Parameters:
650785d1243SStefano Zampini .  NeumannBoundaries - index set defining the Neumann boundaries
651785d1243SStefano Zampini 
652785d1243SStefano Zampini    Level: intermediate
653785d1243SStefano Zampini 
654785d1243SStefano Zampini    Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
655785d1243SStefano Zampini 
656785d1243SStefano Zampini .seealso: PCBDDC
657785d1243SStefano Zampini @*/
658785d1243SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
659785d1243SStefano Zampini {
660785d1243SStefano Zampini   PetscErrorCode ierr;
661785d1243SStefano Zampini 
662785d1243SStefano Zampini   PetscFunctionBegin;
663785d1243SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
664785d1243SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
665785d1243SStefano Zampini   PetscFunctionReturn(0);
666785d1243SStefano Zampini }
667785d1243SStefano Zampini /* -------------------------------------------------------------------------- */
668785d1243SStefano Zampini 
669785d1243SStefano Zampini #undef __FUNCT__
670785d1243SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC"
671785d1243SStefano Zampini static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
672785d1243SStefano Zampini {
673785d1243SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
674785d1243SStefano Zampini 
675785d1243SStefano Zampini   PetscFunctionBegin;
676785d1243SStefano Zampini   *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
677785d1243SStefano Zampini   PetscFunctionReturn(0);
678785d1243SStefano Zampini }
679785d1243SStefano Zampini 
680785d1243SStefano Zampini #undef __FUNCT__
68182ba6b80SStefano Zampini #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal"
68253cdbc3dSStefano Zampini /*@
68382ba6b80SStefano Zampini  PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
68453cdbc3dSStefano Zampini 
685785d1243SStefano Zampini    Collective
68653cdbc3dSStefano Zampini 
68753cdbc3dSStefano Zampini    Input Parameters:
68828509bceSStefano Zampini .  pc - the preconditioning context
68953cdbc3dSStefano Zampini 
69053cdbc3dSStefano Zampini    Output Parameters:
69128509bceSStefano Zampini .  NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
69253cdbc3dSStefano Zampini 
69353cdbc3dSStefano Zampini    Level: intermediate
69453cdbc3dSStefano Zampini 
69553cdbc3dSStefano Zampini    Notes:
69653cdbc3dSStefano Zampini 
69753cdbc3dSStefano Zampini .seealso: PCBDDC
69853cdbc3dSStefano Zampini @*/
69982ba6b80SStefano Zampini PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
70053cdbc3dSStefano Zampini {
70153cdbc3dSStefano Zampini   PetscErrorCode ierr;
70253cdbc3dSStefano Zampini 
70353cdbc3dSStefano Zampini   PetscFunctionBegin;
70453cdbc3dSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
70582ba6b80SStefano Zampini   ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr);
7060c7d97c5SJed Brown   PetscFunctionReturn(0);
7070c7d97c5SJed Brown }
70836e030ebSStefano Zampini /* -------------------------------------------------------------------------- */
7091e6b0712SBarry Smith 
71036e030ebSStefano Zampini #undef __FUNCT__
711da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC"
7121a83f524SJed Brown static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
71336e030ebSStefano Zampini {
71436e030ebSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
715da1bb401SStefano Zampini   PCBDDCGraph    mat_graph = pcbddc->mat_graph;
716da1bb401SStefano Zampini   PetscErrorCode ierr;
71736e030ebSStefano Zampini 
71836e030ebSStefano Zampini   PetscFunctionBegin;
719674ae819SStefano Zampini   /* free old CSR */
720674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr);
721fb223d50SStefano Zampini   /* TODO: PCBDDCGraphSetAdjacency */
722674ae819SStefano Zampini   /* get CSR into graph structure */
723da1bb401SStefano Zampini   if (copymode == PETSC_COPY_VALUES) {
724785e854fSJed Brown     ierr = PetscMalloc1((nvtxs+1),&mat_graph->xadj);CHKERRQ(ierr);
725785e854fSJed Brown     ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr);
726674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
727674ae819SStefano Zampini     ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr);
728da1bb401SStefano Zampini   } else if (copymode == PETSC_OWN_POINTER) {
7291a83f524SJed Brown     mat_graph->xadj = (PetscInt*)xadj;
7301a83f524SJed Brown     mat_graph->adjncy = (PetscInt*)adjncy;
731674ae819SStefano Zampini   }
732575ad6abSStefano Zampini   mat_graph->nvtxs_csr = nvtxs;
73336e030ebSStefano Zampini   PetscFunctionReturn(0);
73436e030ebSStefano Zampini }
7351e6b0712SBarry Smith 
73636e030ebSStefano Zampini #undef __FUNCT__
737da1bb401SStefano Zampini #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph"
73836e030ebSStefano Zampini /*@
73928509bceSStefano Zampini  PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix
74036e030ebSStefano Zampini 
74136e030ebSStefano Zampini    Not collective
74236e030ebSStefano Zampini 
74336e030ebSStefano Zampini    Input Parameters:
74436e030ebSStefano Zampini +  pc - the preconditioning context
74528509bceSStefano Zampini .  nvtxs - number of local vertices of the graph (i.e., the local size of your problem)
74628509bceSStefano Zampini .  xadj, adjncy - the CSR graph
74728509bceSStefano Zampini -  copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER.
74836e030ebSStefano Zampini 
74936e030ebSStefano Zampini    Level: intermediate
75036e030ebSStefano Zampini 
75136e030ebSStefano Zampini    Notes:
75236e030ebSStefano Zampini 
75328509bceSStefano Zampini .seealso: PCBDDC,PetscCopyMode
75436e030ebSStefano Zampini @*/
7551a83f524SJed Brown PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
75636e030ebSStefano Zampini {
757575ad6abSStefano Zampini   void (*f)(void) = 0;
75836e030ebSStefano Zampini   PetscErrorCode ierr;
75936e030ebSStefano Zampini 
76036e030ebSStefano Zampini   PetscFunctionBegin;
76136e030ebSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
762674ae819SStefano Zampini   PetscValidIntPointer(xadj,3);
763674ae819SStefano Zampini   PetscValidIntPointer(xadj,4);
764674ae819SStefano Zampini   if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) {
765674ae819SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__);
766da1bb401SStefano Zampini   }
76736e030ebSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr);
768575ad6abSStefano Zampini   /* free arrays if PCBDDC is not the PC type */
769575ad6abSStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr);
770575ad6abSStefano Zampini   if (!f && copymode == PETSC_OWN_POINTER) {
771575ad6abSStefano Zampini     ierr = PetscFree(xadj);CHKERRQ(ierr);
772575ad6abSStefano Zampini     ierr = PetscFree(adjncy);CHKERRQ(ierr);
77336e030ebSStefano Zampini   }
77436e030ebSStefano Zampini   PetscFunctionReturn(0);
77536e030ebSStefano Zampini }
7769c0446d6SStefano Zampini /* -------------------------------------------------------------------------- */
7771e6b0712SBarry Smith 
7789c0446d6SStefano Zampini #undef __FUNCT__
77963602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC"
78063602bcaSStefano Zampini static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
78163602bcaSStefano Zampini {
78263602bcaSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
78363602bcaSStefano Zampini   PetscInt i;
78463602bcaSStefano Zampini   PetscErrorCode ierr;
78563602bcaSStefano Zampini 
78663602bcaSStefano Zampini   PetscFunctionBegin;
78763602bcaSStefano Zampini   /* Destroy ISes if they were already set */
78863602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
78963602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
79063602bcaSStefano Zampini   }
79163602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
79263602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
79363602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
79463602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
79563602bcaSStefano Zampini   }
79663602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
79763602bcaSStefano Zampini   pcbddc->n_ISForDofs = 0;
79863602bcaSStefano Zampini   /* allocate space then set */
79963602bcaSStefano Zampini   ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr);
80063602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
80163602bcaSStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
80263602bcaSStefano Zampini     pcbddc->ISForDofsLocal[i]=ISForDofs[i];
80363602bcaSStefano Zampini   }
80463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal=n_is;
80563602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8061a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
80763602bcaSStefano Zampini   PetscFunctionReturn(0);
80863602bcaSStefano Zampini }
80963602bcaSStefano Zampini 
81063602bcaSStefano Zampini #undef __FUNCT__
81163602bcaSStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplittingLocal"
81263602bcaSStefano Zampini /*@
81363602bcaSStefano Zampini  PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
81463602bcaSStefano Zampini 
81563602bcaSStefano Zampini    Collective
81663602bcaSStefano Zampini 
81763602bcaSStefano Zampini    Input Parameters:
81863602bcaSStefano Zampini +  pc - the preconditioning context
81963602bcaSStefano Zampini -  n_is - number of index sets defining the fields
82063602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in local ordering
82163602bcaSStefano Zampini 
82263602bcaSStefano Zampini    Level: intermediate
82363602bcaSStefano Zampini 
82463602bcaSStefano 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.
82563602bcaSStefano Zampini 
82663602bcaSStefano Zampini .seealso: PCBDDC
82763602bcaSStefano Zampini @*/
82863602bcaSStefano Zampini PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
82963602bcaSStefano Zampini {
83063602bcaSStefano Zampini   PetscInt       i;
83163602bcaSStefano Zampini   PetscErrorCode ierr;
83263602bcaSStefano Zampini 
83363602bcaSStefano Zampini   PetscFunctionBegin;
83463602bcaSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
83563602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
83663602bcaSStefano Zampini   for (i=0;i<n_is;i++) {
83763602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
83863602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
83963602bcaSStefano Zampini   }
84063602bcaSStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
84163602bcaSStefano Zampini   PetscFunctionReturn(0);
84263602bcaSStefano Zampini }
84363602bcaSStefano Zampini /* -------------------------------------------------------------------------- */
84463602bcaSStefano Zampini 
84563602bcaSStefano Zampini #undef __FUNCT__
8469c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC"
8479c0446d6SStefano Zampini static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
8489c0446d6SStefano Zampini {
8499c0446d6SStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
8509c0446d6SStefano Zampini   PetscInt i;
8519c0446d6SStefano Zampini   PetscErrorCode ierr;
8529c0446d6SStefano Zampini 
8539c0446d6SStefano Zampini   PetscFunctionBegin;
854da1bb401SStefano Zampini   /* Destroy ISes if they were already set */
8559c0446d6SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
8569c0446d6SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
8579c0446d6SStefano Zampini   }
858d11ae9bbSstefano_zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
85963602bcaSStefano Zampini   /* last user setting takes precendence -> destroy any other customization */
86063602bcaSStefano Zampini   for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
86163602bcaSStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr);
86263602bcaSStefano Zampini   }
86363602bcaSStefano Zampini   ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr);
86463602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal = 0;
865da1bb401SStefano Zampini   /* allocate space then set */
866785e854fSJed Brown   ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr);
8679c0446d6SStefano Zampini   for (i=0;i<n_is;i++) {
868da1bb401SStefano Zampini     ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr);
869da1bb401SStefano Zampini     pcbddc->ISForDofs[i]=ISForDofs[i];
8709c0446d6SStefano Zampini   }
8719c0446d6SStefano Zampini   pcbddc->n_ISForDofs=n_is;
87263602bcaSStefano Zampini   if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
8731a8a16d1SStefano Zampini   pcbddc->recompute_topography = PETSC_TRUE;
8749c0446d6SStefano Zampini   PetscFunctionReturn(0);
8759c0446d6SStefano Zampini }
8761e6b0712SBarry Smith 
8779c0446d6SStefano Zampini #undef __FUNCT__
8789c0446d6SStefano Zampini #define __FUNCT__ "PCBDDCSetDofsSplitting"
8799c0446d6SStefano Zampini /*@
88063602bcaSStefano Zampini  PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
8819c0446d6SStefano Zampini 
88263602bcaSStefano Zampini    Collective
8839c0446d6SStefano Zampini 
8849c0446d6SStefano Zampini    Input Parameters:
8859c0446d6SStefano Zampini +  pc - the preconditioning context
88628509bceSStefano Zampini -  n_is - number of index sets defining the fields
88763602bcaSStefano Zampini .  ISForDofs - array of IS describing the fields in global ordering
8889c0446d6SStefano Zampini 
8899c0446d6SStefano Zampini    Level: intermediate
8909c0446d6SStefano Zampini 
89163602bcaSStefano Zampini    Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field.
8929c0446d6SStefano Zampini 
8939c0446d6SStefano Zampini .seealso: PCBDDC
8949c0446d6SStefano Zampini @*/
8959c0446d6SStefano Zampini PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
8969c0446d6SStefano Zampini {
8972b510759SStefano Zampini   PetscInt       i;
8989c0446d6SStefano Zampini   PetscErrorCode ierr;
8999c0446d6SStefano Zampini 
9009c0446d6SStefano Zampini   PetscFunctionBegin;
9019c0446d6SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
90263602bcaSStefano Zampini   PetscValidLogicalCollectiveInt(pc,n_is,2);
9032b510759SStefano Zampini   for (i=0;i<n_is;i++) {
90463602bcaSStefano Zampini     PetscCheckSameComm(pc,1,ISForDofs[i],3);
90563602bcaSStefano Zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3);
9062b510759SStefano Zampini   }
9079c0446d6SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr);
9089c0446d6SStefano Zampini   PetscFunctionReturn(0);
9099c0446d6SStefano Zampini }
910*906d46d4SStefano Zampini 
911da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
912534831adSStefano Zampini #undef __FUNCT__
913534831adSStefano Zampini #define __FUNCT__ "PCPreSolve_BDDC"
914534831adSStefano Zampini /* -------------------------------------------------------------------------- */
915534831adSStefano Zampini /*
916534831adSStefano Zampini    PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
917534831adSStefano Zampini                      guess if a transformation of basis approach has been selected.
9189c0446d6SStefano Zampini 
919534831adSStefano Zampini    Input Parameter:
920534831adSStefano Zampini +  pc - the preconditioner contex
921534831adSStefano Zampini 
922534831adSStefano Zampini    Application Interface Routine: PCPreSolve()
923534831adSStefano Zampini 
924534831adSStefano Zampini    Notes:
925534831adSStefano Zampini    The interface routine PCPreSolve() is not usually called directly by
926534831adSStefano Zampini    the user, but instead is called by KSPSolve().
927534831adSStefano Zampini */
928534831adSStefano Zampini static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
929534831adSStefano Zampini {
930534831adSStefano Zampini   PetscErrorCode ierr;
931534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
932534831adSStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
9333972b0daSStefano Zampini   IS             dirIS;
9343972b0daSStefano Zampini   Vec            used_vec;
935534831adSStefano Zampini 
936534831adSStefano Zampini   PetscFunctionBegin;
93785c4d303SStefano Zampini   /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */
93885c4d303SStefano Zampini   if (ksp) {
93985c4d303SStefano Zampini     PetscBool iscg;
94085c4d303SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr);
94185c4d303SStefano Zampini     if (!iscg) {
94285c4d303SStefano Zampini       ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr);
94385c4d303SStefano Zampini     }
94485c4d303SStefano Zampini   }
94585c4d303SStefano Zampini   /* Creates parallel work vectors used in presolve */
94662a6ff1dSStefano Zampini   if (!pcbddc->original_rhs) {
94762a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
94862a6ff1dSStefano Zampini   }
94962a6ff1dSStefano Zampini   if (!pcbddc->temp_solution) {
95062a6ff1dSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
95162a6ff1dSStefano Zampini   }
9523972b0daSStefano Zampini   if (x) {
9533972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
9543972b0daSStefano Zampini     used_vec = x;
9553972b0daSStefano Zampini   } else {
9563972b0daSStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr);
9573972b0daSStefano Zampini     used_vec = pcbddc->temp_solution;
9583972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9593972b0daSStefano Zampini   }
9603972b0daSStefano Zampini   /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */
9613972b0daSStefano Zampini   if (ksp) {
962*906d46d4SStefano Zampini     PetscBool guess_nonzero;
9633972b0daSStefano Zampini     ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr);
9643972b0daSStefano Zampini     if (!guess_nonzero) {
9653972b0daSStefano Zampini       ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
9663972b0daSStefano Zampini     }
9673972b0daSStefano Zampini   }
9683308cffdSStefano Zampini 
9693972b0daSStefano Zampini   /* store the original rhs */
9703972b0daSStefano Zampini   ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr);
9713972b0daSStefano Zampini 
9723972b0daSStefano Zampini   /* Take into account zeroed rows -> change rhs and store solution removed */
973785d1243SStefano Zampini   /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */
97482ba6b80SStefano Zampini   ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr);
97582ba6b80SStefano Zampini   if (rhs && dirIS) {
976*906d46d4SStefano Zampini     Mat_IS      *matis = (Mat_IS*)pc->pmat->data;
977785d1243SStefano Zampini     PetscInt    dirsize,i,*is_indices;
978785d1243SStefano Zampini     PetscScalar *array_x,*array_diagonal;
979785d1243SStefano Zampini 
9803972b0daSStefano Zampini     ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr);
9813972b0daSStefano Zampini     ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr);
9823972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9833972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9843972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9853972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
98682ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr);
9873972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9883972b0daSStefano Zampini     ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9893972b0daSStefano Zampini     ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9902fa5cd67SKarl Rupp     for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
9913972b0daSStefano Zampini     ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9923972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr);
9933972b0daSStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr);
9943972b0daSStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9953972b0daSStefano Zampini     ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
996b76ba322SStefano Zampini 
9973972b0daSStefano Zampini     /* remove the computed solution from the rhs */
9983972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
9993972b0daSStefano Zampini     ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr);
10003972b0daSStefano Zampini     ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr);
10013308cffdSStefano Zampini   }
1002b76ba322SStefano Zampini 
1003b76ba322SStefano Zampini   /* store partially computed solution and set initial guess */
10043972b0daSStefano Zampini   if (x) {
10053972b0daSStefano Zampini     ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr);
10063972b0daSStefano Zampini     ierr = VecSet(used_vec,0.0);CHKERRQ(ierr);
100785c4d303SStefano Zampini     if (pcbddc->use_exact_dirichlet_trick) {
1008b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1009b76ba322SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1010b76ba322SStefano Zampini       ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
1011b76ba322SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1012b76ba322SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1013b76ba322SStefano Zampini       if (ksp) {
1014b76ba322SStefano Zampini         ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
1015b76ba322SStefano Zampini       }
1016b76ba322SStefano Zampini     }
10173972b0daSStefano Zampini   }
1018b76ba322SStefano Zampini 
1019b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1020*906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1021*906d46d4SStefano Zampini 
1022*906d46d4SStefano Zampini     /* get change ctx */
1023*906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1024*906d46d4SStefano Zampini 
1025*906d46d4SStefano Zampini     /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */
1026*906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1027*906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1028*906d46d4SStefano Zampini     change_ctx->original_mat = pc->mat;
1029*906d46d4SStefano Zampini 
1030*906d46d4SStefano Zampini     /* change iteration matrix */
1031*906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1032*906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr);
1033*906d46d4SStefano Zampini     pc->mat = pcbddc->new_global_mat;
1034*906d46d4SStefano Zampini 
1035*906d46d4SStefano Zampini     /* change rhs */
1036*906d46d4SStefano Zampini     ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr);
1037*906d46d4SStefano Zampini     ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr);
103892e3dcfbSStefano Zampini   }
103992e3dcfbSStefano Zampini 
104092e3dcfbSStefano Zampini   /* remove nullspace if present */
10410bdf917eSStefano Zampini   if (ksp && pcbddc->NullSpace) {
1042d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr);
1043d0195637SJed Brown     ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr);
1044b76ba322SStefano Zampini   }
10450bdf917eSStefano Zampini   ierr = VecDestroy(&used_vec);CHKERRQ(ierr);
1046534831adSStefano Zampini   PetscFunctionReturn(0);
1047534831adSStefano Zampini }
1048*906d46d4SStefano Zampini 
1049534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1050534831adSStefano Zampini #undef __FUNCT__
1051534831adSStefano Zampini #define __FUNCT__ "PCPostSolve_BDDC"
1052534831adSStefano Zampini /* -------------------------------------------------------------------------- */
1053534831adSStefano Zampini /*
1054534831adSStefano Zampini    PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1055534831adSStefano Zampini                      approach has been selected. Also, restores rhs to its original state.
1056534831adSStefano Zampini 
1057534831adSStefano Zampini    Input Parameter:
1058534831adSStefano Zampini +  pc - the preconditioner contex
1059534831adSStefano Zampini 
1060534831adSStefano Zampini    Application Interface Routine: PCPostSolve()
1061534831adSStefano Zampini 
1062534831adSStefano Zampini    Notes:
1063534831adSStefano Zampini    The interface routine PCPostSolve() is not usually called directly by
1064534831adSStefano Zampini    the user, but instead is called by KSPSolve().
1065534831adSStefano Zampini */
1066534831adSStefano Zampini static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1067534831adSStefano Zampini {
1068534831adSStefano Zampini   PetscErrorCode ierr;
1069534831adSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1070534831adSStefano Zampini 
1071534831adSStefano Zampini   PetscFunctionBegin;
1072b9b85e73SStefano Zampini   if (pcbddc->ChangeOfBasisMatrix) {
1073*906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1074*906d46d4SStefano Zampini 
1075*906d46d4SStefano Zampini     /* get change ctx */
1076*906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1077*906d46d4SStefano Zampini 
1078*906d46d4SStefano Zampini     /* restore iteration matrix */
1079*906d46d4SStefano Zampini     ierr = MatDestroy(&pc->mat);CHKERRQ(ierr);
1080*906d46d4SStefano Zampini     ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr);
1081*906d46d4SStefano Zampini     pc->mat = change_ctx->original_mat;
1082*906d46d4SStefano Zampini 
1083*906d46d4SStefano Zampini     /* get solution in original basis */
1084*906d46d4SStefano Zampini     if (x) {
1085*906d46d4SStefano Zampini       PC_IS *pcis = (PC_IS*)(pc->data);
1086*906d46d4SStefano Zampini       ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr);
1087*906d46d4SStefano Zampini       ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr);
10883425bc38SStefano Zampini     }
1089534831adSStefano Zampini   }
1090*906d46d4SStefano Zampini 
10913972b0daSStefano Zampini   /* add solution removed in presolve */
10923425bc38SStefano Zampini   if (x) {
10933425bc38SStefano Zampini     ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr);
10943425bc38SStefano Zampini   }
1095*906d46d4SStefano Zampini 
1096fb223d50SStefano Zampini   /* restore rhs to its original state */
1097fb223d50SStefano Zampini   if (rhs) {
1098fb223d50SStefano Zampini     ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr);
1099fb223d50SStefano Zampini   }
1100534831adSStefano Zampini   PetscFunctionReturn(0);
1101534831adSStefano Zampini }
1102534831adSStefano Zampini /* -------------------------------------------------------------------------- */
110353cdbc3dSStefano Zampini #undef __FUNCT__
110453cdbc3dSStefano Zampini #define __FUNCT__ "PCSetUp_BDDC"
11050c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
11060c7d97c5SJed Brown /*
11070c7d97c5SJed Brown    PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
11080c7d97c5SJed Brown                   by setting data structures and options.
11090c7d97c5SJed Brown 
11100c7d97c5SJed Brown    Input Parameter:
111153cdbc3dSStefano Zampini +  pc - the preconditioner context
11120c7d97c5SJed Brown 
11130c7d97c5SJed Brown    Application Interface Routine: PCSetUp()
11140c7d97c5SJed Brown 
11150c7d97c5SJed Brown    Notes:
11160c7d97c5SJed Brown    The interface routine PCSetUp() is not usually called directly by
11170c7d97c5SJed Brown    the user, but instead is called by PCApply() if necessary.
11180c7d97c5SJed Brown */
111953cdbc3dSStefano Zampini PetscErrorCode PCSetUp_BDDC(PC pc)
11200c7d97c5SJed Brown {
11210c7d97c5SJed Brown   PetscErrorCode   ierr;
11220c7d97c5SJed Brown   PC_BDDC*         pcbddc = (PC_BDDC*)pc->data;
1123f4ddd8eeSStefano Zampini   MatNullSpace     nearnullspace;
1124674ae819SStefano Zampini   PetscBool        computeis,computetopography,computesolvers;
1125165b64e2SStefano Zampini   PetscBool        new_nearnullspace_provided;
11260c7d97c5SJed Brown 
11270c7d97c5SJed Brown   PetscFunctionBegin;
1128f4ddd8eeSStefano Zampini   /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
11293b03a366Sstefano_zampini   /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1130fcd409f5SStefano Zampini      Also, BDDC directly build the Dirichlet problem */
1131f4ddd8eeSStefano Zampini   /* split work */
1132674ae819SStefano Zampini   if (pc->setupcalled) {
1133674ae819SStefano Zampini     computeis = PETSC_FALSE;
1134d1e9a80fSBarry Smith     if (pc->flag == SAME_NONZERO_PATTERN) {
1135674ae819SStefano Zampini       computetopography = PETSC_FALSE;
1136674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1137674ae819SStefano Zampini     } else { /* DIFFERENT_NONZERO_PATTERN */
1138674ae819SStefano Zampini       computetopography = PETSC_TRUE;
1139674ae819SStefano Zampini       computesolvers = PETSC_TRUE;
1140674ae819SStefano Zampini     }
1141674ae819SStefano Zampini   } else {
1142674ae819SStefano Zampini     computeis = PETSC_TRUE;
1143674ae819SStefano Zampini     computetopography = PETSC_TRUE;
1144674ae819SStefano Zampini     computesolvers = PETSC_TRUE;
1145674ae819SStefano Zampini   }
1146fb180af4SStefano Zampini   if (pcbddc->recompute_topography) {
1147fb180af4SStefano Zampini     computetopography = PETSC_TRUE;
1148fb180af4SStefano Zampini   }
1149f4ddd8eeSStefano Zampini 
1150f4ddd8eeSStefano Zampini   /* Get stdout for dbg */
115170cf5478SStefano Zampini   if (pcbddc->dbg_flag) {
115270cf5478SStefano Zampini     if (!pcbddc->dbg_viewer) {
115358a03d70SStefano Zampini       pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1154f4ddd8eeSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
1155f4ddd8eeSStefano Zampini     }
115658a03d70SStefano Zampini     ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
1157f4ddd8eeSStefano Zampini   }
1158f4ddd8eeSStefano Zampini 
1159fcd409f5SStefano Zampini   /* Set up all the "iterative substructuring" common block without computing solvers */
1160674ae819SStefano Zampini   if (computeis) {
1161fcd409f5SStefano Zampini     /* HACK INTO PCIS */
1162fcd409f5SStefano Zampini     PC_IS* pcis = (PC_IS*)pc->data;
1163fcd409f5SStefano Zampini     pcis->computesolvers = PETSC_FALSE;
1164674ae819SStefano Zampini     ierr = PCISSetUp(pc);CHKERRQ(ierr);
116539e2fb2aSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr);
1166674ae819SStefano Zampini   }
1167f4ddd8eeSStefano Zampini 
1168b9b85e73SStefano Zampini   /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1169*906d46d4SStefano Zampini   if (pcbddc->user_ChangeOfBasisMatrix) {
1170b9b85e73SStefano Zampini     pcbddc->use_change_of_basis = PETSC_FALSE;
1171b9b85e73SStefano Zampini   }
1172*906d46d4SStefano Zampini 
1173f4ddd8eeSStefano Zampini   /* Analyze interface */
1174674ae819SStefano Zampini   if (computetopography) {
1175674ae819SStefano Zampini     ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr);
1176674ae819SStefano Zampini   }
1177fb8d54d4SStefano Zampini 
1178f4ddd8eeSStefano Zampini   /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1179fb8d54d4SStefano Zampini   new_nearnullspace_provided = PETSC_FALSE;
1180f4ddd8eeSStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr);
1181f4ddd8eeSStefano Zampini   if (pcbddc->onearnullspace) { /* already used nearnullspace */
1182f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1183f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1184f4ddd8eeSStefano Zampini     } else {
1185f4ddd8eeSStefano Zampini       /* determine if the two nullspaces are different (should be lightweight) */
1186f4ddd8eeSStefano Zampini       if (nearnullspace != pcbddc->onearnullspace) {
1187f4ddd8eeSStefano Zampini         new_nearnullspace_provided = PETSC_TRUE;
1188165b64e2SStefano Zampini       } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1189f4ddd8eeSStefano Zampini         PetscInt         i;
1190165b64e2SStefano Zampini         const Vec        *nearnullvecs;
1191165b64e2SStefano Zampini         PetscObjectState state;
1192165b64e2SStefano Zampini         PetscInt         nnsp_size;
1193165b64e2SStefano Zampini         ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
1194f4ddd8eeSStefano Zampini         for (i=0;i<nnsp_size;i++) {
1195f4ddd8eeSStefano Zampini           ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr);
1196165b64e2SStefano Zampini           if (pcbddc->onearnullvecs_state[i] != state) {
1197f4ddd8eeSStefano Zampini             new_nearnullspace_provided = PETSC_TRUE;
1198f4ddd8eeSStefano Zampini             break;
1199f4ddd8eeSStefano Zampini           }
1200f4ddd8eeSStefano Zampini         }
1201f4ddd8eeSStefano Zampini       }
1202f4ddd8eeSStefano Zampini     }
1203f4ddd8eeSStefano Zampini   } else {
1204f4ddd8eeSStefano Zampini     if (!nearnullspace) { /* both nearnullspaces are null */
1205f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_FALSE;
1206f4ddd8eeSStefano Zampini     } else { /* nearnullspace attached later */
1207f4ddd8eeSStefano Zampini       new_nearnullspace_provided = PETSC_TRUE;
1208f4ddd8eeSStefano Zampini     }
1209f4ddd8eeSStefano Zampini   }
1210f4ddd8eeSStefano Zampini 
1211f4ddd8eeSStefano Zampini   /* Setup constraints and related work vectors */
1212727cdba6SStefano Zampini   /* reset primal space flags */
1213f4ddd8eeSStefano Zampini   pcbddc->new_primal_space = PETSC_FALSE;
1214727cdba6SStefano Zampini   pcbddc->new_primal_space_local = PETSC_FALSE;
1215fb8d54d4SStefano Zampini   if (computetopography || new_nearnullspace_provided) {
1216727cdba6SStefano Zampini     /* It also sets the primal space flags */
1217674ae819SStefano Zampini     ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr);
1218e7b262bdSStefano Zampini     /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1219f4ddd8eeSStefano Zampini     ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr);
1220674ae819SStefano Zampini   }
1221fb8d54d4SStefano Zampini 
1222f4ddd8eeSStefano Zampini   if (computesolvers || pcbddc->new_primal_space) {
1223674ae819SStefano Zampini     /* reset data */
1224674ae819SStefano Zampini     ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
122599cc7994SStefano Zampini     /* Create coarse and local stuffs */
122699cc7994SStefano Zampini     ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr);
1227674ae819SStefano Zampini     ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr);
12280c7d97c5SJed Brown   }
122970cf5478SStefano Zampini 
123058a03d70SStefano Zampini   if (pcbddc->dbg_flag) {
123158a03d70SStefano Zampini     ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr);
12322b510759SStefano Zampini   }
12330c7d97c5SJed Brown   PetscFunctionReturn(0);
12340c7d97c5SJed Brown }
12350c7d97c5SJed Brown 
12360c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
12370c7d97c5SJed Brown /*
123850efa1b5SStefano Zampini    PCApply_BDDC - Applies the BDDC operator to a vector.
12390c7d97c5SJed Brown 
12400c7d97c5SJed Brown    Input Parameters:
12410c7d97c5SJed Brown .  pc - the preconditioner context
12420c7d97c5SJed Brown .  r - input vector (global)
12430c7d97c5SJed Brown 
12440c7d97c5SJed Brown    Output Parameter:
12450c7d97c5SJed Brown .  z - output vector (global)
12460c7d97c5SJed Brown 
12470c7d97c5SJed Brown    Application Interface Routine: PCApply()
12480c7d97c5SJed Brown  */
12490c7d97c5SJed Brown #undef __FUNCT__
12500c7d97c5SJed Brown #define __FUNCT__ "PCApply_BDDC"
125153cdbc3dSStefano Zampini PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
12520c7d97c5SJed Brown {
12530c7d97c5SJed Brown   PC_IS             *pcis = (PC_IS*)(pc->data);
12540c7d97c5SJed Brown   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
12550c7d97c5SJed Brown   PetscErrorCode    ierr;
12563b03a366Sstefano_zampini   const PetscScalar one = 1.0;
12573b03a366Sstefano_zampini   const PetscScalar m_one = -1.0;
12582617d88aSStefano Zampini   const PetscScalar zero = 0.0;
12590c7d97c5SJed Brown 
12600c7d97c5SJed Brown /* This code is similar to that provided in nn.c for PCNN
12610c7d97c5SJed Brown    NN interface preconditioner changed to BDDC
12628eeda7d8SStefano Zampini    Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */
12630c7d97c5SJed Brown 
12640c7d97c5SJed Brown   PetscFunctionBegin;
126585c4d303SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
12660c7d97c5SJed Brown     /* First Dirichlet solve */
12670c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12680c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126953cdbc3dSStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
12700c7d97c5SJed Brown     /*
12710c7d97c5SJed Brown       Assembling right hand side for BDDC operator
1272674ae819SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
1273674ae819SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
12740c7d97c5SJed Brown     */
12750c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12760c7d97c5SJed Brown     ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
12778eeda7d8SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
12780c7d97c5SJed Brown     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
12790c7d97c5SJed Brown     ierr = VecCopy(r,z);CHKERRQ(ierr);
12800c7d97c5SJed Brown     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12810c7d97c5SJed Brown     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1282674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
1283b76ba322SStefano Zampini   } else {
12840bdf917eSStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
1285b76ba322SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
1286674ae819SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
1287b76ba322SStefano Zampini   }
1288b76ba322SStefano Zampini 
12892617d88aSStefano Zampini   /* Apply interface preconditioner
12902617d88aSStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1291dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr);
12922617d88aSStefano Zampini 
1293674ae819SStefano Zampini   /* Apply transpose of partition of unity operator */
1294674ae819SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
12950c7d97c5SJed Brown 
12963b03a366Sstefano_zampini   /* Second Dirichlet solve and assembling of output */
12970c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12980c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12990c7d97c5SJed Brown   ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
13008eeda7d8SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1301df187020SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1302df187020SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1303df187020SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1304df187020SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
13050c7d97c5SJed Brown   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13060c7d97c5SJed Brown   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13070c7d97c5SJed Brown   PetscFunctionReturn(0);
13080c7d97c5SJed Brown }
130950efa1b5SStefano Zampini 
131050efa1b5SStefano Zampini /* -------------------------------------------------------------------------- */
131150efa1b5SStefano Zampini /*
131250efa1b5SStefano Zampini    PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
131350efa1b5SStefano Zampini 
131450efa1b5SStefano Zampini    Input Parameters:
131550efa1b5SStefano Zampini .  pc - the preconditioner context
131650efa1b5SStefano Zampini .  r - input vector (global)
131750efa1b5SStefano Zampini 
131850efa1b5SStefano Zampini    Output Parameter:
131950efa1b5SStefano Zampini .  z - output vector (global)
132050efa1b5SStefano Zampini 
132150efa1b5SStefano Zampini    Application Interface Routine: PCApplyTranspose()
132250efa1b5SStefano Zampini  */
132350efa1b5SStefano Zampini #undef __FUNCT__
132450efa1b5SStefano Zampini #define __FUNCT__ "PCApplyTranspose_BDDC"
132550efa1b5SStefano Zampini PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
132650efa1b5SStefano Zampini {
132750efa1b5SStefano Zampini   PC_IS             *pcis = (PC_IS*)(pc->data);
132850efa1b5SStefano Zampini   PC_BDDC           *pcbddc = (PC_BDDC*)(pc->data);
132950efa1b5SStefano Zampini   PetscErrorCode    ierr;
133050efa1b5SStefano Zampini   const PetscScalar one = 1.0;
133150efa1b5SStefano Zampini   const PetscScalar m_one = -1.0;
133250efa1b5SStefano Zampini   const PetscScalar zero = 0.0;
133350efa1b5SStefano Zampini 
133450efa1b5SStefano Zampini   PetscFunctionBegin;
133550efa1b5SStefano Zampini   if (!pcbddc->use_exact_dirichlet_trick) {
133650efa1b5SStefano Zampini     /* First Dirichlet solve */
133750efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133850efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133950efa1b5SStefano Zampini     ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr);
134050efa1b5SStefano Zampini     /*
134150efa1b5SStefano Zampini       Assembling right hand side for BDDC operator
134250efa1b5SStefano Zampini       - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE)
134350efa1b5SStefano Zampini       - pcis->vec1_B the interface part of the global vector z
134450efa1b5SStefano Zampini     */
134550efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
134650efa1b5SStefano Zampini     ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr);
134750efa1b5SStefano Zampini     if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
134850efa1b5SStefano Zampini     ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr);
134950efa1b5SStefano Zampini     ierr = VecCopy(r,z);CHKERRQ(ierr);
135050efa1b5SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
135150efa1b5SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
135250efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr);
135350efa1b5SStefano Zampini   } else {
135450efa1b5SStefano Zampini     ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr);
135550efa1b5SStefano Zampini     ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr);
135650efa1b5SStefano Zampini     ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr);
135750efa1b5SStefano Zampini   }
135850efa1b5SStefano Zampini 
135950efa1b5SStefano Zampini   /* Apply interface preconditioner
136050efa1b5SStefano Zampini      input/output vecs: pcis->vec1_B and pcis->vec1_D */
1361dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr);
136250efa1b5SStefano Zampini 
136350efa1b5SStefano Zampini   /* Apply transpose of partition of unity operator */
136450efa1b5SStefano Zampini   ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr);
136550efa1b5SStefano Zampini 
136650efa1b5SStefano Zampini   /* Second Dirichlet solve and assembling of output */
136750efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
136850efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
136950efa1b5SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr);
137050efa1b5SStefano Zampini   if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); }
1371b0147a47SStefano Zampini   ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr);
1372b0147a47SStefano Zampini   ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr);
1373b0147a47SStefano Zampini   if (pcbddc->switch_static) { ierr = VecAXPY(pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); }
1374b0147a47SStefano Zampini   ierr = VecAXPY(pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr);
137550efa1b5SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
137650efa1b5SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
137750efa1b5SStefano Zampini   PetscFunctionReturn(0);
137850efa1b5SStefano Zampini }
1379da1bb401SStefano Zampini /* -------------------------------------------------------------------------- */
1380674ae819SStefano Zampini 
1381da1bb401SStefano Zampini #undef __FUNCT__
1382da1bb401SStefano Zampini #define __FUNCT__ "PCDestroy_BDDC"
1383da1bb401SStefano Zampini PetscErrorCode PCDestroy_BDDC(PC pc)
1384da1bb401SStefano Zampini {
1385da1bb401SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
1386da1bb401SStefano Zampini   PetscErrorCode ierr;
1387da1bb401SStefano Zampini 
1388da1bb401SStefano Zampini   PetscFunctionBegin;
1389da1bb401SStefano Zampini   /* free data created by PCIS */
1390da1bb401SStefano Zampini   ierr = PCISDestroy(pc);CHKERRQ(ierr);
1391674ae819SStefano Zampini   /* free BDDC custom data  */
1392674ae819SStefano Zampini   ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr);
1393674ae819SStefano Zampini   /* destroy objects related to topography */
1394674ae819SStefano Zampini   ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr);
1395674ae819SStefano Zampini   /* free allocated graph structure */
1396da1bb401SStefano Zampini   ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr);
1397674ae819SStefano Zampini   /* free data for scaling operator */
1398674ae819SStefano Zampini   ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr);
1399674ae819SStefano Zampini   /* free solvers stuff */
1400674ae819SStefano Zampini   ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr);
140162a6ff1dSStefano Zampini   /* free global vectors needed in presolve */
140262a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
140362a6ff1dSStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
1404*906d46d4SStefano Zampini   /* free stuff for change of basis hooks */
1405*906d46d4SStefano Zampini   if (pcbddc->new_global_mat) {
1406*906d46d4SStefano Zampini     PCBDDCChange_ctx change_ctx;
1407*906d46d4SStefano Zampini     ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr);
1408*906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr);
1409*906d46d4SStefano Zampini     ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr);
1410*906d46d4SStefano Zampini     ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr);
1411*906d46d4SStefano Zampini     ierr = PetscFree(change_ctx);CHKERRQ(ierr);
1412*906d46d4SStefano Zampini   }
1413*906d46d4SStefano Zampini   ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr);
1414*906d46d4SStefano Zampini   /* remove map from local boundary to local numbering */
141539e2fb2aSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr);
14163425bc38SStefano Zampini   /* remove functions */
1417*906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr);
1418674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr);
1419bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr);
14202b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr);
1421b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr);
14222b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr);
1423bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr);
1424bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
142582ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1426bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
142782ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1428bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr);
142982ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr);
1430bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr);
1431785d1243SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr);
1432bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr);
143363602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr);
1434bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr);
1435bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr);
1436bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr);
1437bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr);
1438674ae819SStefano Zampini   /* Free the private data structure */
1439674ae819SStefano Zampini   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1440da1bb401SStefano Zampini   PetscFunctionReturn(0);
1441da1bb401SStefano Zampini }
14423425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
14431e6b0712SBarry Smith 
14443425bc38SStefano Zampini #undef __FUNCT__
14453425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC"
14463425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
14473425bc38SStefano Zampini {
1448674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
14493425bc38SStefano Zampini   PC_IS*         pcis;
14503425bc38SStefano Zampini   PC_BDDC*       pcbddc;
14513425bc38SStefano Zampini   PetscErrorCode ierr;
14520c7d97c5SJed Brown 
14533425bc38SStefano Zampini   PetscFunctionBegin;
14543425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
14553425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
14563425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
14573425bc38SStefano Zampini 
14583425bc38SStefano Zampini   /* change of basis for physical rhs if needed
14593425bc38SStefano Zampini      It also changes the rhs in case of dirichlet boundaries */
14603308cffdSStefano Zampini   ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr);
14613425bc38SStefano Zampini   /* store vectors for computation of fetidp final solution */
14623425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14633425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1464fb223d50SStefano Zampini   /* scale rhs since it should be unassembled */
1465fb223d50SStefano Zampini   /* TODO use counter scaling? (also below) */
14663425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14673425bc38SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1468674ae819SStefano Zampini   /* Apply partition of unity */
14693425bc38SStefano Zampini   ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
1470674ae819SStefano Zampini   /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
14718eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
14723425bc38SStefano Zampini     /* compute partially subassembled Schur complement right-hand side */
14733425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14743425bc38SStefano Zampini     ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr);
14753425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr);
14763425bc38SStefano Zampini     ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr);
14773425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14783425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1479674ae819SStefano Zampini     /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */
14803425bc38SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14813425bc38SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14823425bc38SStefano Zampini     ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr);
14833425bc38SStefano Zampini   }
14843425bc38SStefano Zampini   /* BDDC rhs */
14853425bc38SStefano Zampini   ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr);
14868eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
14873425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
14883425bc38SStefano Zampini   }
14893425bc38SStefano Zampini   /* apply BDDC */
1490dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
14913425bc38SStefano Zampini   /* Application of B_delta and assembling of rhs for fetidp fluxes */
14923425bc38SStefano Zampini   ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr);
14933425bc38SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
14943425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14953425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14963425bc38SStefano Zampini   /* restore original rhs */
14973425bc38SStefano Zampini   ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr);
14983425bc38SStefano Zampini   PetscFunctionReturn(0);
14993425bc38SStefano Zampini }
15001e6b0712SBarry Smith 
15013425bc38SStefano Zampini #undef __FUNCT__
15023425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetRHS"
15033425bc38SStefano Zampini /*@
150428509bceSStefano Zampini  PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system
15053425bc38SStefano Zampini 
15063425bc38SStefano Zampini    Collective
15073425bc38SStefano Zampini 
15083425bc38SStefano Zampini    Input Parameters:
150928509bceSStefano Zampini +  fetidp_mat   - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators
151028509bceSStefano Zampini .  standard_rhs - the right-hand side for your linear system
15113425bc38SStefano Zampini 
15123425bc38SStefano Zampini    Output Parameters:
151328509bceSStefano Zampini -  fetidp_flux_rhs   - the right-hand side for the FETIDP linear system
15143425bc38SStefano Zampini 
15153425bc38SStefano Zampini    Level: developer
15163425bc38SStefano Zampini 
15173425bc38SStefano Zampini    Notes:
15183425bc38SStefano Zampini 
151928509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
15203425bc38SStefano Zampini @*/
15213425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
15223425bc38SStefano Zampini {
1523674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15243425bc38SStefano Zampini   PetscErrorCode ierr;
15253425bc38SStefano Zampini 
15263425bc38SStefano Zampini   PetscFunctionBegin;
15273425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15283425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr);
15293425bc38SStefano Zampini   PetscFunctionReturn(0);
15303425bc38SStefano Zampini }
15313425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
15321e6b0712SBarry Smith 
15333425bc38SStefano Zampini #undef __FUNCT__
15343425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC"
15353425bc38SStefano Zampini static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
15363425bc38SStefano Zampini {
1537674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15383425bc38SStefano Zampini   PC_IS*         pcis;
15393425bc38SStefano Zampini   PC_BDDC*       pcbddc;
15403425bc38SStefano Zampini   PetscErrorCode ierr;
15413425bc38SStefano Zampini 
15423425bc38SStefano Zampini   PetscFunctionBegin;
15433425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
15443425bc38SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
15453425bc38SStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
15463425bc38SStefano Zampini 
15473425bc38SStefano Zampini   /* apply B_delta^T */
15483425bc38SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15493425bc38SStefano Zampini   ierr = VecScatterEnd  (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15503425bc38SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
15513425bc38SStefano Zampini   /* compute rhs for BDDC application */
15523425bc38SStefano Zampini   ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr);
15538eeda7d8SStefano Zampini   if (pcbddc->switch_static) {
15543425bc38SStefano Zampini     ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15553425bc38SStefano Zampini   }
15563425bc38SStefano Zampini   /* apply BDDC */
1557dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
15583425bc38SStefano Zampini   /* put values into standard global vector */
15593425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15603425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15618eeda7d8SStefano Zampini   if (!pcbddc->switch_static) {
15623425bc38SStefano Zampini     /* compute values into the interior if solved for the partially subassembled Schur complement */
15633425bc38SStefano Zampini     ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr);
15643425bc38SStefano Zampini     ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr);
15653425bc38SStefano Zampini     ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr);
15663425bc38SStefano Zampini   }
15673425bc38SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15683425bc38SStefano Zampini   ierr = VecScatterEnd  (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15693425bc38SStefano Zampini   /* final change of basis if needed
15703425bc38SStefano Zampini      Is also sums the dirichlet part removed during RHS assembling */
15713308cffdSStefano Zampini   ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr);
15723425bc38SStefano Zampini   PetscFunctionReturn(0);
15733425bc38SStefano Zampini }
15741e6b0712SBarry Smith 
15753425bc38SStefano Zampini #undef __FUNCT__
15763425bc38SStefano Zampini #define __FUNCT__ "PCBDDCMatFETIDPGetSolution"
15773425bc38SStefano Zampini /*@
157828509bceSStefano Zampini  PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system
15793425bc38SStefano Zampini 
15803425bc38SStefano Zampini    Collective
15813425bc38SStefano Zampini 
15823425bc38SStefano Zampini    Input Parameters:
158328509bceSStefano Zampini +  fetidp_mat        - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators
158428509bceSStefano Zampini .  fetidp_flux_sol - the solution of the FETIDP linear system
15853425bc38SStefano Zampini 
15863425bc38SStefano Zampini    Output Parameters:
158728509bceSStefano Zampini -  standard_sol      - the solution defined on the physical domain
15883425bc38SStefano Zampini 
15893425bc38SStefano Zampini    Level: developer
15903425bc38SStefano Zampini 
15913425bc38SStefano Zampini    Notes:
15923425bc38SStefano Zampini 
159328509bceSStefano Zampini .seealso: PCBDDC,PCBDDCCreateFETIDPOperators
15943425bc38SStefano Zampini @*/
15953425bc38SStefano Zampini PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
15963425bc38SStefano Zampini {
1597674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
15983425bc38SStefano Zampini   PetscErrorCode ierr;
15993425bc38SStefano Zampini 
16003425bc38SStefano Zampini   PetscFunctionBegin;
16013425bc38SStefano Zampini   ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr);
16023425bc38SStefano Zampini   ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr);
16033425bc38SStefano Zampini   PetscFunctionReturn(0);
16043425bc38SStefano Zampini }
16053425bc38SStefano Zampini /* -------------------------------------------------------------------------- */
16061e6b0712SBarry Smith 
1607f23aa3ddSBarry Smith extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec);
1608edf7251bSStefano Zampini extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec);
1609f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat);
1610f23aa3ddSBarry Smith extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec);
1611edf7251bSStefano Zampini extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec);
1612f23aa3ddSBarry Smith extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC);
1613674ae819SStefano Zampini 
16143425bc38SStefano Zampini #undef __FUNCT__
16153425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC"
16163425bc38SStefano Zampini static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
16173425bc38SStefano Zampini {
1618674ae819SStefano Zampini 
1619674ae819SStefano Zampini   FETIDPMat_ctx  fetidpmat_ctx;
16203425bc38SStefano Zampini   Mat            newmat;
1621674ae819SStefano Zampini   FETIDPPC_ctx   fetidppc_ctx;
16223425bc38SStefano Zampini   PC             newpc;
1623ce94432eSBarry Smith   MPI_Comm       comm;
16243425bc38SStefano Zampini   PetscErrorCode ierr;
16253425bc38SStefano Zampini 
16263425bc38SStefano Zampini   PetscFunctionBegin;
1627ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
16283425bc38SStefano Zampini   /* FETIDP linear matrix */
16293425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr);
16303425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr);
16313425bc38SStefano Zampini   ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr);
16323425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr);
1633edf7251bSStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr);
16343425bc38SStefano Zampini   ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr);
16353425bc38SStefano Zampini   ierr = MatSetUp(newmat);CHKERRQ(ierr);
16363425bc38SStefano Zampini   /* FETIDP preconditioner */
16373425bc38SStefano Zampini   ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr);
16383425bc38SStefano Zampini   ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr);
16393425bc38SStefano Zampini   ierr = PCCreate(comm,&newpc);CHKERRQ(ierr);
16403425bc38SStefano Zampini   ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr);
16413425bc38SStefano Zampini   ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr);
16423425bc38SStefano Zampini   ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr);
1643edf7251bSStefano Zampini   ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr);
16443425bc38SStefano Zampini   ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr);
164523ee1639SBarry Smith   ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr);
16463425bc38SStefano Zampini   ierr = PCSetUp(newpc);CHKERRQ(ierr);
16473425bc38SStefano Zampini   /* return pointers for objects created */
16483425bc38SStefano Zampini   *fetidp_mat=newmat;
16493425bc38SStefano Zampini   *fetidp_pc=newpc;
16503425bc38SStefano Zampini   PetscFunctionReturn(0);
16513425bc38SStefano Zampini }
16521e6b0712SBarry Smith 
16533425bc38SStefano Zampini #undef __FUNCT__
16543425bc38SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPOperators"
16553425bc38SStefano Zampini /*@
165628509bceSStefano Zampini  PCBDDCCreateFETIDPOperators - Create operators for FETIDP
16573425bc38SStefano Zampini 
16583425bc38SStefano Zampini    Collective
16593425bc38SStefano Zampini 
16603425bc38SStefano Zampini    Input Parameters:
166128509bceSStefano Zampini +  pc - the BDDC preconditioning context already setup
166228509bceSStefano Zampini 
166328509bceSStefano Zampini    Output Parameters:
166428509bceSStefano Zampini .  fetidp_mat - shell FETIDP matrix object
166528509bceSStefano Zampini .  fetidp_pc  - shell Dirichlet preconditioner for FETIDP matrix
166628509bceSStefano Zampini 
166728509bceSStefano Zampini    Options Database Keys:
166828509bceSStefano Zampini -    -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers
16693425bc38SStefano Zampini 
16703425bc38SStefano Zampini    Level: developer
16713425bc38SStefano Zampini 
16723425bc38SStefano Zampini    Notes:
167328509bceSStefano Zampini      Currently the only operation provided for FETIDP matrix is MatMult
16743425bc38SStefano Zampini 
16753425bc38SStefano Zampini .seealso: PCBDDC
16763425bc38SStefano Zampini @*/
16773425bc38SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc)
16783425bc38SStefano Zampini {
16793425bc38SStefano Zampini   PetscErrorCode ierr;
16803425bc38SStefano Zampini 
16813425bc38SStefano Zampini   PetscFunctionBegin;
16823425bc38SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16833425bc38SStefano Zampini   if (pc->setupcalled) {
1684516d51a7SStefano Zampini     ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr);
1685f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n");
16863425bc38SStefano Zampini   PetscFunctionReturn(0);
16873425bc38SStefano Zampini }
16880c7d97c5SJed Brown /* -------------------------------------------------------------------------- */
1689da1bb401SStefano Zampini /*MC
1690da1bb401SStefano Zampini    PCBDDC - Balancing Domain Decomposition by Constraints.
16910c7d97c5SJed Brown 
169228509bceSStefano Zampini    An implementation of the BDDC preconditioner based on
169328509bceSStefano Zampini 
169428509bceSStefano Zampini .vb
169528509bceSStefano Zampini    [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
169628509bceSStefano 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
169728509bceSStefano Zampini    [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
169828509bceSStefano Zampini .ve
169928509bceSStefano Zampini 
170028509bceSStefano Zampini    The matrix to be preconditioned (Pmat) must be of type MATIS.
170128509bceSStefano Zampini 
1702b6fdb6dfSStefano Zampini    Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
170328509bceSStefano Zampini 
170428509bceSStefano Zampini    It also works with unsymmetric and indefinite problems.
170528509bceSStefano Zampini 
1706b6fdb6dfSStefano 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.
1707b6fdb6dfSStefano Zampini 
170828509bceSStefano Zampini    Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace
170928509bceSStefano Zampini 
17106a818285SBarry 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()
171128509bceSStefano Zampini 
17126a818285SBarry Smith    Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace().
171328509bceSStefano Zampini 
1714b6fdb6dfSStefano 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.
171528509bceSStefano Zampini 
171628509bceSStefano Zampini    The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object.
171728509bceSStefano Zampini 
1718da1bb401SStefano Zampini    Options Database Keys:
171928509bceSStefano Zampini 
172028509bceSStefano Zampini .    -pc_bddc_use_vertices <1> - use or not vertices in primal space
172128509bceSStefano Zampini .    -pc_bddc_use_edges <1> - use or not edges in primal space
1722b6fdb6dfSStefano Zampini .    -pc_bddc_use_faces <0> - use or not faces in primal space
172328509bceSStefano Zampini .    -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only)
172428509bceSStefano Zampini .    -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested
172528509bceSStefano Zampini .    -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1])
172628509bceSStefano Zampini .    -pc_bddc_levels <0> - maximum number of levels for multilevel
172728509bceSStefano Zampini .    -pc_bddc_coarsening_ratio - H/h ratio at the coarser level
172828509bceSStefano Zampini -    -pc_bddc_check_level <0> - set verbosity level of debugging output
172928509bceSStefano Zampini 
173028509bceSStefano Zampini    Options for Dirichlet, Neumann or coarse solver can be set with
173128509bceSStefano Zampini .vb
173228509bceSStefano Zampini       -pc_bddc_dirichlet_
173328509bceSStefano Zampini       -pc_bddc_neumann_
173428509bceSStefano Zampini       -pc_bddc_coarse_
173528509bceSStefano Zampini .ve
173628509bceSStefano Zampini    e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg
173728509bceSStefano Zampini 
173828509bceSStefano Zampini    When using a multilevel approach, solvers' options at the N-th level can be specified as
173928509bceSStefano Zampini .vb
1740312be037SStefano Zampini       -pc_bddc_dirichlet_lN_
1741312be037SStefano Zampini       -pc_bddc_neumann_lN_
1742312be037SStefano Zampini       -pc_bddc_coarse_lN_
174328509bceSStefano Zampini .ve
1744312be037SStefano Zampini    Note that level number ranges from the finest 0 to the coarsest N.
1745da1bb401SStefano Zampini 
1746da1bb401SStefano Zampini    Level: intermediate
1747da1bb401SStefano Zampini 
1748b6fdb6dfSStefano Zampini    Developer notes:
1749da1bb401SStefano Zampini 
175028509bceSStefano Zampini    New deluxe scaling operator will be available soon.
1751da1bb401SStefano Zampini 
1752da1bb401SStefano Zampini    Contributed by Stefano Zampini
1753da1bb401SStefano Zampini 
1754da1bb401SStefano Zampini .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
1755da1bb401SStefano Zampini M*/
1756b2573a8aSBarry Smith 
1757da1bb401SStefano Zampini #undef __FUNCT__
1758da1bb401SStefano Zampini #define __FUNCT__ "PCCreate_BDDC"
17598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
1760da1bb401SStefano Zampini {
1761da1bb401SStefano Zampini   PetscErrorCode      ierr;
1762da1bb401SStefano Zampini   PC_BDDC             *pcbddc;
1763da1bb401SStefano Zampini 
1764da1bb401SStefano Zampini   PetscFunctionBegin;
1765da1bb401SStefano Zampini   /* Creates the private data structure for this preconditioner and attach it to the PC object. */
1766b00a9115SJed Brown   ierr      = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr);
1767da1bb401SStefano Zampini   pc->data  = (void*)pcbddc;
1768da1bb401SStefano Zampini 
1769da1bb401SStefano Zampini   /* create PCIS data structure */
1770da1bb401SStefano Zampini   ierr = PCISCreate(pc);CHKERRQ(ierr);
1771da1bb401SStefano Zampini 
177247d04d0dSStefano Zampini   /* BDDC customization */
177308a5cf49SStefano Zampini   pcbddc->use_local_adj       = PETSC_TRUE;
177447d04d0dSStefano Zampini   pcbddc->use_vertices        = PETSC_TRUE;
177547d04d0dSStefano Zampini   pcbddc->use_edges           = PETSC_TRUE;
177647d04d0dSStefano Zampini   pcbddc->use_faces           = PETSC_FALSE;
177747d04d0dSStefano Zampini   pcbddc->use_change_of_basis = PETSC_FALSE;
177847d04d0dSStefano Zampini   pcbddc->use_change_on_faces = PETSC_FALSE;
177947d04d0dSStefano Zampini   pcbddc->switch_static       = PETSC_FALSE;
178047d04d0dSStefano Zampini   pcbddc->use_nnsp_true       = PETSC_FALSE; /* not yet exposed */
178147d04d0dSStefano Zampini   pcbddc->dbg_flag            = 0;
1782b9d89cd5SStefano Zampini   /* private */
1783b9d89cd5SStefano Zampini   pcbddc->issym                      = PETSC_FALSE;
178439e2fb2aSStefano Zampini   pcbddc->BtoNmap                    = 0;
1785727cdba6SStefano Zampini   pcbddc->local_primal_size          = 0;
1786e9189074SStefano Zampini   pcbddc->n_vertices                 = 0;
1787e9189074SStefano Zampini   pcbddc->n_actual_vertices          = 0;
1788e9189074SStefano Zampini   pcbddc->n_constraints              = 0;
1789727cdba6SStefano Zampini   pcbddc->primal_indices_local_idxs  = 0;
1790fb180af4SStefano Zampini   pcbddc->recompute_topography       = PETSC_FALSE;
179168457ee5SStefano Zampini   pcbddc->coarse_size                = -1;
1792f4ddd8eeSStefano Zampini   pcbddc->new_primal_space           = PETSC_FALSE;
1793727cdba6SStefano Zampini   pcbddc->new_primal_space_local     = PETSC_FALSE;
1794f4ddd8eeSStefano Zampini   pcbddc->global_primal_indices      = 0;
1795f4ddd8eeSStefano Zampini   pcbddc->onearnullspace             = 0;
1796f4ddd8eeSStefano Zampini   pcbddc->onearnullvecs_state        = 0;
1797674ae819SStefano Zampini   pcbddc->user_primal_vertices       = 0;
17980bdf917eSStefano Zampini   pcbddc->NullSpace                  = 0;
17993972b0daSStefano Zampini   pcbddc->temp_solution              = 0;
1800534831adSStefano Zampini   pcbddc->original_rhs               = 0;
1801534831adSStefano Zampini   pcbddc->local_mat                  = 0;
1802534831adSStefano Zampini   pcbddc->ChangeOfBasisMatrix        = 0;
1803b9b85e73SStefano Zampini   pcbddc->user_ChangeOfBasisMatrix   = 0;
1804*906d46d4SStefano Zampini   pcbddc->new_global_mat             = 0;
1805da1bb401SStefano Zampini   pcbddc->coarse_vec                 = 0;
1806da1bb401SStefano Zampini   pcbddc->coarse_rhs                 = 0;
1807da1bb401SStefano Zampini   pcbddc->coarse_ksp                 = 0;
1808da1bb401SStefano Zampini   pcbddc->coarse_phi_B               = 0;
1809da1bb401SStefano Zampini   pcbddc->coarse_phi_D               = 0;
181015aaf578SStefano Zampini   pcbddc->coarse_psi_B               = 0;
181115aaf578SStefano Zampini   pcbddc->coarse_psi_D               = 0;
1812da1bb401SStefano Zampini   pcbddc->vec1_P                     = 0;
1813da1bb401SStefano Zampini   pcbddc->vec1_R                     = 0;
1814da1bb401SStefano Zampini   pcbddc->vec2_R                     = 0;
1815da1bb401SStefano Zampini   pcbddc->local_auxmat1              = 0;
1816da1bb401SStefano Zampini   pcbddc->local_auxmat2              = 0;
1817da1bb401SStefano Zampini   pcbddc->R_to_B                     = 0;
1818da1bb401SStefano Zampini   pcbddc->R_to_D                     = 0;
1819da1bb401SStefano Zampini   pcbddc->ksp_D                      = 0;
1820da1bb401SStefano Zampini   pcbddc->ksp_R                      = 0;
1821da1bb401SStefano Zampini   pcbddc->NeumannBoundaries          = 0;
1822785d1243SStefano Zampini   pcbddc->NeumannBoundariesLocal     = 0;
1823785d1243SStefano Zampini   pcbddc->DirichletBoundaries        = 0;
1824785d1243SStefano Zampini   pcbddc->DirichletBoundariesLocal   = 0;
182560d44989SStefano Zampini   pcbddc->user_provided_isfordofs    = PETSC_FALSE;
182660d44989SStefano Zampini   pcbddc->n_ISForDofs                = 0;
182763602bcaSStefano Zampini   pcbddc->n_ISForDofsLocal           = 0;
1828da1bb401SStefano Zampini   pcbddc->ISForDofs                  = 0;
182963602bcaSStefano Zampini   pcbddc->ISForDofsLocal             = 0;
1830da1bb401SStefano Zampini   pcbddc->ConstraintMatrix           = 0;
183185c4d303SStefano Zampini   pcbddc->use_exact_dirichlet_trick  = PETSC_TRUE;
183247d04d0dSStefano Zampini   pcbddc->coarse_loc_to_glob         = 0;
183347d04d0dSStefano Zampini   pcbddc->coarsening_ratio           = 8;
18344fad6a16SStefano Zampini   pcbddc->current_level              = 0;
18352b510759SStefano Zampini   pcbddc->max_levels                 = 0;
1836323d395dSStefano Zampini   pcbddc->use_coarse_estimates       = PETSC_FALSE;
1837f3bde8b3SStefano Zampini   pcbddc->coarse_subassembling       = 0;
1838323d395dSStefano Zampini   pcbddc->coarse_subassembling_init  = 0;
1839da1bb401SStefano Zampini 
1840674ae819SStefano Zampini   /* create local graph structure */
1841674ae819SStefano Zampini   ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr);
1842674ae819SStefano Zampini 
1843674ae819SStefano Zampini   /* scaling */
1844674ae819SStefano Zampini   pcbddc->use_deluxe_scaling         = PETSC_FALSE;
184508a5cf49SStefano Zampini   pcbddc->deluxe_ctx                 = 0;
1846674ae819SStefano Zampini   pcbddc->work_scaling               = 0;
1847da1bb401SStefano Zampini 
1848da1bb401SStefano Zampini   /* function pointers */
1849da1bb401SStefano Zampini   pc->ops->apply               = PCApply_BDDC;
185093bd9ae7SStefano Zampini   pc->ops->applytranspose      = PCApplyTranspose_BDDC;
1851da1bb401SStefano Zampini   pc->ops->setup               = PCSetUp_BDDC;
1852da1bb401SStefano Zampini   pc->ops->destroy             = PCDestroy_BDDC;
1853da1bb401SStefano Zampini   pc->ops->setfromoptions      = PCSetFromOptions_BDDC;
1854da1bb401SStefano Zampini   pc->ops->view                = 0;
1855da1bb401SStefano Zampini   pc->ops->applyrichardson     = 0;
1856da1bb401SStefano Zampini   pc->ops->applysymmetricleft  = 0;
1857da1bb401SStefano Zampini   pc->ops->applysymmetricright = 0;
1858534831adSStefano Zampini   pc->ops->presolve            = PCPreSolve_BDDC;
1859534831adSStefano Zampini   pc->ops->postsolve           = PCPostSolve_BDDC;
1860da1bb401SStefano Zampini 
1861da1bb401SStefano Zampini   /* composing function */
1862*906d46d4SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr);
1863674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr);
1864bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr);
18652b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr);
1866b8ffe317SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr);
18672b510759SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr);
1868bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr);
1869bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr);
187082ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1871bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr);
187282ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1873bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr);
187482ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr);
1875bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr);
187682ba6b80SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr);
1877bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr);
187863602bcaSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr);
1879bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr);
1880bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr);
1881bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr);
1882bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr);
1883da1bb401SStefano Zampini   PetscFunctionReturn(0);
1884da1bb401SStefano Zampini }
18853425bc38SStefano Zampini 
1886