Actual source code: bddc.c
petsc-3.11.4 2019-09-28
1: /* TODOLIST
3: Solvers
4: - Add support for cholesky for coarse solver (similar to local solvers)
5: - Propagate ksp prefixes for solvers to mat objects?
7: User interface
8: - ** DM attached to pc?
10: Debugging output
11: - * Better management of verbosity levels of debugging output
13: Extra
14: - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)?
15: - BDDC with MG framework?
17: MATIS related operations contained in BDDC code
18: - Provide general case for subassembling
20: */
22: #include <../src/ksp/pc/impls/bddc/bddc.h>
23: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
24: #include <petscblaslapack.h>
26: static PetscBool PCBDDCPackageInitialized = PETSC_FALSE;
28: static PetscBool cited = PETSC_FALSE;
29: static const char citation[] =
30: "@article{ZampiniPCBDDC,\n"
31: "author = {Stefano Zampini},\n"
32: "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n"
33: "journal = {SIAM Journal on Scientific Computing},\n"
34: "volume = {38},\n"
35: "number = {5},\n"
36: "pages = {S282-S306},\n"
37: "year = {2016},\n"
38: "doi = {10.1137/15M1025785},\n"
39: "URL = {http://dx.doi.org/10.1137/15M1025785},\n"
40: "eprint = {http://dx.doi.org/10.1137/15M1025785}\n"
41: "}\n";
43: PetscLogEvent PC_BDDC_Topology[PETSC_PCBDDC_MAXLEVELS];
44: PetscLogEvent PC_BDDC_LocalSolvers[PETSC_PCBDDC_MAXLEVELS];
45: PetscLogEvent PC_BDDC_LocalWork[PETSC_PCBDDC_MAXLEVELS];
46: PetscLogEvent PC_BDDC_CorrectionSetUp[PETSC_PCBDDC_MAXLEVELS];
47: PetscLogEvent PC_BDDC_CoarseSetUp[PETSC_PCBDDC_MAXLEVELS];
48: PetscLogEvent PC_BDDC_CoarseSolver[PETSC_PCBDDC_MAXLEVELS];
49: PetscLogEvent PC_BDDC_AdaptiveSetUp[PETSC_PCBDDC_MAXLEVELS];
50: PetscLogEvent PC_BDDC_Scaling[PETSC_PCBDDC_MAXLEVELS];
51: PetscLogEvent PC_BDDC_Schurs[PETSC_PCBDDC_MAXLEVELS];
53: PetscErrorCode PCApply_BDDC(PC,Vec,Vec);
55: PetscErrorCode PCSetFromOptions_BDDC(PetscOptionItems *PetscOptionsObject,PC pc)
56: {
57: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
58: PetscInt nt,i;
62: PetscOptionsHead(PetscOptionsObject,"BDDC options");
63: /* Verbose debugging */
64: PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);
65: /* Approximate solvers */
66: PetscOptionsBool("-pc_bddc_dirichlet_approximate","Inform PCBDDC that we are using approximate Dirichlet solvers","none",pcbddc->NullSpace_corr[0],&pcbddc->NullSpace_corr[0],NULL);
67: PetscOptionsBool("-pc_bddc_dirichlet_approximate_scale","Inform PCBDDC that we need to scale the Dirichlet solve","none",pcbddc->NullSpace_corr[1],&pcbddc->NullSpace_corr[1],NULL);
68: PetscOptionsBool("-pc_bddc_neumann_approximate","Inform PCBDDC that we are using approximate Neumann solvers","none",pcbddc->NullSpace_corr[2],&pcbddc->NullSpace_corr[2],NULL);
69: PetscOptionsBool("-pc_bddc_neumann_approximate_scale","Inform PCBDDC that we need to scale the Neumann solve","none",pcbddc->NullSpace_corr[3],&pcbddc->NullSpace_corr[3],NULL);
70: /* Primal space customization */
71: 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);
72: PetscOptionsInt("-pc_bddc_graph_maxcount","Maximum number of shared subdomains for a connected component","none",pcbddc->graphmaxcount,&pcbddc->graphmaxcount,NULL);
73: PetscOptionsBool("-pc_bddc_corner_selection","Activates face-based corner selection","none",pcbddc->corner_selection,&pcbddc->corner_selection,NULL);
74: PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);
75: PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);
76: PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);
77: PetscOptionsInt("-pc_bddc_vertex_size","Connected components smaller or equal to vertex size will be considered as primal vertices","none",pcbddc->vertex_size,&pcbddc->vertex_size,NULL);
78: PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);
79: PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is always used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);
80: /* Change of basis */
81: 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);
82: 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);
83: if (!pcbddc->use_change_of_basis) {
84: pcbddc->use_change_on_faces = PETSC_FALSE;
85: }
86: /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */
87: PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);
88: PetscOptionsInt("-pc_bddc_coarse_eqs_per_proc","Target number of equations per process for coarse problem redistribution (significant only at the coarsest level)","none",pcbddc->coarse_eqs_per_proc,&pcbddc->coarse_eqs_per_proc,NULL);
89: i = pcbddc->coarsening_ratio;
90: PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","PCBDDCSetCoarseningRatio",i,&i,NULL);
91: PCBDDCSetCoarseningRatio(pc,i);
92: i = pcbddc->max_levels;
93: PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","PCBDDCSetLevels",i,&i,NULL);
94: PCBDDCSetLevels(pc,i);
95: PetscOptionsInt("-pc_bddc_coarse_eqs_limit","Set maximum number of equations on coarsest grid to aim for","none",pcbddc->coarse_eqs_limit,&pcbddc->coarse_eqs_limit,NULL);
96: PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);
97: PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);
98: PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);
99: PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);
100: PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);
101: PetscOptionsBool("-pc_bddc_schur_exact","Whether or not to use the exact Schur complement instead of the reduced one (which excludes size 1 cc)","none",pcbddc->sub_schurs_exact_schur,&pcbddc->sub_schurs_exact_schur,NULL);
102: PetscOptionsBool("-pc_bddc_deluxe_zerorows","Zero rows and columns of deluxe operators associated with primal dofs","none",pcbddc->deluxe_zerorows,&pcbddc->deluxe_zerorows,NULL);
103: PetscOptionsBool("-pc_bddc_deluxe_singlemat","Collapse deluxe operators","none",pcbddc->deluxe_singlemat,&pcbddc->deluxe_singlemat,NULL);
104: PetscOptionsBool("-pc_bddc_adaptive_userdefined","Use user-defined constraints (should be attached via MatSetNearNullSpace to pmat) in addition to those adaptively generated","none",pcbddc->adaptive_userdefined,&pcbddc->adaptive_userdefined,NULL);
105: nt = 2;
106: PetscOptionsRealArray("-pc_bddc_adaptive_threshold","Thresholds to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&nt,NULL);
107: if (nt == 1) pcbddc->adaptive_threshold[1] = pcbddc->adaptive_threshold[0];
108: PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);
109: PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);
110: PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);
111: PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);
112: PetscOptionsBool("-pc_bddc_benign_trick","Apply the benign subspace trick to saddle point problems with discontinuous pressures","none",pcbddc->benign_saddle_point,&pcbddc->benign_saddle_point,NULL);
113: PetscOptionsBool("-pc_bddc_benign_change","Compute the pressure change of basis explicitly","none",pcbddc->benign_change_explicit,&pcbddc->benign_change_explicit,NULL);
114: PetscOptionsBool("-pc_bddc_benign_compute_correction","Compute the benign correction during PreSolve","none",pcbddc->benign_compute_correction,&pcbddc->benign_compute_correction,NULL);
115: PetscOptionsBool("-pc_bddc_nonetflux","Automatic computation of no-net-flux quadrature weights","none",pcbddc->compute_nonetflux,&pcbddc->compute_nonetflux,NULL);
116: PetscOptionsBool("-pc_bddc_detect_disconnected","Detects disconnected subdomains","none",pcbddc->detect_disconnected,&pcbddc->detect_disconnected,NULL);
117: PetscOptionsBool("-pc_bddc_detect_disconnected_filter","Filters out small entries in the local matrix when detecting disconnected subdomains","none",pcbddc->detect_disconnected_filter,&pcbddc->detect_disconnected_filter,NULL);
118: PetscOptionsBool("-pc_bddc_eliminate_dirichlet","Whether or not we want to eliminate dirichlet dofs during presolve","none",pcbddc->eliminate_dirdofs,&pcbddc->eliminate_dirdofs,NULL);
119: PetscOptionsTail();
120: return(0);
121: }
123: static PetscErrorCode PCView_BDDC(PC pc,PetscViewer viewer)
124: {
125: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
126: PC_IS *pcis = (PC_IS*)pc->data;
127: PetscErrorCode ierr;
128: PetscBool isascii;
129: PetscSubcomm subcomm;
130: PetscViewer subviewer;
133: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
134: /* ASCII viewer */
135: if (isascii) {
136: PetscMPIInt color,rank,size;
137: PetscInt64 loc[7],gsum[6],gmax[6],gmin[6],totbenign;
138: PetscScalar interface_size;
139: PetscReal ratio1=0.,ratio2=0.;
140: Vec counter;
142: if (!pc->setupcalled) {
143: PetscViewerASCIIPrintf(viewer," Partial information available: preconditioner has not been setup yet\n");
144: }
145: PetscViewerASCIIPrintf(viewer," Use verbose output: %D\n",pcbddc->dbg_flag);
146: PetscViewerASCIIPrintf(viewer," Use user-defined CSR: %d\n",!!pcbddc->mat_graph->nvtxs_csr);
147: PetscViewerASCIIPrintf(viewer," Use local mat graph: %d\n",pcbddc->use_local_adj && !pcbddc->mat_graph->nvtxs_csr);
148: if (pcbddc->mat_graph->twodim) {
149: PetscViewerASCIIPrintf(viewer," Connectivity graph topological dimension: 2\n");
150: } else {
151: PetscViewerASCIIPrintf(viewer," Connectivity graph topological dimension: 3\n");
152: }
153: if (pcbddc->graphmaxcount != PETSC_MAX_INT) {
154: PetscViewerASCIIPrintf(viewer," Graph max count: %D\n",pcbddc->graphmaxcount);
155: }
156: PetscViewerASCIIPrintf(viewer," Use vertices: %d (vertex size %D)\n",pcbddc->use_vertices,pcbddc->vertex_size);
157: PetscViewerASCIIPrintf(viewer," Use edges: %d\n",pcbddc->use_edges);
158: PetscViewerASCIIPrintf(viewer," Use faces: %d\n",pcbddc->use_faces);
159: PetscViewerASCIIPrintf(viewer," Use true near null space: %d\n",pcbddc->use_nnsp_true);
160: PetscViewerASCIIPrintf(viewer," Use QR for single constraints on cc: %d\n",pcbddc->use_qr_single);
161: PetscViewerASCIIPrintf(viewer," Use change of basis on local edge nodes: %d\n",pcbddc->use_change_of_basis);
162: PetscViewerASCIIPrintf(viewer," Use change of basis on local face nodes: %d\n",pcbddc->use_change_on_faces);
163: PetscViewerASCIIPrintf(viewer," User defined change of basis matrix: %d\n",!!pcbddc->user_ChangeOfBasisMatrix);
164: PetscViewerASCIIPrintf(viewer," Has change of basis matrix: %d\n",!!pcbddc->ChangeOfBasisMatrix);
165: PetscViewerASCIIPrintf(viewer," Eliminate dirichlet boundary dofs: %d\n",pcbddc->eliminate_dirdofs);
166: PetscViewerASCIIPrintf(viewer," Switch on static condensation ops around the interface preconditioner: %d\n",pcbddc->switch_static);
167: PetscViewerASCIIPrintf(viewer," Use exact dirichlet trick: %d\n",pcbddc->use_exact_dirichlet_trick);
168: PetscViewerASCIIPrintf(viewer," Multilevel max levels: %D\n",pcbddc->max_levels);
169: PetscViewerASCIIPrintf(viewer," Multilevel coarsening ratio: %D\n",pcbddc->coarsening_ratio);
170: PetscViewerASCIIPrintf(viewer," Use estimated eigs for coarse problem: %d\n",pcbddc->use_coarse_estimates);
171: PetscViewerASCIIPrintf(viewer," Use deluxe scaling: %d\n",pcbddc->use_deluxe_scaling);
172: PetscViewerASCIIPrintf(viewer," Use deluxe zerorows: %d\n",pcbddc->deluxe_zerorows);
173: PetscViewerASCIIPrintf(viewer," Use deluxe singlemat: %d\n",pcbddc->deluxe_singlemat);
174: PetscViewerASCIIPrintf(viewer," Rebuild interface graph for Schur principal minors: %d\n",pcbddc->sub_schurs_rebuild);
175: PetscViewerASCIIPrintf(viewer," Number of dofs' layers for the computation of principal minors: %D\n",pcbddc->sub_schurs_layers);
176: PetscViewerASCIIPrintf(viewer," Use user CSR graph to compute successive layers: %d\n",pcbddc->sub_schurs_use_useradj);
177: if (pcbddc->adaptive_threshold[1] != pcbddc->adaptive_threshold[0]) {
178: PetscViewerASCIIPrintf(viewer," Adaptive constraint selection thresholds (active %d, userdefined %d): %g,%g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0],pcbddc->adaptive_threshold[1]);
179: } else {
180: PetscViewerASCIIPrintf(viewer," Adaptive constraint selection threshold (active %d, userdefined %d): %g\n",pcbddc->adaptive_selection,pcbddc->adaptive_userdefined,pcbddc->adaptive_threshold[0]);
181: }
182: PetscViewerASCIIPrintf(viewer," Min constraints / connected component: %D\n",pcbddc->adaptive_nmin);
183: PetscViewerASCIIPrintf(viewer," Max constraints / connected component: %D\n",pcbddc->adaptive_nmax);
184: PetscViewerASCIIPrintf(viewer," Invert exact Schur complement for adaptive selection: %d\n",pcbddc->sub_schurs_exact_schur);
185: PetscViewerASCIIPrintf(viewer," Symmetric computation of primal basis functions: %d\n",pcbddc->symmetric_primal);
186: PetscViewerASCIIPrintf(viewer," Num. Procs. to map coarse adjacency list: %D\n",pcbddc->coarse_adj_red);
187: PetscViewerASCIIPrintf(viewer," Coarse eqs per proc (significant at the coarsest level): %D\n",pcbddc->coarse_eqs_per_proc);
188: PetscViewerASCIIPrintf(viewer," Detect disconnected: %d (filter %d)\n",pcbddc->detect_disconnected,pcbddc->detect_disconnected_filter);
189: PetscViewerASCIIPrintf(viewer," Benign subspace trick: %d (change explicit %d)\n",pcbddc->benign_saddle_point,pcbddc->benign_change_explicit);
190: PetscViewerASCIIPrintf(viewer," Benign subspace trick is active: %d\n",pcbddc->benign_have_null);
191: PetscViewerASCIIPrintf(viewer," Algebraic computation of no-net-flux: %d\n",pcbddc->compute_nonetflux);
192: if (!pc->setupcalled) return(0);
194: /* compute interface size */
195: VecSet(pcis->vec1_B,1.0);
196: MatCreateVecs(pc->pmat,&counter,0);
197: VecSet(counter,0.0);
198: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);
199: VecScatterEnd(pcis->global_to_B,pcis->vec1_B,counter,INSERT_VALUES,SCATTER_REVERSE);
200: VecSum(counter,&interface_size);
201: VecDestroy(&counter);
203: /* compute some statistics on the domain decomposition */
204: gsum[0] = 1;
205: gsum[1] = gsum[2] = gsum[3] = gsum[4] = gsum[5] = 0;
206: loc[0] = !!pcis->n;
207: loc[1] = pcis->n - pcis->n_B;
208: loc[2] = pcis->n_B;
209: loc[3] = pcbddc->local_primal_size;
210: loc[4] = pcis->n;
211: loc[5] = pcbddc->n_local_subs > 0 ? pcbddc->n_local_subs : (pcis->n ? 1 : 0);
212: loc[6] = pcbddc->benign_n;
213: MPI_Reduce(loc,gsum,6,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));
214: if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = -1;
215: MPI_Reduce(loc,gmax,6,MPIU_INT64,MPI_MAX,0,PetscObjectComm((PetscObject)pc));
216: if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = PETSC_MAX_INT;
217: MPI_Reduce(loc,gmin,6,MPIU_INT64,MPI_MIN,0,PetscObjectComm((PetscObject)pc));
218: MPI_Reduce(&loc[6],&totbenign,1,MPIU_INT64,MPI_SUM,0,PetscObjectComm((PetscObject)pc));
219: if (pcbddc->coarse_size) {
220: ratio1 = pc->pmat->rmap->N/(1.*pcbddc->coarse_size);
221: ratio2 = PetscRealPart(interface_size)/pcbddc->coarse_size;
222: }
223: PetscViewerASCIIPrintf(viewer,"********************************** STATISTICS AT LEVEL %d **********************************\n",pcbddc->current_level);
224: PetscViewerASCIIPrintf(viewer," Global dofs sizes: all %D interface %D coarse %D\n",pc->pmat->rmap->N,(PetscInt)PetscRealPart(interface_size),pcbddc->coarse_size);
225: PetscViewerASCIIPrintf(viewer," Coarsening ratios: all/coarse %D interface/coarse %D\n",(PetscInt)ratio1,(PetscInt)ratio2);
226: PetscViewerASCIIPrintf(viewer," Active processes : %D\n",(PetscInt)gsum[0]);
227: PetscViewerASCIIPrintf(viewer," Total subdomains : %D\n",(PetscInt)gsum[5]);
228: if (pcbddc->benign_have_null) {
229: PetscViewerASCIIPrintf(viewer," Benign subs : %D\n",(PetscInt)totbenign);
230: }
231: PetscViewerASCIIPrintf(viewer," Dofs type :\tMIN\tMAX\tMEAN\n");
232: PetscViewerASCIIPrintf(viewer," Interior dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[1],(PetscInt)gmax[1],(PetscInt)(gsum[1]/gsum[0]));
233: PetscViewerASCIIPrintf(viewer," Interface dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[2],(PetscInt)gmax[2],(PetscInt)(gsum[2]/gsum[0]));
234: PetscViewerASCIIPrintf(viewer," Primal dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[3],(PetscInt)gmax[3],(PetscInt)(gsum[3]/gsum[0]));
235: PetscViewerASCIIPrintf(viewer," Local dofs :\t%D\t%D\t%D\n",(PetscInt)gmin[4],(PetscInt)gmax[4],(PetscInt)(gsum[4]/gsum[0]));
236: PetscViewerASCIIPrintf(viewer," Local subs :\t%D\t%D\n" ,(PetscInt)gmin[5],(PetscInt)gmax[5]);
237: PetscViewerFlush(viewer);
239: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
241: /* local solvers */
242: PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pcbddc->ksp_D),&subviewer);
243: if (!rank) {
244: PetscViewerASCIIPrintf(subviewer,"--- Interior solver (rank 0)\n");
245: PetscViewerASCIIPushTab(subviewer);
246: KSPView(pcbddc->ksp_D,subviewer);
247: PetscViewerASCIIPopTab(subviewer);
248: PetscViewerASCIIPrintf(subviewer,"--- Correction solver (rank 0)\n");
249: PetscViewerASCIIPushTab(subviewer);
250: KSPView(pcbddc->ksp_R,subviewer);
251: PetscViewerASCIIPopTab(subviewer);
252: PetscViewerFlush(subviewer);
253: }
254: PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pcbddc->ksp_D),&subviewer);
255: PetscViewerFlush(viewer);
257: /* the coarse problem can be handled by a different communicator */
258: if (pcbddc->coarse_ksp) color = 1;
259: else color = 0;
260: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
261: PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&subcomm);
262: PetscSubcommSetNumber(subcomm,PetscMin(size,2));
263: PetscSubcommSetTypeGeneral(subcomm,color,rank);
264: PetscViewerGetSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);
265: if (color == 1) {
266: PetscViewerASCIIPrintf(subviewer,"--- Coarse solver\n");
267: PetscViewerASCIIPushTab(subviewer);
268: KSPView(pcbddc->coarse_ksp,subviewer);
269: PetscViewerASCIIPopTab(subviewer);
270: PetscViewerFlush(subviewer);
271: }
272: PetscViewerRestoreSubViewer(viewer,PetscSubcommChild(subcomm),&subviewer);
273: PetscSubcommDestroy(&subcomm);
274: PetscViewerFlush(viewer);
275: }
276: return(0);
277: }
279: static PetscErrorCode PCBDDCSetDiscreteGradient_BDDC(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming)
280: {
281: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
285: PetscObjectReference((PetscObject)G);
286: MatDestroy(&pcbddc->discretegradient);
287: pcbddc->discretegradient = G;
288: pcbddc->nedorder = order > 0 ? order : -order;
289: pcbddc->nedfield = field;
290: pcbddc->nedglobal = global;
291: pcbddc->conforming = conforming;
292: return(0);
293: }
295: /*@
296: PCBDDCSetDiscreteGradient - Sets the discrete gradient
298: Collective on PC
300: Input Parameters:
301: + pc - the preconditioning context
302: . G - the discrete gradient matrix (should be in AIJ format)
303: . order - the order of the Nedelec space (1 for the lowest order)
304: . field - the field id of the Nedelec dofs (not used if the fields have not been specified)
305: . global - the type of global ordering for the rows of G
306: - conforming - whether the mesh is conforming or not
308: Level: advanced
310: Notes:
311: The discrete gradient matrix G is used to analyze the subdomain edges, and it should not contain any zero entry.
312: For variable order spaces, the order should be set to zero.
313: If global is true, the rows of G should be given in global ordering for the whole dofs;
314: if false, the ordering should be global for the Nedelec field.
315: In the latter case, it should hold gid[i] < gid[j] iff geid[i] < geid[j], with gid the global orderding for all the dofs
316: and geid the one for the Nedelec field.
318: .seealso: PCBDDC,PCBDDCSetDofsSplitting(),PCBDDCSetDofsSplittingLocal()
319: @*/
320: PetscErrorCode PCBDDCSetDiscreteGradient(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming)
321: {
332: PetscTryMethod(pc,"PCBDDCSetDiscreteGradient_C",(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool),(pc,G,order,field,global,conforming));
333: return(0);
334: }
336: static PetscErrorCode PCBDDCSetDivergenceMat_BDDC(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
337: {
338: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
342: PetscObjectReference((PetscObject)divudotp);
343: MatDestroy(&pcbddc->divudotp);
344: pcbddc->divudotp = divudotp;
345: pcbddc->divudotp_trans = trans;
346: pcbddc->compute_nonetflux = PETSC_TRUE;
347: if (vl2l) {
348: PetscObjectReference((PetscObject)vl2l);
349: ISDestroy(&pcbddc->divudotp_vl2l);
350: pcbddc->divudotp_vl2l = vl2l;
351: }
352: return(0);
353: }
355: /*@
356: PCBDDCSetDivergenceMat - Sets the linear operator representing \int_\Omega \div {\bf u} \cdot p dx
358: Collective on PC
360: Input Parameters:
361: + pc - the preconditioning context
362: . divudotp - the matrix (must be of type MATIS)
363: . trans - if trans if false (resp. true), then pressures are in the test (trial) space and velocities are in the trial (test) space.
364: - vl2l - optional index set describing the local (wrt the local matrix in divudotp) to local (wrt the local matrix in the preconditioning matrix) map for the velocities
366: Level: advanced
368: Notes:
369: This auxiliary matrix is used to compute quadrature weights representing the net-flux across subdomain boundaries
370: If vl2l is NULL, the local ordering for velocities in divudotp should match that of the preconditioning matrix
372: .seealso: PCBDDC
373: @*/
374: PetscErrorCode PCBDDCSetDivergenceMat(PC pc, Mat divudotp, PetscBool trans, IS vl2l)
375: {
376: PetscBool ismatis;
385: PetscObjectTypeCompare((PetscObject)divudotp,MATIS,&ismatis);
386: if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Divergence matrix needs to be of type MATIS");
387: PetscTryMethod(pc,"PCBDDCSetDivergenceMat_C",(PC,Mat,PetscBool,IS),(pc,divudotp,trans,vl2l));
388: return(0);
389: }
391: static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change, PetscBool interior)
392: {
393: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
397: PetscObjectReference((PetscObject)change);
398: MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);
399: pcbddc->user_ChangeOfBasisMatrix = change;
400: pcbddc->change_interior = interior;
401: return(0);
402: }
403: /*@
404: PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs
406: Collective on PC
408: Input Parameters:
409: + pc - the preconditioning context
410: . change - the change of basis matrix
411: - interior - whether or not the change of basis modifies interior dofs
413: Level: intermediate
415: Notes:
417: .seealso: PCBDDC
418: @*/
419: PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change, PetscBool interior)
420: {
427: if (pc->mat) {
428: PetscInt rows_c,cols_c,rows,cols;
429: MatGetSize(pc->mat,&rows,&cols);
430: MatGetSize(change,&rows_c,&cols_c);
431: if (rows_c != rows) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %D != %D",rows_c,rows);
432: if (cols_c != cols) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %D != %D",cols_c,cols);
433: MatGetLocalSize(pc->mat,&rows,&cols);
434: MatGetLocalSize(change,&rows_c,&cols_c);
435: if (rows_c != rows) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %D != %D",rows_c,rows);
436: if (cols_c != cols) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %D != %D",cols_c,cols);
437: }
438: PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat,PetscBool),(pc,change,interior));
439: return(0);
440: }
442: static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices)
443: {
444: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
445: PetscBool isequal = PETSC_FALSE;
449: PetscObjectReference((PetscObject)PrimalVertices);
450: if (pcbddc->user_primal_vertices) {
451: ISEqual(PrimalVertices,pcbddc->user_primal_vertices,&isequal);
452: }
453: ISDestroy(&pcbddc->user_primal_vertices);
454: ISDestroy(&pcbddc->user_primal_vertices_local);
455: pcbddc->user_primal_vertices = PrimalVertices;
456: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
457: return(0);
458: }
460: /*@
461: PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in PCBDDC
463: Collective
465: Input Parameters:
466: + pc - the preconditioning context
467: - PrimalVertices - index set of primal vertices in global numbering (can be empty)
469: Level: intermediate
471: Notes:
472: Any process can list any global node
474: .seealso: PCBDDC, PCBDDCGetPrimalVerticesIS(), PCBDDCSetPrimalVerticesLocalIS(), PCBDDCGetPrimalVerticesLocalIS()
475: @*/
476: PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices)
477: {
484: PetscTryMethod(pc,"PCBDDCSetPrimalVerticesIS_C",(PC,IS),(pc,PrimalVertices));
485: return(0);
486: }
488: static PetscErrorCode PCBDDCGetPrimalVerticesIS_BDDC(PC pc, IS *is)
489: {
490: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
493: *is = pcbddc->user_primal_vertices;
494: return(0);
495: }
497: /*@
498: PCBDDCGetPrimalVerticesIS - Get user defined primal vertices set with PCBDDCSetPrimalVerticesIS()
500: Collective
502: Input Parameters:
503: . pc - the preconditioning context
505: Output Parameters:
506: . is - index set of primal vertices in global numbering (NULL if not set)
508: Level: intermediate
510: Notes:
512: .seealso: PCBDDC, PCBDDCSetPrimalVerticesIS(), PCBDDCSetPrimalVerticesLocalIS(), PCBDDCGetPrimalVerticesLocalIS()
513: @*/
514: PetscErrorCode PCBDDCGetPrimalVerticesIS(PC pc, IS *is)
515: {
521: PetscUseMethod(pc,"PCBDDCGetPrimalVerticesIS_C",(PC,IS*),(pc,is));
522: return(0);
523: }
525: static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices)
526: {
527: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
528: PetscBool isequal = PETSC_FALSE;
532: PetscObjectReference((PetscObject)PrimalVertices);
533: if (pcbddc->user_primal_vertices_local) {
534: ISEqual(PrimalVertices,pcbddc->user_primal_vertices_local,&isequal);
535: }
536: ISDestroy(&pcbddc->user_primal_vertices);
537: ISDestroy(&pcbddc->user_primal_vertices_local);
538: pcbddc->user_primal_vertices_local = PrimalVertices;
539: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
540: return(0);
541: }
543: /*@
544: PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC
546: Collective
548: Input Parameters:
549: + pc - the preconditioning context
550: - PrimalVertices - index set of primal vertices in local numbering (can be empty)
552: Level: intermediate
554: Notes:
556: .seealso: PCBDDC, PCBDDCSetPrimalVerticesIS(), PCBDDCGetPrimalVerticesIS(), PCBDDCGetPrimalVerticesLocalIS()
557: @*/
558: PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices)
559: {
566: PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));
567: return(0);
568: }
570: static PetscErrorCode PCBDDCGetPrimalVerticesLocalIS_BDDC(PC pc, IS *is)
571: {
572: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
575: *is = pcbddc->user_primal_vertices_local;
576: return(0);
577: }
579: /*@
580: PCBDDCGetPrimalVerticesLocalIS - Get user defined primal vertices set with PCBDDCSetPrimalVerticesLocalIS()
582: Collective
584: Input Parameters:
585: . pc - the preconditioning context
587: Output Parameters:
588: . is - index set of primal vertices in local numbering (NULL if not set)
590: Level: intermediate
592: Notes:
594: .seealso: PCBDDC, PCBDDCSetPrimalVerticesIS(), PCBDDCGetPrimalVerticesIS(), PCBDDCSetPrimalVerticesLocalIS()
595: @*/
596: PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC pc, IS *is)
597: {
603: PetscUseMethod(pc,"PCBDDCGetPrimalVerticesLocalIS_C",(PC,IS*),(pc,is));
604: return(0);
605: }
607: static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k)
608: {
609: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
612: pcbddc->coarsening_ratio = k;
613: return(0);
614: }
616: /*@
617: PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel
619: Logically collective on PC
621: Input Parameters:
622: + pc - the preconditioning context
623: - k - coarsening ratio (H/h at the coarser level)
625: Options Database Keys:
626: . -pc_bddc_coarsening_ratio
628: Level: intermediate
630: Notes:
631: Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level
633: .seealso: PCBDDC, PCBDDCSetLevels()
634: @*/
635: PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k)
636: {
642: PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));
643: return(0);
644: }
646: /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */
647: static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg)
648: {
649: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
652: pcbddc->use_exact_dirichlet_trick = flg;
653: return(0);
654: }
656: PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg)
657: {
663: PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));
664: return(0);
665: }
667: static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level)
668: {
669: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
672: pcbddc->current_level = level;
673: return(0);
674: }
676: PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level)
677: {
683: PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));
684: return(0);
685: }
687: static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels)
688: {
689: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
692: if (levels > PETSC_PCBDDC_MAXLEVELS-1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of additional levels for BDDC is %d",PETSC_PCBDDC_MAXLEVELS-1);
693: pcbddc->max_levels = levels;
694: return(0);
695: }
697: /*@
698: PCBDDCSetLevels - Sets the maximum number of additional levels allowed for multilevel BDDC
700: Logically collective on PC
702: Input Parameters:
703: + pc - the preconditioning context
704: - levels - the maximum number of levels
706: Options Database Keys:
707: . -pc_bddc_levels
709: Level: intermediate
711: Notes:
712: The default value is 0, that gives the classical two-levels BDDC
714: .seealso: PCBDDC, PCBDDCSetCoarseningRatio()
715: @*/
716: PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels)
717: {
723: PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));
724: return(0);
725: }
727: static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries)
728: {
729: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
730: PetscBool isequal = PETSC_FALSE;
734: PetscObjectReference((PetscObject)DirichletBoundaries);
735: if (pcbddc->DirichletBoundaries) {
736: ISEqual(DirichletBoundaries,pcbddc->DirichletBoundaries,&isequal);
737: }
738: /* last user setting takes precendence -> destroy any other customization */
739: ISDestroy(&pcbddc->DirichletBoundariesLocal);
740: ISDestroy(&pcbddc->DirichletBoundaries);
741: pcbddc->DirichletBoundaries = DirichletBoundaries;
742: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
743: return(0);
744: }
746: /*@
747: PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem.
749: Collective
751: Input Parameters:
752: + pc - the preconditioning context
753: - DirichletBoundaries - parallel IS defining the Dirichlet boundaries
755: Level: intermediate
757: Notes:
758: Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node
760: .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns()
761: @*/
762: PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries)
763: {
770: PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));
771: return(0);
772: }
774: static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries)
775: {
776: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
777: PetscBool isequal = PETSC_FALSE;
781: PetscObjectReference((PetscObject)DirichletBoundaries);
782: if (pcbddc->DirichletBoundariesLocal) {
783: ISEqual(DirichletBoundaries,pcbddc->DirichletBoundariesLocal,&isequal);
784: }
785: /* last user setting takes precendence -> destroy any other customization */
786: ISDestroy(&pcbddc->DirichletBoundariesLocal);
787: ISDestroy(&pcbddc->DirichletBoundaries);
788: pcbddc->DirichletBoundariesLocal = DirichletBoundaries;
789: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
790: return(0);
791: }
793: /*@
794: PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering.
796: Collective
798: Input Parameters:
799: + pc - the preconditioning context
800: - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering)
802: Level: intermediate
804: Notes:
806: .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns()
807: @*/
808: PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries)
809: {
816: PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));
817: return(0);
818: }
820: static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries)
821: {
822: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
823: PetscBool isequal = PETSC_FALSE;
827: PetscObjectReference((PetscObject)NeumannBoundaries);
828: if (pcbddc->NeumannBoundaries) {
829: ISEqual(NeumannBoundaries,pcbddc->NeumannBoundaries,&isequal);
830: }
831: /* last user setting takes precendence -> destroy any other customization */
832: ISDestroy(&pcbddc->NeumannBoundariesLocal);
833: ISDestroy(&pcbddc->NeumannBoundaries);
834: pcbddc->NeumannBoundaries = NeumannBoundaries;
835: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
836: return(0);
837: }
839: /*@
840: PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem.
842: Collective
844: Input Parameters:
845: + pc - the preconditioning context
846: - NeumannBoundaries - parallel IS defining the Neumann boundaries
848: Level: intermediate
850: Notes:
851: Any process can list any global node
853: .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal()
854: @*/
855: PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries)
856: {
863: PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));
864: return(0);
865: }
867: static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries)
868: {
869: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
870: PetscBool isequal = PETSC_FALSE;
874: PetscObjectReference((PetscObject)NeumannBoundaries);
875: if (pcbddc->NeumannBoundariesLocal) {
876: ISEqual(NeumannBoundaries,pcbddc->NeumannBoundariesLocal,&isequal);
877: }
878: /* last user setting takes precendence -> destroy any other customization */
879: ISDestroy(&pcbddc->NeumannBoundariesLocal);
880: ISDestroy(&pcbddc->NeumannBoundaries);
881: pcbddc->NeumannBoundariesLocal = NeumannBoundaries;
882: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
883: return(0);
884: }
886: /*@
887: PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering.
889: Collective
891: Input Parameters:
892: + pc - the preconditioning context
893: - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering)
895: Level: intermediate
897: Notes:
899: .seealso: PCBDDC, PCBDDCSetNeumannBoundaries()
900: @*/
901: PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries)
902: {
909: PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));
910: return(0);
911: }
913: static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries)
914: {
915: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
918: *DirichletBoundaries = pcbddc->DirichletBoundaries;
919: return(0);
920: }
922: /*@
923: PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries
925: Collective
927: Input Parameters:
928: . pc - the preconditioning context
930: Output Parameters:
931: . DirichletBoundaries - index set defining the Dirichlet boundaries
933: Level: intermediate
935: Notes:
936: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries
938: .seealso: PCBDDC
939: @*/
940: PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries)
941: {
946: PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));
947: return(0);
948: }
950: static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries)
951: {
952: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
955: *DirichletBoundaries = pcbddc->DirichletBoundariesLocal;
956: return(0);
957: }
959: /*@
960: PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering)
962: Collective
964: Input Parameters:
965: . pc - the preconditioning context
967: Output Parameters:
968: . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries
970: Level: intermediate
972: Notes:
973: The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries).
974: In the latter case, the IS will be available after PCSetUp.
976: .seealso: PCBDDC
977: @*/
978: PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries)
979: {
984: PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));
985: return(0);
986: }
988: static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries)
989: {
990: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
993: *NeumannBoundaries = pcbddc->NeumannBoundaries;
994: return(0);
995: }
997: /*@
998: PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries
1000: Collective
1002: Input Parameters:
1003: . pc - the preconditioning context
1005: Output Parameters:
1006: . NeumannBoundaries - index set defining the Neumann boundaries
1008: Level: intermediate
1010: Notes:
1011: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries
1013: .seealso: PCBDDC
1014: @*/
1015: PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries)
1016: {
1021: PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));
1022: return(0);
1023: }
1025: static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries)
1026: {
1027: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1030: *NeumannBoundaries = pcbddc->NeumannBoundariesLocal;
1031: return(0);
1032: }
1034: /*@
1035: PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering)
1037: Collective
1039: Input Parameters:
1040: . pc - the preconditioning context
1042: Output Parameters:
1043: . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries
1045: Level: intermediate
1047: Notes:
1048: The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries).
1049: In the latter case, the IS will be available after PCSetUp.
1051: .seealso: PCBDDC
1052: @*/
1053: PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries)
1054: {
1059: PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));
1060: return(0);
1061: }
1063: static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
1064: {
1065: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1066: PCBDDCGraph mat_graph = pcbddc->mat_graph;
1067: PetscBool same_data = PETSC_FALSE;
1071: if (!nvtxs) {
1072: if (copymode == PETSC_OWN_POINTER) {
1073: PetscFree(xadj);
1074: PetscFree(adjncy);
1075: }
1076: PCBDDCGraphResetCSR(mat_graph);
1077: return(0);
1078: }
1079: if (mat_graph->nvtxs == nvtxs && mat_graph->freecsr) { /* we own the data */
1080: if (mat_graph->xadj == xadj && mat_graph->adjncy == adjncy) same_data = PETSC_TRUE;
1081: if (!same_data && mat_graph->xadj[nvtxs] == xadj[nvtxs]) {
1082: PetscMemcmp(xadj,mat_graph->xadj,(nvtxs+1)*sizeof(PetscInt),&same_data);
1083: if (same_data) {
1084: PetscMemcmp(adjncy,mat_graph->adjncy,xadj[nvtxs]*sizeof(PetscInt),&same_data);
1085: }
1086: }
1087: }
1088: if (!same_data) {
1089: /* free old CSR */
1090: PCBDDCGraphResetCSR(mat_graph);
1091: /* get CSR into graph structure */
1092: if (copymode == PETSC_COPY_VALUES) {
1093: PetscMalloc1(nvtxs+1,&mat_graph->xadj);
1094: PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);
1095: PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));
1096: PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));
1097: mat_graph->freecsr = PETSC_TRUE;
1098: } else if (copymode == PETSC_OWN_POINTER) {
1099: mat_graph->xadj = (PetscInt*)xadj;
1100: mat_graph->adjncy = (PetscInt*)adjncy;
1101: mat_graph->freecsr = PETSC_TRUE;
1102: } else if (copymode == PETSC_USE_POINTER) {
1103: mat_graph->xadj = (PetscInt*)xadj;
1104: mat_graph->adjncy = (PetscInt*)adjncy;
1105: mat_graph->freecsr = PETSC_FALSE;
1106: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %D",copymode);
1107: mat_graph->nvtxs_csr = nvtxs;
1108: pcbddc->recompute_topography = PETSC_TRUE;
1109: }
1110: return(0);
1111: }
1113: /*@
1114: PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local degrees of freedom.
1116: Not collective
1118: Input Parameters:
1119: + pc - the preconditioning context.
1120: . nvtxs - number of local vertices of the graph (i.e., the number of local dofs).
1121: . xadj, adjncy - the connectivity of the dofs in CSR format.
1122: - copymode - supported modes are PETSC_COPY_VALUES, PETSC_USE_POINTER or PETSC_OWN_POINTER.
1124: Level: intermediate
1126: Notes:
1127: A dof is considered connected with all local dofs if xadj[dof+1]-xadj[dof] == 1 and adjncy[xadj[dof]] is negative.
1129: .seealso: PCBDDC,PetscCopyMode
1130: @*/
1131: PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode)
1132: {
1133: void (*f)(void) = 0;
1138: if (nvtxs) {
1141: }
1142: PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));
1143: /* free arrays if PCBDDC is not the PC type */
1144: PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);
1145: if (!f && copymode == PETSC_OWN_POINTER) {
1146: PetscFree(xadj);
1147: PetscFree(adjncy);
1148: }
1149: return(0);
1150: }
1152: static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
1153: {
1154: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1155: PetscInt i;
1156: PetscBool isequal = PETSC_FALSE;
1160: if (pcbddc->n_ISForDofsLocal == n_is) {
1161: for (i=0;i<n_is;i++) {
1162: PetscBool isequalt;
1163: ISEqual(ISForDofs[i],pcbddc->ISForDofsLocal[i],&isequalt);
1164: if (!isequalt) break;
1165: }
1166: if (i == n_is) isequal = PETSC_TRUE;
1167: }
1168: for (i=0;i<n_is;i++) {
1169: PetscObjectReference((PetscObject)ISForDofs[i]);
1170: }
1171: /* Destroy ISes if they were already set */
1172: for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
1173: ISDestroy(&pcbddc->ISForDofsLocal[i]);
1174: }
1175: PetscFree(pcbddc->ISForDofsLocal);
1176: /* last user setting takes precendence -> destroy any other customization */
1177: for (i=0;i<pcbddc->n_ISForDofs;i++) {
1178: ISDestroy(&pcbddc->ISForDofs[i]);
1179: }
1180: PetscFree(pcbddc->ISForDofs);
1181: pcbddc->n_ISForDofs = 0;
1182: /* allocate space then set */
1183: if (n_is) {
1184: PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);
1185: }
1186: for (i=0;i<n_is;i++) {
1187: pcbddc->ISForDofsLocal[i] = ISForDofs[i];
1188: }
1189: pcbddc->n_ISForDofsLocal = n_is;
1190: if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
1191: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
1192: return(0);
1193: }
1195: /*@
1196: PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix
1198: Collective
1200: Input Parameters:
1201: + pc - the preconditioning context
1202: . n_is - number of index sets defining the fields
1203: - ISForDofs - array of IS describing the fields in local ordering
1205: Level: intermediate
1207: Notes:
1208: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
1210: .seealso: PCBDDC
1211: @*/
1212: PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[])
1213: {
1214: PetscInt i;
1220: for (i=0;i<n_is;i++) {
1223: }
1224: PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
1225: return(0);
1226: }
1228: static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[])
1229: {
1230: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1231: PetscInt i;
1232: PetscBool isequal = PETSC_FALSE;
1236: if (pcbddc->n_ISForDofs == n_is) {
1237: for (i=0;i<n_is;i++) {
1238: PetscBool isequalt;
1239: ISEqual(ISForDofs[i],pcbddc->ISForDofs[i],&isequalt);
1240: if (!isequalt) break;
1241: }
1242: if (i == n_is) isequal = PETSC_TRUE;
1243: }
1244: for (i=0;i<n_is;i++) {
1245: PetscObjectReference((PetscObject)ISForDofs[i]);
1246: }
1247: /* Destroy ISes if they were already set */
1248: for (i=0;i<pcbddc->n_ISForDofs;i++) {
1249: ISDestroy(&pcbddc->ISForDofs[i]);
1250: }
1251: PetscFree(pcbddc->ISForDofs);
1252: /* last user setting takes precendence -> destroy any other customization */
1253: for (i=0;i<pcbddc->n_ISForDofsLocal;i++) {
1254: ISDestroy(&pcbddc->ISForDofsLocal[i]);
1255: }
1256: PetscFree(pcbddc->ISForDofsLocal);
1257: pcbddc->n_ISForDofsLocal = 0;
1258: /* allocate space then set */
1259: if (n_is) {
1260: PetscMalloc1(n_is,&pcbddc->ISForDofs);
1261: }
1262: for (i=0;i<n_is;i++) {
1263: pcbddc->ISForDofs[i] = ISForDofs[i];
1264: }
1265: pcbddc->n_ISForDofs = n_is;
1266: if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE;
1267: if (!isequal) pcbddc->recompute_topography = PETSC_TRUE;
1268: return(0);
1269: }
1271: /*@
1272: PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix
1274: Collective
1276: Input Parameters:
1277: + pc - the preconditioning context
1278: . n_is - number of index sets defining the fields
1279: - ISForDofs - array of IS describing the fields in global ordering
1281: Level: intermediate
1283: Notes:
1284: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field.
1286: .seealso: PCBDDC
1287: @*/
1288: PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[])
1289: {
1290: PetscInt i;
1296: for (i=0;i<n_is;i++) {
1299: }
1300: PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));
1301: return(0);
1302: }
1304: /*
1305: PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial
1306: guess if a transformation of basis approach has been selected.
1308: Input Parameter:
1309: + pc - the preconditioner contex
1311: Application Interface Routine: PCPreSolve()
1313: Notes:
1314: The interface routine PCPreSolve() is not usually called directly by
1315: the user, but instead is called by KSPSolve().
1316: */
1317: static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1318: {
1320: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1321: PC_IS *pcis = (PC_IS*)(pc->data);
1322: Vec used_vec;
1323: PetscBool iscg = PETSC_FALSE, save_rhs = PETSC_TRUE, benign_correction_computed;
1326: /* if we are working with CG, one dirichlet solve can be avoided during Krylov iterations */
1327: if (ksp) {
1328: PetscBool isgroppcg, ispipecg, ispipelcg, ispipecgrr;
1330: PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);
1331: PetscObjectTypeCompare((PetscObject)ksp,KSPGROPPCG,&isgroppcg);
1332: PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECG,&ispipecg);
1333: PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECG,&ispipelcg);
1334: PetscObjectTypeCompare((PetscObject)ksp,KSPPIPECGRR,&ispipecgrr);
1335: iscg = (PetscBool)(iscg || isgroppcg || ispipecg || ispipelcg || ispipecgrr);
1336: if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || !iscg || pc->mat != pc->pmat) {
1337: PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
1338: }
1339: }
1340: if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || pc->mat != pc->pmat) {
1341: PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);
1342: }
1344: /* Creates parallel work vectors used in presolve */
1345: if (!pcbddc->original_rhs) {
1346: VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
1347: }
1348: if (!pcbddc->temp_solution) {
1349: VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);
1350: }
1352: pcbddc->temp_solution_used = PETSC_FALSE;
1353: if (x) {
1354: PetscObjectReference((PetscObject)x);
1355: used_vec = x;
1356: } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */
1357: PetscObjectReference((PetscObject)pcbddc->temp_solution);
1358: used_vec = pcbddc->temp_solution;
1359: VecSet(used_vec,0.0);
1360: pcbddc->temp_solution_used = PETSC_TRUE;
1361: VecCopy(rhs,pcbddc->original_rhs);
1362: save_rhs = PETSC_FALSE;
1363: pcbddc->eliminate_dirdofs = PETSC_TRUE;
1364: }
1366: /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */
1367: if (ksp) {
1368: /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */
1369: KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);
1370: if (!pcbddc->ksp_guess_nonzero) {
1371: VecSet(used_vec,0.0);
1372: }
1373: }
1375: pcbddc->rhs_change = PETSC_FALSE;
1376: /* Take into account zeroed rows -> change rhs and store solution removed */
1377: if (rhs && pcbddc->eliminate_dirdofs) {
1378: IS dirIS = NULL;
1380: /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */
1381: PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);
1382: if (dirIS) {
1383: Mat_IS *matis = (Mat_IS*)pc->pmat->data;
1384: PetscInt dirsize,i,*is_indices;
1385: PetscScalar *array_x;
1386: const PetscScalar *array_diagonal;
1388: MatGetDiagonal(pc->pmat,pcis->vec1_global);
1389: VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);
1390: VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1391: VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);
1392: VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1393: VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
1394: ISGetLocalSize(dirIS,&dirsize);
1395: VecGetArray(pcis->vec1_N,&array_x);
1396: VecGetArrayRead(pcis->vec2_N,&array_diagonal);
1397: ISGetIndices(dirIS,(const PetscInt**)&is_indices);
1398: for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]];
1399: ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);
1400: VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);
1401: VecRestoreArray(pcis->vec1_N,&array_x);
1402: VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1403: VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);
1404: pcbddc->rhs_change = PETSC_TRUE;
1405: ISDestroy(&dirIS);
1406: }
1407: }
1409: /* remove the computed solution or the initial guess from the rhs */
1410: if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) {
1411: /* save the original rhs */
1412: if (save_rhs) {
1413: VecSwap(rhs,pcbddc->original_rhs);
1414: save_rhs = PETSC_FALSE;
1415: }
1416: pcbddc->rhs_change = PETSC_TRUE;
1417: VecScale(used_vec,-1.0);
1418: MatMultAdd(pc->mat,used_vec,pcbddc->original_rhs,rhs);
1419: VecScale(used_vec,-1.0);
1420: VecCopy(used_vec,pcbddc->temp_solution);
1421: pcbddc->temp_solution_used = PETSC_TRUE;
1422: if (ksp) {
1423: KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);
1424: }
1425: }
1426: VecDestroy(&used_vec);
1428: /* compute initial vector in benign space if needed
1429: and remove non-benign solution from the rhs */
1430: benign_correction_computed = PETSC_FALSE;
1431: if (rhs && pcbddc->benign_compute_correction && (pcbddc->benign_have_null || pcbddc->benign_apply_coarse_only)) {
1432: /* compute u^*_h using ideas similar to those in Xuemin Tu's PhD thesis (see Section 4.8.1)
1433: Recursively apply BDDC in the multilevel case */
1434: if (!pcbddc->benign_vec) {
1435: VecDuplicate(rhs,&pcbddc->benign_vec);
1436: }
1437: /* keep applying coarse solver unless we no longer have benign subdomains */
1438: pcbddc->benign_apply_coarse_only = pcbddc->benign_have_null ? PETSC_TRUE : PETSC_FALSE;
1439: if (!pcbddc->benign_skip_correction) {
1440: PCApply_BDDC(pc,rhs,pcbddc->benign_vec);
1441: benign_correction_computed = PETSC_TRUE;
1442: if (pcbddc->temp_solution_used) {
1443: VecAXPY(pcbddc->temp_solution,1.0,pcbddc->benign_vec);
1444: }
1445: VecScale(pcbddc->benign_vec,-1.0);
1446: /* store the original rhs if not done earlier */
1447: if (save_rhs) {
1448: VecSwap(rhs,pcbddc->original_rhs);
1449: }
1450: if (pcbddc->rhs_change) {
1451: MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);
1452: } else {
1453: MatMultAdd(pc->mat,pcbddc->benign_vec,pcbddc->original_rhs,rhs);
1454: }
1455: pcbddc->rhs_change = PETSC_TRUE;
1456: }
1457: pcbddc->benign_apply_coarse_only = PETSC_FALSE;
1458: } else {
1459: VecDestroy(&pcbddc->benign_vec);
1460: }
1462: /* dbg output */
1463: if (pcbddc->dbg_flag && benign_correction_computed) {
1464: Vec v;
1466: VecDuplicate(pcis->vec1_global,&v);
1467: if (pcbddc->ChangeOfBasisMatrix) {
1468: MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,v);
1469: } else {
1470: VecCopy(rhs,v);
1471: }
1472: PCBDDCBenignGetOrSetP0(pc,v,PETSC_TRUE);
1473: PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"LEVEL %D: is the correction benign?\n",pcbddc->current_level);
1474: PetscScalarView(pcbddc->benign_n,pcbddc->benign_p0,pcbddc->dbg_viewer);
1475: PetscViewerFlush(pcbddc->dbg_viewer);
1476: VecDestroy(&v);
1477: }
1479: /* set initial guess if using PCG */
1480: pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1481: if (x && pcbddc->use_exact_dirichlet_trick) {
1482: VecSet(x,0.0);
1483: if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) {
1484: if (benign_correction_computed) { /* we have already saved the changed rhs */
1485: VecLockReadPop(pcis->vec1_global);
1486: } else {
1487: MatMultTranspose(pcbddc->ChangeOfBasisMatrix,rhs,pcis->vec1_global);
1488: }
1489: VecScatterBegin(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1490: VecScatterEnd(pcis->global_to_D,pcis->vec1_global,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1491: } else {
1492: VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1493: VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1494: }
1495: KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1496: KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec2_D);
1497: if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) {
1498: VecSet(pcis->vec1_global,0.);
1499: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);
1500: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);
1501: MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_global,x);
1502: } else {
1503: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1504: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);
1505: }
1506: if (ksp) {
1507: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
1508: }
1509: pcbddc->exact_dirichlet_trick_app = PETSC_TRUE;
1510: } else if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior && benign_correction_computed && pcbddc->use_exact_dirichlet_trick) {
1511: VecLockReadPop(pcis->vec1_global);
1512: }
1513: return(0);
1514: }
1516: /*
1517: PCPostSolve_BDDC - Changes the computed solution if a transformation of basis
1518: approach has been selected. Also, restores rhs to its original state.
1520: Input Parameter:
1521: + pc - the preconditioner contex
1523: Application Interface Routine: PCPostSolve()
1525: Notes:
1526: The interface routine PCPostSolve() is not usually called directly by
1527: the user, but instead is called by KSPSolve().
1528: */
1529: static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x)
1530: {
1532: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
1535: /* add solution removed in presolve */
1536: if (x && pcbddc->rhs_change) {
1537: if (pcbddc->temp_solution_used) {
1538: VecAXPY(x,1.0,pcbddc->temp_solution);
1539: } else if (pcbddc->benign_compute_correction && pcbddc->benign_vec) {
1540: VecAXPY(x,-1.0,pcbddc->benign_vec);
1541: }
1542: /* restore to original state (not for FETI-DP) */
1543: if (ksp) pcbddc->temp_solution_used = PETSC_FALSE;
1544: }
1546: /* restore rhs to its original state (not needed for FETI-DP) */
1547: if (rhs && pcbddc->rhs_change) {
1548: VecSwap(rhs,pcbddc->original_rhs);
1549: pcbddc->rhs_change = PETSC_FALSE;
1550: }
1551: /* restore ksp guess state */
1552: if (ksp) {
1553: KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);
1554: /* reset flag for exact dirichlet trick */
1555: pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1556: }
1557: return(0);
1558: }
1560: /*
1561: PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner
1562: by setting data structures and options.
1564: Input Parameter:
1565: + pc - the preconditioner context
1567: Application Interface Routine: PCSetUp()
1569: Notes:
1570: The interface routine PCSetUp() is not usually called directly by
1571: the user, but instead is called by PCApply() if necessary.
1572: */
1573: PetscErrorCode PCSetUp_BDDC(PC pc)
1574: {
1575: PC_BDDC* pcbddc = (PC_BDDC*)pc->data;
1576: PCBDDCSubSchurs sub_schurs;
1577: Mat_IS* matis;
1578: MatNullSpace nearnullspace;
1579: Mat lA;
1580: IS lP,zerodiag = NULL;
1581: PetscInt nrows,ncols;
1582: PetscMPIInt size;
1583: PetscBool computesubschurs;
1584: PetscBool computeconstraintsmatrix;
1585: PetscBool new_nearnullspace_provided,ismatis,rl;
1586: PetscErrorCode ierr;
1589: PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);
1590: if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS");
1591: MatGetSize(pc->pmat,&nrows,&ncols);
1592: if (nrows != ncols) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix");
1593: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
1595: matis = (Mat_IS*)pc->pmat->data;
1596: /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */
1597: /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup
1598: Also, BDDC builds its own KSP for the Dirichlet problem */
1599: rl = pcbddc->recompute_topography;
1600: if (!pc->setupcalled || pc->flag == DIFFERENT_NONZERO_PATTERN) rl = PETSC_TRUE;
1601: MPIU_Allreduce(&rl,&pcbddc->recompute_topography,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));
1602: if (pcbddc->recompute_topography) {
1603: pcbddc->graphanalyzed = PETSC_FALSE;
1604: computeconstraintsmatrix = PETSC_TRUE;
1605: } else {
1606: computeconstraintsmatrix = PETSC_FALSE;
1607: }
1609: /* check parameters' compatibility */
1610: if (!pcbddc->use_deluxe_scaling) pcbddc->deluxe_zerorows = PETSC_FALSE;
1611: pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold[0] != 0.0 || pcbddc->adaptive_threshold[1] != 0.0);
1612: pcbddc->use_deluxe_scaling = (PetscBool)(pcbddc->use_deluxe_scaling && size > 1);
1613: pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_selection && size > 1);
1614: pcbddc->adaptive_userdefined = (PetscBool)(pcbddc->adaptive_selection && pcbddc->adaptive_userdefined);
1615: if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE;
1617: computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling);
1619: /* activate all connected components if the netflux has been requested */
1620: if (pcbddc->compute_nonetflux) {
1621: pcbddc->use_vertices = PETSC_TRUE;
1622: pcbddc->use_edges = PETSC_TRUE;
1623: pcbddc->use_faces = PETSC_TRUE;
1624: }
1626: /* Get stdout for dbg */
1627: if (pcbddc->dbg_flag) {
1628: if (!pcbddc->dbg_viewer) {
1629: pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc));
1630: }
1631: PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer);
1632: PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1633: }
1635: /* process topology information */
1636: PetscLogEventBegin(PC_BDDC_Topology[pcbddc->current_level],pc,0,0,0);
1637: if (pcbddc->recompute_topography) {
1638: PCBDDCComputeLocalTopologyInfo(pc);
1639: if (pcbddc->discretegradient) {
1640: PCBDDCNedelecSupport(pc);
1641: }
1642: }
1643: if (pcbddc->corner_selected) pcbddc->use_vertices = PETSC_TRUE;
1645: /* change basis if requested by the user */
1646: if (pcbddc->user_ChangeOfBasisMatrix) {
1647: /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */
1648: pcbddc->use_change_of_basis = PETSC_FALSE;
1649: PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);
1650: } else {
1651: MatDestroy(&pcbddc->local_mat);
1652: PetscObjectReference((PetscObject)matis->A);
1653: pcbddc->local_mat = matis->A;
1654: }
1656: /*
1657: Compute change of basis on local pressures (aka zerodiag dofs) with the benign trick
1658: This should come earlier then PCISSetUp for extracting the correct subdomain matrices
1659: */
1660: PCBDDCBenignShellMat(pc,PETSC_TRUE);
1661: if (pcbddc->benign_saddle_point) {
1662: PC_IS* pcis = (PC_IS*)pc->data;
1664: if (pcbddc->user_ChangeOfBasisMatrix || pcbddc->use_change_of_basis || !computesubschurs) pcbddc->benign_change_explicit = PETSC_TRUE;
1665: /* detect local saddle point and change the basis in pcbddc->local_mat */
1666: PCBDDCBenignDetectSaddlePoint(pc,(PetscBool)(!pcbddc->recompute_topography),&zerodiag);
1667: /* pop B0 mat from local mat */
1668: PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);
1669: /* give pcis a hint to not reuse submatrices during PCISCreate */
1670: if (pc->flag == SAME_NONZERO_PATTERN && pcis->reusesubmatrices == PETSC_TRUE) {
1671: if (pcbddc->benign_n && (pcbddc->benign_change_explicit || pcbddc->dbg_flag)) {
1672: pcis->reusesubmatrices = PETSC_FALSE;
1673: } else {
1674: pcis->reusesubmatrices = PETSC_TRUE;
1675: }
1676: } else {
1677: pcis->reusesubmatrices = PETSC_FALSE;
1678: }
1679: }
1681: /* propagate relevant information */
1682: if (matis->A->symmetric_set) {
1683: MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);
1684: }
1685: if (matis->A->spd_set) {
1686: MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);
1687: }
1689: /* Set up all the "iterative substructuring" common block without computing solvers */
1690: {
1691: Mat temp_mat;
1693: temp_mat = matis->A;
1694: matis->A = pcbddc->local_mat;
1695: PCISSetUp(pc,PETSC_TRUE,PETSC_FALSE);
1696: pcbddc->local_mat = matis->A;
1697: matis->A = temp_mat;
1698: }
1700: /* Analyze interface */
1701: if (!pcbddc->graphanalyzed) {
1702: PCBDDCAnalyzeInterface(pc);
1703: computeconstraintsmatrix = PETSC_TRUE;
1704: if (pcbddc->adaptive_selection && !pcbddc->use_deluxe_scaling && !pcbddc->mat_graph->twodim) {
1705: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute the adaptive primal space for a problem with 3D edges without deluxe scaling");
1706: }
1707: if (pcbddc->compute_nonetflux) {
1708: MatNullSpace nnfnnsp;
1710: if (!pcbddc->divudotp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Missing divudotp operator");
1711: PCBDDCComputeNoNetFlux(pc->pmat,pcbddc->divudotp,pcbddc->divudotp_trans,pcbddc->divudotp_vl2l,pcbddc->mat_graph,&nnfnnsp);
1712: /* TODO what if a nearnullspace is already attached? */
1713: if (nnfnnsp) {
1714: MatSetNearNullSpace(pc->pmat,nnfnnsp);
1715: MatNullSpaceDestroy(&nnfnnsp);
1716: }
1717: }
1718: }
1719: PetscLogEventEnd(PC_BDDC_Topology[pcbddc->current_level],pc,0,0,0);
1721: /* check existence of a divergence free extension, i.e.
1722: b(v_I,p_0) = 0 for all v_I (raise error if not).
1723: Also, check that PCBDDCBenignGetOrSetP0 works */
1724: if (pcbddc->benign_saddle_point && pcbddc->dbg_flag > 1) {
1725: PCBDDCBenignCheck(pc,zerodiag);
1726: }
1727: ISDestroy(&zerodiag);
1729: /* Setup local dirichlet solver ksp_D and sub_schurs solvers */
1730: if (computesubschurs && pcbddc->recompute_topography) {
1731: PCBDDCInitSubSchurs(pc);
1732: }
1733: /* SetUp Scaling operator (scaling matrices could be needed in SubSchursSetUp)*/
1734: if (!pcbddc->use_deluxe_scaling) {
1735: PCBDDCScalingSetUp(pc);
1736: }
1738: /* finish setup solvers and do adaptive selection of constraints */
1739: sub_schurs = pcbddc->sub_schurs;
1740: if (sub_schurs && sub_schurs->schur_explicit) {
1741: if (computesubschurs) {
1742: PCBDDCSetUpSubSchurs(pc);
1743: }
1744: PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1745: } else {
1746: PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);
1747: if (computesubschurs) {
1748: PCBDDCSetUpSubSchurs(pc);
1749: }
1750: }
1751: if (pcbddc->adaptive_selection) {
1752: PCBDDCAdaptiveSelection(pc);
1753: computeconstraintsmatrix = PETSC_TRUE;
1754: }
1756: /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */
1757: new_nearnullspace_provided = PETSC_FALSE;
1758: MatGetNearNullSpace(pc->pmat,&nearnullspace);
1759: if (pcbddc->onearnullspace) { /* already used nearnullspace */
1760: if (!nearnullspace) { /* near null space attached to mat has been destroyed */
1761: new_nearnullspace_provided = PETSC_TRUE;
1762: } else {
1763: /* determine if the two nullspaces are different (should be lightweight) */
1764: if (nearnullspace != pcbddc->onearnullspace) {
1765: new_nearnullspace_provided = PETSC_TRUE;
1766: } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */
1767: PetscInt i;
1768: const Vec *nearnullvecs;
1769: PetscObjectState state;
1770: PetscInt nnsp_size;
1771: MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);
1772: for (i=0;i<nnsp_size;i++) {
1773: PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);
1774: if (pcbddc->onearnullvecs_state[i] != state) {
1775: new_nearnullspace_provided = PETSC_TRUE;
1776: break;
1777: }
1778: }
1779: }
1780: }
1781: } else {
1782: if (!nearnullspace) { /* both nearnullspaces are null */
1783: new_nearnullspace_provided = PETSC_FALSE;
1784: } else { /* nearnullspace attached later */
1785: new_nearnullspace_provided = PETSC_TRUE;
1786: }
1787: }
1789: /* Setup constraints and related work vectors */
1790: /* reset primal space flags */
1791: PetscLogEventBegin(PC_BDDC_LocalWork[pcbddc->current_level],pc,0,0,0);
1792: pcbddc->new_primal_space = PETSC_FALSE;
1793: pcbddc->new_primal_space_local = PETSC_FALSE;
1794: if (computeconstraintsmatrix || new_nearnullspace_provided) {
1795: /* It also sets the primal space flags */
1796: PCBDDCConstraintsSetUp(pc);
1797: }
1798: /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */
1799: PCBDDCSetUpLocalWorkVectors(pc);
1801: if (pcbddc->use_change_of_basis) {
1802: PC_IS *pcis = (PC_IS*)(pc->data);
1804: PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);
1805: if (pcbddc->benign_change) {
1806: MatDestroy(&pcbddc->benign_B0);
1807: /* pop B0 from pcbddc->local_mat */
1808: PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);
1809: }
1810: /* get submatrices */
1811: MatDestroy(&pcis->A_IB);
1812: MatDestroy(&pcis->A_BI);
1813: MatDestroy(&pcis->A_BB);
1814: MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);
1815: MatCreateSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);
1816: MatCreateSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);
1817: /* set flag in pcis to not reuse submatrices during PCISCreate */
1818: pcis->reusesubmatrices = PETSC_FALSE;
1819: } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) {
1820: MatDestroy(&pcbddc->local_mat);
1821: PetscObjectReference((PetscObject)matis->A);
1822: pcbddc->local_mat = matis->A;
1823: }
1825: /* interface pressure block row for B_C */
1826: PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lP" ,(PetscObject*)&lP);
1827: PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_lA" ,(PetscObject*)&lA);
1828: if (lA && lP) {
1829: PC_IS* pcis = (PC_IS*)pc->data;
1830: Mat B_BI,B_BB,Bt_BI,Bt_BB;
1831: PetscBool issym;
1832: MatIsSymmetric(lA,PETSC_SMALL,&issym);
1833: if (issym) {
1834: MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI);
1835: MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB);
1836: MatCreateTranspose(B_BI,&Bt_BI);
1837: MatCreateTranspose(B_BB,&Bt_BB);
1838: } else {
1839: MatCreateSubMatrix(lA,lP,pcis->is_I_local,MAT_INITIAL_MATRIX,&B_BI);
1840: MatCreateSubMatrix(lA,lP,pcis->is_B_local,MAT_INITIAL_MATRIX,&B_BB);
1841: MatCreateSubMatrix(lA,pcis->is_I_local,lP,MAT_INITIAL_MATRIX,&Bt_BI);
1842: MatCreateSubMatrix(lA,pcis->is_B_local,lP,MAT_INITIAL_MATRIX,&Bt_BB);
1843: }
1844: PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BI",(PetscObject)B_BI);
1845: PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_B_BB",(PetscObject)B_BB);
1846: PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BI",(PetscObject)Bt_BI);
1847: PetscObjectCompose((PetscObject)pc,"__KSPFETIDP_Bt_BB",(PetscObject)Bt_BB);
1848: MatDestroy(&B_BI);
1849: MatDestroy(&B_BB);
1850: MatDestroy(&Bt_BI);
1851: MatDestroy(&Bt_BB);
1852: }
1853: PetscLogEventEnd(PC_BDDC_LocalWork[pcbddc->current_level],pc,0,0,0);
1855: /* SetUp coarse and local Neumann solvers */
1856: PCBDDCSetUpSolvers(pc);
1857: /* SetUp Scaling operator */
1858: if (pcbddc->use_deluxe_scaling) {
1859: PCBDDCScalingSetUp(pc);
1860: }
1862: /* mark topography as done */
1863: pcbddc->recompute_topography = PETSC_FALSE;
1865: /* wrap pcis->A_IB and pcis->A_BI if we did not change explicitly the variables on the pressures */
1866: PCBDDCBenignShellMat(pc,PETSC_FALSE);
1868: if (pcbddc->dbg_flag) {
1869: PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);
1870: PetscViewerASCIIPopSynchronized(pcbddc->dbg_viewer);
1871: }
1872: return(0);
1873: }
1875: /*
1876: PCApply_BDDC - Applies the BDDC operator to a vector.
1878: Input Parameters:
1879: + pc - the preconditioner context
1880: - r - input vector (global)
1882: Output Parameter:
1883: . z - output vector (global)
1885: Application Interface Routine: PCApply()
1886: */
1887: PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z)
1888: {
1889: PC_IS *pcis = (PC_IS*)(pc->data);
1890: PC_BDDC *pcbddc = (PC_BDDC*)(pc->data);
1891: Mat lA = NULL;
1892: PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B;
1893: PetscErrorCode ierr;
1894: const PetscScalar one = 1.0;
1895: const PetscScalar m_one = -1.0;
1896: const PetscScalar zero = 0.0;
1897: /* This code is similar to that provided in nn.c for PCNN
1898: NN interface preconditioner changed to BDDC
1899: Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */
1902: PetscCitationsRegister(citation,&cited);
1903: if (pcbddc->switch_static) {
1904: MatISGetLocalMat(pc->useAmat ? pc->mat : pc->pmat,&lA);
1905: }
1907: if (pcbddc->ChangeOfBasisMatrix) {
1908: Vec swap;
1910: MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change);
1911: swap = pcbddc->work_change;
1912: pcbddc->work_change = r;
1913: r = swap;
1914: /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */
1915: if (pcbddc->benign_apply_coarse_only && pcbddc->use_exact_dirichlet_trick && pcbddc->change_interior) {
1916: VecCopy(r,pcis->vec1_global);
1917: VecLockReadPush(pcis->vec1_global);
1918: }
1919: }
1920: if (pcbddc->benign_have_null) { /* get p0 from r */
1921: PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);
1922: }
1923: if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
1924: VecCopy(r,z);
1925: /* First Dirichlet solve */
1926: VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1927: VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
1928: /*
1929: Assembling right hand side for BDDC operator
1930: - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
1931: - pcis->vec1_B the interface part of the global vector z
1932: */
1933: if (n_D) {
1934: KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
1935: KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec2_D);
1936: VecScale(pcis->vec2_D,m_one);
1937: if (pcbddc->switch_static) {
1938: VecSet(pcis->vec1_N,0.);
1939: VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1940: VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1941: if (!pcbddc->switch_static_change) {
1942: MatMult(lA,pcis->vec1_N,pcis->vec2_N);
1943: } else {
1944: MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1945: MatMult(lA,pcis->vec2_N,pcis->vec1_N);
1946: MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1947: }
1948: VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
1949: VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
1950: VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1951: VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1952: } else {
1953: MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);
1954: }
1955: } else {
1956: VecSet(pcis->vec1_B,zero);
1957: }
1958: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1959: VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
1960: PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
1961: } else {
1962: if (!pcbddc->benign_apply_coarse_only) {
1963: PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
1964: }
1965: }
1967: /* Apply interface preconditioner
1968: input/output vecs: pcis->vec1_B and pcis->vec1_D */
1969: PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);
1971: /* Apply transpose of partition of unity operator */
1972: PCBDDCScalingExtension(pc,pcis->vec1_B,z);
1974: /* Second Dirichlet solve and assembling of output */
1975: VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1976: VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
1977: if (n_B) {
1978: if (pcbddc->switch_static) {
1979: VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1980: VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1981: VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1982: VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
1983: if (!pcbddc->switch_static_change) {
1984: MatMult(lA,pcis->vec1_N,pcis->vec2_N);
1985: } else {
1986: MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1987: MatMult(lA,pcis->vec2_N,pcis->vec1_N);
1988: MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
1989: }
1990: VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
1991: VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
1992: } else {
1993: MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);
1994: }
1995: } else if (pcbddc->switch_static) { /* n_B is zero */
1996: if (!pcbddc->switch_static_change) {
1997: MatMult(lA,pcis->vec1_D,pcis->vec3_D);
1998: } else {
1999: MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);
2000: MatMult(lA,pcis->vec1_N,pcis->vec2_N);
2001: MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);
2002: }
2003: }
2004: KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
2005: KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec4_D);
2007: if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
2008: if (pcbddc->switch_static) {
2009: VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
2010: } else {
2011: VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
2012: }
2013: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
2014: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
2015: } else {
2016: if (pcbddc->switch_static) {
2017: VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
2018: } else {
2019: VecScale(pcis->vec4_D,m_one);
2020: }
2021: VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
2022: VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
2023: }
2024: if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */
2025: if (pcbddc->benign_apply_coarse_only) {
2026: PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
2027: }
2028: PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);
2029: }
2031: if (pcbddc->ChangeOfBasisMatrix) {
2032: pcbddc->work_change = r;
2033: VecCopy(z,pcbddc->work_change);
2034: MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);
2035: }
2036: return(0);
2037: }
2039: /*
2040: PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector.
2042: Input Parameters:
2043: + pc - the preconditioner context
2044: - r - input vector (global)
2046: Output Parameter:
2047: . z - output vector (global)
2049: Application Interface Routine: PCApplyTranspose()
2050: */
2051: PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z)
2052: {
2053: PC_IS *pcis = (PC_IS*)(pc->data);
2054: PC_BDDC *pcbddc = (PC_BDDC*)(pc->data);
2055: Mat lA = NULL;
2056: PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B;
2057: PetscErrorCode ierr;
2058: const PetscScalar one = 1.0;
2059: const PetscScalar m_one = -1.0;
2060: const PetscScalar zero = 0.0;
2063: PetscCitationsRegister(citation,&cited);
2064: if (pcbddc->switch_static) {
2065: MatISGetLocalMat(pc->useAmat ? pc->mat : pc->pmat,&lA);
2066: }
2067: if (pcbddc->ChangeOfBasisMatrix) {
2068: Vec swap;
2070: MatMultTranspose(pcbddc->ChangeOfBasisMatrix,r,pcbddc->work_change);
2071: swap = pcbddc->work_change;
2072: pcbddc->work_change = r;
2073: r = swap;
2074: /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */
2075: if (pcbddc->benign_apply_coarse_only && pcbddc->exact_dirichlet_trick_app && pcbddc->change_interior) {
2076: VecCopy(r,pcis->vec1_global);
2077: VecLockReadPush(pcis->vec1_global);
2078: }
2079: }
2080: if (pcbddc->benign_have_null) { /* get p0 from r */
2081: PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);
2082: }
2083: if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
2084: VecCopy(r,z);
2085: /* First Dirichlet solve */
2086: VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
2087: VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
2088: /*
2089: Assembling right hand side for BDDC operator
2090: - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE)
2091: - pcis->vec1_B the interface part of the global vector z
2092: */
2093: if (n_D) {
2094: KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);
2095: KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec2_D);
2096: VecScale(pcis->vec2_D,m_one);
2097: if (pcbddc->switch_static) {
2098: VecSet(pcis->vec1_N,0.);
2099: VecScatterBegin(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2100: VecScatterEnd(pcis->N_to_D,pcis->vec2_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2101: if (!pcbddc->switch_static_change) {
2102: MatMultTranspose(lA,pcis->vec1_N,pcis->vec2_N);
2103: } else {
2104: MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
2105: MatMultTranspose(lA,pcis->vec2_N,pcis->vec1_N);
2106: MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
2107: }
2108: VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
2109: VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec1_D,ADD_VALUES,SCATTER_FORWARD);
2110: VecScatterBegin(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
2111: VecScatterEnd(pcis->N_to_B,pcis->vec2_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
2112: } else {
2113: MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);
2114: }
2115: } else {
2116: VecSet(pcis->vec1_B,zero);
2117: }
2118: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
2119: VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
2120: PCBDDCScalingRestriction(pc,z,pcis->vec1_B);
2121: } else {
2122: PCBDDCScalingRestriction(pc,r,pcis->vec1_B);
2123: }
2125: /* Apply interface preconditioner
2126: input/output vecs: pcis->vec1_B and pcis->vec1_D */
2127: PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);
2129: /* Apply transpose of partition of unity operator */
2130: PCBDDCScalingExtension(pc,pcis->vec1_B,z);
2132: /* Second Dirichlet solve and assembling of output */
2133: VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
2134: VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
2135: if (n_B) {
2136: if (pcbddc->switch_static) {
2137: VecScatterBegin(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2138: VecScatterEnd(pcis->N_to_D,pcis->vec1_D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2139: VecScatterBegin(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2140: VecScatterEnd(pcis->N_to_B,pcis->vec1_B,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
2141: if (!pcbddc->switch_static_change) {
2142: MatMultTranspose(lA,pcis->vec1_N,pcis->vec2_N);
2143: } else {
2144: MatMult(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
2145: MatMultTranspose(lA,pcis->vec2_N,pcis->vec1_N);
2146: MatMultTranspose(pcbddc->switch_static_change,pcis->vec1_N,pcis->vec2_N);
2147: }
2148: VecScatterBegin(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
2149: VecScatterEnd(pcis->N_to_D,pcis->vec2_N,pcis->vec3_D,INSERT_VALUES,SCATTER_FORWARD);
2150: } else {
2151: MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);
2152: }
2153: } else if (pcbddc->switch_static) { /* n_B is zero */
2154: if (!pcbddc->switch_static_change) {
2155: MatMultTranspose(lA,pcis->vec1_D,pcis->vec3_D);
2156: } else {
2157: MatMult(pcbddc->switch_static_change,pcis->vec1_D,pcis->vec1_N);
2158: MatMultTranspose(lA,pcis->vec1_N,pcis->vec2_N);
2159: MatMultTranspose(pcbddc->switch_static_change,pcis->vec2_N,pcis->vec3_D);
2160: }
2161: }
2162: KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);
2163: KSPCheckSolve(pcbddc->ksp_D,pc,pcis->vec4_D);
2164: if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) {
2165: if (pcbddc->switch_static) {
2166: VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);
2167: } else {
2168: VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);
2169: }
2170: VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
2171: VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
2172: } else {
2173: if (pcbddc->switch_static) {
2174: VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);
2175: } else {
2176: VecScale(pcis->vec4_D,m_one);
2177: }
2178: VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
2179: VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);
2180: }
2181: if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */
2182: PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);
2183: }
2184: if (pcbddc->ChangeOfBasisMatrix) {
2185: pcbddc->work_change = r;
2186: VecCopy(z,pcbddc->work_change);
2187: MatMult(pcbddc->ChangeOfBasisMatrix,pcbddc->work_change,z);
2188: }
2189: return(0);
2190: }
2192: PetscErrorCode PCReset_BDDC(PC pc)
2193: {
2194: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
2195: PC_IS *pcis = (PC_IS*)pc->data;
2196: KSP kspD,kspR,kspC;
2200: /* free BDDC custom data */
2201: PCBDDCResetCustomization(pc);
2202: /* destroy objects related to topography */
2203: PCBDDCResetTopography(pc);
2204: /* destroy objects for scaling operator */
2205: PCBDDCScalingDestroy(pc);
2206: /* free solvers stuff */
2207: PCBDDCResetSolvers(pc);
2208: /* free global vectors needed in presolve */
2209: VecDestroy(&pcbddc->temp_solution);
2210: VecDestroy(&pcbddc->original_rhs);
2211: /* free data created by PCIS */
2212: PCISDestroy(pc);
2214: /* restore defaults */
2215: kspD = pcbddc->ksp_D;
2216: kspR = pcbddc->ksp_R;
2217: kspC = pcbddc->coarse_ksp;
2218: PetscMemzero(pc->data,sizeof(*pcbddc));
2219: pcis->n_neigh = -1;
2220: pcis->scaling_factor = 1.0;
2221: pcis->reusesubmatrices = PETSC_TRUE;
2222: pcbddc->use_local_adj = PETSC_TRUE;
2223: pcbddc->use_vertices = PETSC_TRUE;
2224: pcbddc->use_edges = PETSC_TRUE;
2225: pcbddc->symmetric_primal = PETSC_TRUE;
2226: pcbddc->vertex_size = 1;
2227: pcbddc->recompute_topography = PETSC_TRUE;
2228: pcbddc->coarse_size = -1;
2229: pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
2230: pcbddc->coarsening_ratio = 8;
2231: pcbddc->coarse_eqs_per_proc = 1;
2232: pcbddc->benign_compute_correction = PETSC_TRUE;
2233: pcbddc->nedfield = -1;
2234: pcbddc->nedglobal = PETSC_TRUE;
2235: pcbddc->graphmaxcount = PETSC_MAX_INT;
2236: pcbddc->sub_schurs_layers = -1;
2237: pcbddc->ksp_D = kspD;
2238: pcbddc->ksp_R = kspR;
2239: pcbddc->coarse_ksp = kspC;
2240: return(0);
2241: }
2243: PetscErrorCode PCDestroy_BDDC(PC pc)
2244: {
2245: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
2249: PCReset_BDDC(pc);
2250: KSPDestroy(&pcbddc->ksp_D);
2251: KSPDestroy(&pcbddc->ksp_R);
2252: KSPDestroy(&pcbddc->coarse_ksp);
2253: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",NULL);
2254: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",NULL);
2255: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);
2256: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);
2257: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",NULL);
2258: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);
2259: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);
2260: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);
2261: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);
2262: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);
2263: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);
2264: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);
2265: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);
2266: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);
2267: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);
2268: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);
2269: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);
2270: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);
2271: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);
2272: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);
2273: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);
2274: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);
2275: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);
2276: PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL);
2277: PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);
2278: PetscFree(pc->data);
2279: return(0);
2280: }
2282: static PetscErrorCode PCSetCoordinates_BDDC(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
2283: {
2284: PC_BDDC *pcbddc = (PC_BDDC*)pc->data;
2285: PCBDDCGraph mat_graph = pcbddc->mat_graph;
2289: PetscFree(mat_graph->coords);
2290: PetscMalloc1(nloc*dim,&mat_graph->coords);
2291: PetscMemcpy(mat_graph->coords,coords,nloc*dim*sizeof(PetscReal));
2292: mat_graph->cnloc = nloc;
2293: mat_graph->cdim = dim;
2294: mat_graph->cloc = PETSC_FALSE;
2295: /* flg setup */
2296: pcbddc->recompute_topography = PETSC_TRUE;
2297: pcbddc->corner_selected = PETSC_FALSE;
2298: return(0);
2299: }
2301: static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool* change)
2302: {
2304: *change = PETSC_TRUE;
2305: return(0);
2306: }
2308: static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
2309: {
2310: FETIDPMat_ctx mat_ctx;
2311: Vec work;
2312: PC_IS* pcis;
2313: PC_BDDC* pcbddc;
2317: MatShellGetContext(fetidp_mat,&mat_ctx);
2318: pcis = (PC_IS*)mat_ctx->pc->data;
2319: pcbddc = (PC_BDDC*)mat_ctx->pc->data;
2321: VecSet(fetidp_flux_rhs,0.0);
2322: /* copy rhs since we may change it during PCPreSolve_BDDC */
2323: if (!pcbddc->original_rhs) {
2324: VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);
2325: }
2326: if (mat_ctx->rhs_flip) {
2327: VecPointwiseMult(pcbddc->original_rhs,standard_rhs,mat_ctx->rhs_flip);
2328: } else {
2329: VecCopy(standard_rhs,pcbddc->original_rhs);
2330: }
2331: if (mat_ctx->g2g_p) {
2332: /* interface pressure rhs */
2333: VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE);
2334: VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_rhs,pcbddc->original_rhs,INSERT_VALUES,SCATTER_REVERSE);
2335: VecScatterBegin(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD);
2336: VecScatterEnd(mat_ctx->g2g_p,standard_rhs,fetidp_flux_rhs,INSERT_VALUES,SCATTER_FORWARD);
2337: if (!mat_ctx->rhs_flip) {
2338: VecScale(fetidp_flux_rhs,-1.);
2339: }
2340: }
2341: /*
2342: change of basis for physical rhs if needed
2343: It also changes the rhs in case of dirichlet boundaries
2344: */
2345: PCPreSolve_BDDC(mat_ctx->pc,NULL,pcbddc->original_rhs,NULL);
2346: if (pcbddc->ChangeOfBasisMatrix) {
2347: MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcbddc->original_rhs,pcbddc->work_change);
2348: work = pcbddc->work_change;
2349: } else {
2350: work = pcbddc->original_rhs;
2351: }
2352: /* store vectors for computation of fetidp final solution */
2353: VecScatterBegin(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
2354: VecScatterEnd(pcis->global_to_D,work,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);
2355: /* scale rhs since it should be unassembled */
2356: /* TODO use counter scaling? (also below) */
2357: VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2358: VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2359: /* Apply partition of unity */
2360: VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
2361: /* PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B); */
2362: if (!pcbddc->switch_static) {
2363: /* compute partially subassembled Schur complement right-hand side */
2364: KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);
2365: /* Cannot propagate up error in KSPSolve() because there is no access to the PC */
2366: MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);
2367: VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);
2368: VecSet(work,0.0);
2369: VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE);
2370: VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,work,ADD_VALUES,SCATTER_REVERSE);
2371: /* PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B); */
2372: VecScatterBegin(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2373: VecScatterEnd(pcis->global_to_B,work,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);
2374: VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);
2375: }
2376: /* BDDC rhs */
2377: VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);
2378: if (pcbddc->switch_static) {
2379: VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
2380: }
2381: /* apply BDDC */
2382: PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
2383: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
2384: PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
2386: /* Application of B_delta and assembling of rhs for fetidp fluxes */
2387: MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
2388: VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2389: VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2390: /* Add contribution to interface pressures */
2391: if (mat_ctx->l2g_p) {
2392: MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);
2393: if (pcbddc->switch_static) {
2394: MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);
2395: }
2396: VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2397: VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);
2398: }
2399: return(0);
2400: }
2402: /*@
2403: PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side
2405: Collective
2407: Input Parameters:
2408: + fetidp_mat - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators
2409: - standard_rhs - the right-hand side of the original linear system
2411: Output Parameters:
2412: . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system
2414: Level: developer
2416: Notes:
2418: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution
2419: @*/
2420: PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs)
2421: {
2422: FETIDPMat_ctx mat_ctx;
2429: MatShellGetContext(fetidp_mat,&mat_ctx);
2430: PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));
2431: return(0);
2432: }
2434: static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2435: {
2436: FETIDPMat_ctx mat_ctx;
2437: PC_IS* pcis;
2438: PC_BDDC* pcbddc;
2440: Vec work;
2443: MatShellGetContext(fetidp_mat,&mat_ctx);
2444: pcis = (PC_IS*)mat_ctx->pc->data;
2445: pcbddc = (PC_BDDC*)mat_ctx->pc->data;
2447: /* apply B_delta^T */
2448: VecSet(pcis->vec1_B,0.);
2449: VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
2450: VecScatterEnd(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
2451: MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
2452: if (mat_ctx->l2g_p) {
2453: VecScatterBegin(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
2454: VecScatterEnd(mat_ctx->l2g_p,fetidp_flux_sol,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
2455: MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);
2456: }
2458: /* compute rhs for BDDC application */
2459: VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);
2460: if (pcbddc->switch_static) {
2461: VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);
2462: if (mat_ctx->l2g_p) {
2463: VecScale(mat_ctx->vP,-1.);
2464: MatMultAdd(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D,pcis->vec1_D);
2465: }
2466: }
2468: /* apply BDDC */
2469: PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
2470: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
2472: /* put values into global vector */
2473: if (pcbddc->ChangeOfBasisMatrix) work = pcbddc->work_change;
2474: else work = standard_sol;
2475: VecScatterBegin(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE);
2476: VecScatterEnd(pcis->global_to_B,pcis->vec1_B,work,INSERT_VALUES,SCATTER_REVERSE);
2477: if (!pcbddc->switch_static) {
2478: /* compute values into the interior if solved for the partially subassembled Schur complement */
2479: MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);
2480: VecAYPX(pcis->vec1_D,-1.0,mat_ctx->temp_solution_D);
2481: KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec1_D);
2482: /* Cannot propagate up error in KSPSolve() because there is no access to the PC */
2483: }
2485: VecScatterBegin(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE);
2486: VecScatterEnd(pcis->global_to_D,pcis->vec1_D,work,INSERT_VALUES,SCATTER_REVERSE);
2487: /* add p0 solution to final solution */
2488: PCBDDCBenignGetOrSetP0(mat_ctx->pc,work,PETSC_FALSE);
2489: if (pcbddc->ChangeOfBasisMatrix) {
2490: MatMult(pcbddc->ChangeOfBasisMatrix,work,standard_sol);
2491: }
2492: PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);
2493: if (mat_ctx->g2g_p) {
2494: VecScatterBegin(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
2495: VecScatterEnd(mat_ctx->g2g_p,fetidp_flux_sol,standard_sol,INSERT_VALUES,SCATTER_REVERSE);
2496: }
2497: return(0);
2498: }
2500: static PetscErrorCode PCView_BDDCIPC(PC pc, PetscViewer viewer)
2501: {
2503: BDDCIPC_ctx bddcipc_ctx;
2504: PetscBool isascii;
2507: PCShellGetContext(pc,(void **)&bddcipc_ctx);
2508: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
2509: if (isascii) {
2510: PetscViewerASCIIPrintf(viewer,"BDDC interface preconditioner\n");
2511: }
2512: PetscViewerASCIIPushTab(viewer);
2513: PCView(bddcipc_ctx->bddc,viewer);
2514: PetscViewerASCIIPopTab(viewer);
2515: return(0);
2516: }
2518: static PetscErrorCode PCSetUp_BDDCIPC(PC pc)
2519: {
2521: BDDCIPC_ctx bddcipc_ctx;
2522: PetscBool isbddc;
2523: Vec vv;
2524: IS is;
2525: PC_IS *pcis;
2528: PCShellGetContext(pc,(void **)&bddcipc_ctx);
2529: PetscObjectTypeCompare((PetscObject)bddcipc_ctx->bddc,PCBDDC,&isbddc);
2530: if (!isbddc) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid type %s. Must be of type bddc",((PetscObject)bddcipc_ctx->bddc)->type_name);
2531: PCSetUp(bddcipc_ctx->bddc);
2533: /* create interface scatter */
2534: pcis = (PC_IS*)(bddcipc_ctx->bddc->data);
2535: VecScatterDestroy(&bddcipc_ctx->g2l);
2536: MatCreateVecs(pc->pmat,&vv,NULL);
2537: ISRenumber(pcis->is_B_global,NULL,NULL,&is);
2538: VecScatterCreate(vv,is,pcis->vec1_B,NULL,&bddcipc_ctx->g2l);
2539: ISDestroy(&is);
2540: VecDestroy(&vv);
2541: return(0);
2542: }
2544: static PetscErrorCode PCApply_BDDCIPC(PC pc, Vec r, Vec x)
2545: {
2547: BDDCIPC_ctx bddcipc_ctx;
2548: PC_IS *pcis;
2549: VecScatter tmps;
2552: PCShellGetContext(pc,(void **)&bddcipc_ctx);
2553: pcis = (PC_IS*)(bddcipc_ctx->bddc->data);
2554: tmps = pcis->global_to_B;
2555: pcis->global_to_B = bddcipc_ctx->g2l;
2556: PCBDDCScalingRestriction(bddcipc_ctx->bddc,r,pcis->vec1_B);
2557: PCBDDCApplyInterfacePreconditioner(bddcipc_ctx->bddc,PETSC_FALSE);
2558: PCBDDCScalingExtension(bddcipc_ctx->bddc,pcis->vec1_B,x);
2559: pcis->global_to_B = tmps;
2560: return(0);
2561: }
2563: static PetscErrorCode PCApplyTranspose_BDDCIPC(PC pc, Vec r, Vec x)
2564: {
2566: BDDCIPC_ctx bddcipc_ctx;
2567: PC_IS *pcis;
2568: VecScatter tmps;
2571: PCShellGetContext(pc,(void **)&bddcipc_ctx);
2572: pcis = (PC_IS*)(bddcipc_ctx->bddc->data);
2573: tmps = pcis->global_to_B;
2574: pcis->global_to_B = bddcipc_ctx->g2l;
2575: PCBDDCScalingRestriction(bddcipc_ctx->bddc,r,pcis->vec1_B);
2576: PCBDDCApplyInterfacePreconditioner(bddcipc_ctx->bddc,PETSC_TRUE);
2577: PCBDDCScalingExtension(bddcipc_ctx->bddc,pcis->vec1_B,x);
2578: pcis->global_to_B = tmps;
2579: return(0);
2580: }
2582: static PetscErrorCode PCDestroy_BDDCIPC(PC pc)
2583: {
2585: BDDCIPC_ctx bddcipc_ctx;
2588: PCShellGetContext(pc,(void **)&bddcipc_ctx);
2589: PCDestroy(&bddcipc_ctx->bddc);
2590: VecScatterDestroy(&bddcipc_ctx->g2l);
2591: PetscFree(bddcipc_ctx);
2592: return(0);
2593: }
2595: /*@
2596: PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system
2598: Collective
2600: Input Parameters:
2601: + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators
2602: - fetidp_flux_sol - the solution of the FETI-DP linear system
2604: Output Parameters:
2605: . standard_sol - the solution defined on the physical domain
2607: Level: developer
2609: Notes:
2611: .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS
2612: @*/
2613: PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol)
2614: {
2615: FETIDPMat_ctx mat_ctx;
2622: MatShellGetContext(fetidp_mat,&mat_ctx);
2623: PetscUseMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));
2624: return(0);
2625: }
2627: static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, const char* prefix, Mat *fetidp_mat, PC *fetidp_pc)
2628: {
2630: FETIDPMat_ctx fetidpmat_ctx;
2631: Mat newmat;
2632: FETIDPPC_ctx fetidppc_ctx;
2633: PC newpc;
2634: MPI_Comm comm;
2638: PetscObjectGetComm((PetscObject)pc,&comm);
2639: /* FETI-DP matrix */
2640: PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);
2641: fetidpmat_ctx->fully_redundant = fully_redundant;
2642: PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);
2643: MatCreateShell(comm,fetidpmat_ctx->n,fetidpmat_ctx->n,fetidpmat_ctx->N,fetidpmat_ctx->N,fetidpmat_ctx,&newmat);
2644: PetscObjectSetName((PetscObject)newmat,!fetidpmat_ctx->l2g_lambda_only ? "F" : "G");
2645: MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);
2646: MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);
2647: MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);
2648: /* propagate MatOptions */
2649: {
2650: PC_BDDC *pcbddc = (PC_BDDC*)fetidpmat_ctx->pc->data;
2651: PetscBool issym;
2653: MatGetOption(pc->mat,MAT_SYMMETRIC,&issym);
2654: if (issym || pcbddc->symmetric_primal) {
2655: MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);
2656: }
2657: }
2658: MatSetOptionsPrefix(newmat,prefix);
2659: MatAppendOptionsPrefix(newmat,"fetidp_");
2660: MatSetUp(newmat);
2661: /* FETI-DP preconditioner */
2662: PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);
2663: PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);
2664: PCCreate(comm,&newpc);
2665: PCSetOperators(newpc,newmat,newmat);
2666: PCSetOptionsPrefix(newpc,prefix);
2667: PCAppendOptionsPrefix(newpc,"fetidp_");
2668: PCSetErrorIfFailure(newpc,pc->erroriffailure);
2669: if (!fetidpmat_ctx->l2g_lambda_only) { /* standard FETI-DP */
2670: PCSetType(newpc,PCSHELL);
2671: PCShellSetName(newpc,"FETI-DP multipliers");
2672: PCShellSetContext(newpc,fetidppc_ctx);
2673: PCShellSetApply(newpc,FETIDPPCApply);
2674: PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);
2675: PCShellSetView(newpc,FETIDPPCView);
2676: PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);
2677: } else { /* saddle-point FETI-DP */
2678: Mat M;
2679: PetscInt psize;
2680: PetscBool fake = PETSC_FALSE, isfieldsplit;
2682: ISViewFromOptions(fetidpmat_ctx->lagrange,NULL,"-lag_view");
2683: ISViewFromOptions(fetidpmat_ctx->pressure,NULL,"-press_view");
2684: PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_PPmat",(PetscObject*)&M);
2685: PCSetType(newpc,PCFIELDSPLIT);
2686: PCFieldSplitSetIS(newpc,"lag",fetidpmat_ctx->lagrange);
2687: PCFieldSplitSetIS(newpc,"p",fetidpmat_ctx->pressure);
2688: PCFieldSplitSetType(newpc,PC_COMPOSITE_SCHUR);
2689: PCFieldSplitSetSchurFactType(newpc,PC_FIELDSPLIT_SCHUR_FACT_DIAG);
2690: ISGetSize(fetidpmat_ctx->pressure,&psize);
2691: if (psize != M->rmap->N) {
2692: Mat M2;
2693: PetscInt lpsize;
2695: fake = PETSC_TRUE;
2696: ISGetLocalSize(fetidpmat_ctx->pressure,&lpsize);
2697: MatCreate(comm,&M2);
2698: MatSetType(M2,MATAIJ);
2699: MatSetSizes(M2,lpsize,lpsize,psize,psize);
2700: MatSetUp(M2);
2701: MatAssemblyBegin(M2,MAT_FINAL_ASSEMBLY);
2702: MatAssemblyEnd(M2,MAT_FINAL_ASSEMBLY);
2703: PCFieldSplitSetSchurPre(newpc,PC_FIELDSPLIT_SCHUR_PRE_USER,M2);
2704: MatDestroy(&M2);
2705: } else {
2706: PCFieldSplitSetSchurPre(newpc,PC_FIELDSPLIT_SCHUR_PRE_USER,M);
2707: }
2708: PCFieldSplitSetSchurScale(newpc,1.0);
2710: /* we need to setfromoptions and setup here to access the blocks */
2711: PCSetFromOptions(newpc);
2712: PCSetUp(newpc);
2714: /* user may have changed the type (e.g. -fetidp_pc_type none) */
2715: PetscObjectTypeCompare((PetscObject)newpc,PCFIELDSPLIT,&isfieldsplit);
2716: if (isfieldsplit) {
2717: KSP *ksps;
2718: PC ppc,lagpc;
2719: PetscInt nn;
2720: PetscBool ismatis,matisok = PETSC_FALSE,check = PETSC_FALSE;
2722: /* set the solver for the (0,0) block */
2723: PCFieldSplitSchurGetSubKSP(newpc,&nn,&ksps);
2724: if (!nn) { /* not of type PC_COMPOSITE_SCHUR */
2725: PCFieldSplitGetSubKSP(newpc,&nn,&ksps);
2726: if (!fake) { /* pass pmat to the pressure solver */
2727: Mat F;
2729: KSPGetOperators(ksps[1],&F,NULL);
2730: KSPSetOperators(ksps[1],F,M);
2731: }
2732: } else {
2733: PetscBool issym;
2734: Mat S;
2736: PCFieldSplitSchurGetS(newpc,&S);
2738: MatGetOption(newmat,MAT_SYMMETRIC,&issym);
2739: if (issym) {
2740: MatSetOption(S,MAT_SYMMETRIC,PETSC_TRUE);
2741: }
2742: }
2743: KSPGetPC(ksps[0],&lagpc);
2744: PCSetType(lagpc,PCSHELL);
2745: PCShellSetName(lagpc,"FETI-DP multipliers");
2746: PCShellSetContext(lagpc,fetidppc_ctx);
2747: PCShellSetApply(lagpc,FETIDPPCApply);
2748: PCShellSetApplyTranspose(lagpc,FETIDPPCApplyTranspose);
2749: PCShellSetView(lagpc,FETIDPPCView);
2750: PCShellSetDestroy(lagpc,PCBDDCDestroyFETIDPPC);
2752: /* Olof's idea: interface Schur complement preconditioner for the mass matrix */
2753: KSPGetPC(ksps[1],&ppc);
2754: if (fake) {
2755: BDDCIPC_ctx bddcipc_ctx;
2756: PetscContainer c;
2758: matisok = PETSC_TRUE;
2760: /* create inner BDDC solver */
2761: PetscNew(&bddcipc_ctx);
2762: PCCreate(comm,&bddcipc_ctx->bddc);
2763: PCSetType(bddcipc_ctx->bddc,PCBDDC);
2764: PCSetOperators(bddcipc_ctx->bddc,M,M);
2765: PetscObjectQuery((PetscObject)pc,"__KSPFETIDP_pCSR",(PetscObject*)&c);
2766: PetscObjectTypeCompare((PetscObject)M,MATIS,&ismatis);
2767: if (c && ismatis) {
2768: Mat lM;
2769: PetscInt *csr,n;
2771: MatISGetLocalMat(M,&lM);
2772: MatGetSize(lM,&n,NULL);
2773: PetscContainerGetPointer(c,(void**)&csr);
2774: PCBDDCSetLocalAdjacencyGraph(bddcipc_ctx->bddc,n,csr,csr + (n + 1),PETSC_COPY_VALUES);
2775: MatISRestoreLocalMat(M,&lM);
2776: }
2777: PCSetOptionsPrefix(bddcipc_ctx->bddc,((PetscObject)ksps[1])->prefix);
2778: PCSetErrorIfFailure(bddcipc_ctx->bddc,pc->erroriffailure);
2779: PCSetFromOptions(bddcipc_ctx->bddc);
2781: /* wrap the interface application */
2782: PCSetType(ppc,PCSHELL);
2783: PCShellSetName(ppc,"FETI-DP pressure");
2784: PCShellSetContext(ppc,bddcipc_ctx);
2785: PCShellSetSetUp(ppc,PCSetUp_BDDCIPC);
2786: PCShellSetApply(ppc,PCApply_BDDCIPC);
2787: PCShellSetApplyTranspose(ppc,PCApplyTranspose_BDDCIPC);
2788: PCShellSetView(ppc,PCView_BDDCIPC);
2789: PCShellSetDestroy(ppc,PCDestroy_BDDCIPC);
2790: }
2792: /* determine if we need to assemble M to construct a preconditioner */
2793: if (!matisok) {
2794: PetscObjectTypeCompare((PetscObject)M,MATIS,&ismatis);
2795: PetscObjectTypeCompareAny((PetscObject)ppc,&matisok,PCBDDC,PCJACOBI,PCNONE,PCMG,"");
2796: if (ismatis && !matisok) {
2797: MatConvert(M,MATAIJ,MAT_INPLACE_MATRIX,&M);
2798: }
2799: }
2801: /* run the subproblems to check convergence */
2802: PetscOptionsGetBool(NULL,((PetscObject)newmat)->prefix,"-check_saddlepoint",&check,NULL);
2803: if (check) {
2804: PetscInt i;
2806: for (i=0;i<nn;i++) {
2807: KSP kspC;
2808: PC pc;
2809: Mat F,pF;
2810: Vec x,y;
2811: PetscBool isschur,prec = PETSC_TRUE;
2813: KSPCreate(PetscObjectComm((PetscObject)ksps[i]),&kspC);
2814: KSPSetOptionsPrefix(kspC,((PetscObject)ksps[i])->prefix);
2815: KSPAppendOptionsPrefix(kspC,"check_");
2816: KSPGetOperators(ksps[i],&F,&pF);
2817: PetscObjectTypeCompare((PetscObject)F,MATSCHURCOMPLEMENT,&isschur);
2818: if (isschur) {
2819: KSP kspS,kspS2;
2820: Mat A00,pA00,A10,A01,A11;
2821: char prefix[256];
2823: MatSchurComplementGetKSP(F,&kspS);
2824: MatSchurComplementGetSubMatrices(F,&A00,&pA00,&A01,&A10,&A11);
2825: MatCreateSchurComplement(A00,pA00,A01,A10,A11,&F);
2826: MatSchurComplementGetKSP(F,&kspS2);
2827: PetscSNPrintf(prefix,sizeof(prefix),"%sschur_",((PetscObject)kspC)->prefix);
2828: KSPSetOptionsPrefix(kspS2,prefix);
2829: KSPGetPC(kspS2,&pc);
2830: PCSetType(pc,PCKSP);
2831: PCKSPSetKSP(pc,kspS);
2832: KSPSetFromOptions(kspS2);
2833: KSPGetPC(kspS2,&pc);
2834: PCSetUseAmat(pc,PETSC_TRUE);
2835: } else {
2836: PetscObjectReference((PetscObject)F);
2837: }
2838: KSPSetFromOptions(kspC);
2839: PetscOptionsGetBool(NULL,((PetscObject)kspC)->prefix,"-preconditioned",&prec,NULL);
2840: if (prec) {
2841: KSPGetPC(ksps[i],&pc);
2842: KSPSetPC(kspC,pc);
2843: }
2844: KSPSetOperators(kspC,F,pF);
2845: MatCreateVecs(F,&x,&y);
2846: VecSetRandom(x,NULL);
2847: MatMult(F,x,y);
2848: KSPSolve(kspC,y,x);
2849: KSPCheckSolve(kspC,pc,x);
2850: KSPDestroy(&kspC);
2851: MatDestroy(&F);
2852: VecDestroy(&x);
2853: VecDestroy(&y);
2854: }
2855: }
2856: PetscFree(ksps);
2857: }
2858: }
2859: /* return pointers for objects created */
2860: *fetidp_mat = newmat;
2861: *fetidp_pc = newpc;
2862: return(0);
2863: }
2865: /*@C
2866: PCBDDCCreateFETIDPOperators - Create FETI-DP operators
2868: Collective
2870: Input Parameters:
2871: + pc - the BDDC preconditioning context (setup should have been called before)
2872: . fully_redundant - true for a fully redundant set of Lagrange multipliers
2873: - prefix - optional options database prefix for the objects to be created (can be NULL)
2875: Output Parameters:
2876: + fetidp_mat - shell FETI-DP matrix object
2877: - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix
2879: Level: developer
2881: Notes:
2882: Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose
2884: .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution
2885: @*/
2886: PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc)
2887: {
2892: if (pc->setupcalled) {
2893: PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,PetscBool,const char*,Mat*,PC*),(pc,fully_redundant,prefix,fetidp_mat,fetidp_pc));
2894: } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"You must call PCSetup_BDDC() first");
2895: return(0);
2896: }
2897: /* -------------------------------------------------------------------------- */
2898: /*MC
2899: PCBDDC - Balancing Domain Decomposition by Constraints.
2901: An implementation of the BDDC preconditioner based on
2903: .vb
2904: [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007
2905: [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
2906: [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977
2907: [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf
2908: .ve
2910: The matrix to be preconditioned (Pmat) must be of type MATIS.
2912: Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers.
2914: It also works with unsymmetric and indefinite problems.
2916: 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.
2918: Approximate local solvers are automatically adapted (see [1]) if the user has attached a nullspace object to the subdomain matrices, and informed BDDC of using approximate solvers (via the command line).
2920: Boundary nodes are split in vertices, edges and faces classes using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph()
2921: Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesIS() and their local counterparts.
2923: Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD.
2925: Change of basis is performed similarly to [2] when requested. When more than one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used.
2926: User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat()
2928: The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object.
2930: Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS or MKL_PARDISO are present. Future versions of the code will also consider using PASTIX.
2932: An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases.
2933: Deluxe scaling is not supported yet for FETI-DP.
2935: Options Database Keys (some of them, run with -h for a complete list):
2937: . -pc_bddc_use_vertices <true> - use or not vertices in primal space
2938: . -pc_bddc_use_edges <true> - use or not edges in primal space
2939: . -pc_bddc_use_faces <false> - use or not faces in primal space
2940: . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems
2941: . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only)
2942: . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested
2943: . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1])
2944: . -pc_bddc_levels <0> - maximum number of levels for multilevel
2945: . -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case)
2946: . -pc_bddc_coarse_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level)
2947: . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling
2948: . -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling)
2949: . -pc_bddc_adaptive_threshold <0.0> - when a value different than zero is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS or MKL_PARDISO installed)
2950: - -pc_bddc_check_level <0> - set verbosity level of debugging output
2952: Options for Dirichlet, Neumann or coarse solver can be set with
2953: .vb
2954: -pc_bddc_dirichlet_
2955: -pc_bddc_neumann_
2956: -pc_bddc_coarse_
2957: .ve
2958: e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU.
2960: When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as
2961: .vb
2962: -pc_bddc_dirichlet_lN_
2963: -pc_bddc_neumann_lN_
2964: -pc_bddc_coarse_lN_
2965: .ve
2966: Note that level number ranges from the finest (0) to the coarsest (N).
2967: In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g.
2968: .vb
2969: -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3
2970: .ve
2971: will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors
2973: Level: intermediate
2975: Developer Notes:
2977: Contributed by Stefano Zampini
2979: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS
2980: M*/
2982: PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc)
2983: {
2984: PetscErrorCode ierr;
2985: PC_BDDC *pcbddc;
2988: PetscNewLog(pc,&pcbddc);
2989: pc->data = (void*)pcbddc;
2991: /* create PCIS data structure */
2992: PCISCreate(pc);
2994: /* create local graph structure */
2995: PCBDDCGraphCreate(&pcbddc->mat_graph);
2997: /* BDDC nonzero defaults */
2998: pcbddc->use_local_adj = PETSC_TRUE;
2999: pcbddc->use_vertices = PETSC_TRUE;
3000: pcbddc->use_edges = PETSC_TRUE;
3001: pcbddc->symmetric_primal = PETSC_TRUE;
3002: pcbddc->vertex_size = 1;
3003: pcbddc->recompute_topography = PETSC_TRUE;
3004: pcbddc->coarse_size = -1;
3005: pcbddc->use_exact_dirichlet_trick = PETSC_TRUE;
3006: pcbddc->coarsening_ratio = 8;
3007: pcbddc->coarse_eqs_per_proc = 1;
3008: pcbddc->benign_compute_correction = PETSC_TRUE;
3009: pcbddc->nedfield = -1;
3010: pcbddc->nedglobal = PETSC_TRUE;
3011: pcbddc->graphmaxcount = PETSC_MAX_INT;
3012: pcbddc->sub_schurs_layers = -1;
3013: pcbddc->adaptive_threshold[0] = 0.0;
3014: pcbddc->adaptive_threshold[1] = 0.0;
3016: /* function pointers */
3017: pc->ops->apply = PCApply_BDDC;
3018: pc->ops->applytranspose = PCApplyTranspose_BDDC;
3019: pc->ops->setup = PCSetUp_BDDC;
3020: pc->ops->destroy = PCDestroy_BDDC;
3021: pc->ops->setfromoptions = PCSetFromOptions_BDDC;
3022: pc->ops->view = PCView_BDDC;
3023: pc->ops->applyrichardson = 0;
3024: pc->ops->applysymmetricleft = 0;
3025: pc->ops->applysymmetricright = 0;
3026: pc->ops->presolve = PCPreSolve_BDDC;
3027: pc->ops->postsolve = PCPostSolve_BDDC;
3028: pc->ops->reset = PCReset_BDDC;
3030: /* composing function */
3031: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDiscreteGradient_C",PCBDDCSetDiscreteGradient_BDDC);
3032: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDivergenceMat_C",PCBDDCSetDivergenceMat_BDDC);
3033: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);
3034: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);
3035: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesIS_C",PCBDDCSetPrimalVerticesIS_BDDC);
3036: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetPrimalVerticesLocalIS_C",PCBDDCGetPrimalVerticesLocalIS_BDDC);
3037: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetPrimalVerticesIS_C",PCBDDCGetPrimalVerticesIS_BDDC);
3038: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);
3039: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);
3040: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);
3041: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);
3042: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);
3043: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);
3044: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);
3045: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);
3046: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);
3047: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);
3048: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);
3049: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);
3050: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);
3051: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);
3052: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);
3053: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);
3054: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);
3055: PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);
3056: PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_BDDC);
3057: PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_BDDC);
3058: return(0);
3059: }
3061: /*@C
3062: PCBDDCInitializePackage - This function initializes everything in the PCBDDC package. It is called
3063: from PCInitializePackage().
3065: Level: developer
3067: .keywords: PC, PCBDDC, initialize, package
3068: .seealso: PetscInitialize()
3069: @*/
3070: PetscErrorCode PCBDDCInitializePackage(void)
3071: {
3073: int i;
3076: if (PCBDDCPackageInitialized) return(0);
3077: PCBDDCPackageInitialized = PETSC_TRUE;
3078: PetscRegisterFinalize(PCBDDCFinalizePackage);
3080: /* general events */
3081: PetscLogEventRegister("PCBDDCTopo",PC_CLASSID,&PC_BDDC_Topology[0]);
3082: PetscLogEventRegister("PCBDDCLKSP",PC_CLASSID,&PC_BDDC_LocalSolvers[0]);
3083: PetscLogEventRegister("PCBDDCLWor",PC_CLASSID,&PC_BDDC_LocalWork[0]);
3084: PetscLogEventRegister("PCBDDCCorr",PC_CLASSID,&PC_BDDC_CorrectionSetUp[0]);
3085: PetscLogEventRegister("PCBDDCCSet",PC_CLASSID,&PC_BDDC_CoarseSetUp[0]);
3086: PetscLogEventRegister("PCBDDCCKSP",PC_CLASSID,&PC_BDDC_CoarseSolver[0]);
3087: PetscLogEventRegister("PCBDDCAdap",PC_CLASSID,&PC_BDDC_AdaptiveSetUp[0]);
3088: PetscLogEventRegister("PCBDDCScal",PC_CLASSID,&PC_BDDC_Scaling[0]);
3089: PetscLogEventRegister("PCBDDCSchr",PC_CLASSID,&PC_BDDC_Schurs[0]);
3090: for (i=1;i<PETSC_PCBDDC_MAXLEVELS;i++) {
3091: char ename[32];
3093: PetscSNPrintf(ename,sizeof(ename),"PCBDDCTopo l%02d",i);
3094: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Topology[i]);
3095: PetscSNPrintf(ename,sizeof(ename),"PCBDDCLKSP l%02d",i);
3096: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_LocalSolvers[i]);
3097: PetscSNPrintf(ename,sizeof(ename),"PCBDDCLWor l%02d",i);
3098: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_LocalWork[i]);
3099: PetscSNPrintf(ename,sizeof(ename),"PCBDDCCorr l%02d",i);
3100: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_CorrectionSetUp[i]);
3101: PetscSNPrintf(ename,sizeof(ename),"PCBDDCCSet l%02d",i);
3102: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_CoarseSetUp[i]);
3103: PetscSNPrintf(ename,sizeof(ename),"PCBDDCCKSP l%02d",i);
3104: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_CoarseSolver[i]);
3105: PetscSNPrintf(ename,sizeof(ename),"PCBDDCAdap l%02d",i);
3106: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_AdaptiveSetUp[i]);
3107: PetscSNPrintf(ename,sizeof(ename),"PCBDDCScal l%02d",i);
3108: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Scaling[i]);
3109: PetscSNPrintf(ename,sizeof(ename),"PCBDDCSchr l%02d",i);
3110: PetscLogEventRegister(ename,PC_CLASSID,&PC_BDDC_Schurs[i]);
3111: }
3112: return(0);
3113: }
3115: /*@C
3116: PCBDDCFinalizePackage - This function frees everything from the PCBDDC package. It is
3117: called from PetscFinalize() automatically.
3119: Level: developer
3121: .keywords: Petsc, destroy, package
3122: .seealso: PetscFinalize()
3123: @*/
3124: PetscErrorCode PCBDDCFinalizePackage(void)
3125: {
3127: PCBDDCPackageInitialized = PETSC_FALSE;
3128: return(0);
3129: }