xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision 026de3107f93ffb842e26d1e8c9e3f422c3d1ec0)
1 #include "bddc.h"
2 #include "bddcprivate.h"
3 #include <petscblaslapack.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "PCBDDCResetCustomization"
7 PetscErrorCode PCBDDCResetCustomization(PC pc)
8 {
9   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
10   PetscInt       i;
11   PetscErrorCode ierr;
12 
13   PetscFunctionBegin;
14   ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
15   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
16   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
17   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
18   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
19   for (i=0;i<pcbddc->n_ISForDofs;i++) {
20     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
21   }
22   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "PCBDDCResetTopography"
28 PetscErrorCode PCBDDCResetTopography(PC pc)
29 {
30   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
31   PetscErrorCode ierr;
32 
33   PetscFunctionBegin;
34   ierr = MatDestroy(&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
35   ierr = MatDestroy(&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
36   ierr = PCBDDCGraphReset(pcbddc->mat_graph);CHKERRQ(ierr);
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "PCBDDCResetSolvers"
42 PetscErrorCode PCBDDCResetSolvers(PC pc)
43 {
44   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
49   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
50   ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
51   ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr);
52   ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr);
53   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
54   ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr);
55   ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr);
56   ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr);
57   ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr);
58   ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr);
59   ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr);
60   ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr);
61   ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr);
62   ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr);
63   ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr);
64   ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr);
65   ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr);
66   ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr);
67   ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr);
68   ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr);
69   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
70   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
71   ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr);
72   ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr);
73   ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr);
74   ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr);
75   ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "PCBDDCSolveSaddlePoint"
81 static PetscErrorCode  PCBDDCSolveSaddlePoint(PC pc)
82 {
83   PetscErrorCode ierr;
84   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
85 
86   PetscFunctionBegin;
87   ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr);
88   if (pcbddc->local_auxmat1) {
89     ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr);
90     ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr);
91   }
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner"
97 PetscErrorCode  PCBDDCApplyInterfacePreconditioner(PC pc)
98 {
99   PetscErrorCode ierr;
100   PC_BDDC*        pcbddc = (PC_BDDC*)(pc->data);
101   PC_IS*            pcis = (PC_IS*)  (pc->data);
102   const PetscScalar zero = 0.0;
103 
104   PetscFunctionBegin;
105   /* Application of PHI^T (or PSI^T)  */
106   if (pcbddc->coarse_psi_B) {
107     ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
108     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
109   } else {
110     ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
111     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
112   }
113   /* Scatter data of coarse_rhs */
114   if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); }
115   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
116 
117   /* Local solution on R nodes */
118   ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr);
119   ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
120   ierr = VecScatterEnd  (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
121   if (pcbddc->inexact_prec_type) {
122     ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
123     ierr = VecScatterEnd  (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
124   }
125   ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr);
126   ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
127   ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
128   ierr = VecScatterEnd  (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
129   if (pcbddc->inexact_prec_type) {
130     ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
131     ierr = VecScatterEnd  (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132   }
133 
134   /* Coarse solution */
135   ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
136   if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */
137     ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr);
138   }
139   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
140   ierr = PCBDDCScatterCoarseDataEnd  (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
141 
142   /* Sum contributions from two levels */
143   ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
144   if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "PCBDDCScatterCoarseDataBegin"
150 PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
151 {
152   PetscErrorCode ierr;
153   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
154 
155   PetscFunctionBegin;
156   switch (pcbddc->coarse_communications_type) {
157     case SCATTERS_BDDC:
158       ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
159       break;
160     case GATHERS_BDDC:
161       break;
162   }
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "PCBDDCScatterCoarseDataEnd"
168 PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
169 {
170   PetscErrorCode ierr;
171   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
172   PetscScalar*   array_to;
173   PetscScalar*   array_from;
174   MPI_Comm       comm;
175   PetscInt       i;
176 
177   PetscFunctionBegin;
178   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
179   switch (pcbddc->coarse_communications_type) {
180     case SCATTERS_BDDC:
181       ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
182       break;
183     case GATHERS_BDDC:
184       if (vec_from) {
185         ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr);
186       }
187       if (vec_to) {
188         ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr);
189       }
190       switch(pcbddc->coarse_problem_type){
191         case SEQUENTIAL_BDDC:
192           if (smode == SCATTER_FORWARD) {
193             ierr = MPI_Gatherv(&array_from[0],pcbddc->local_primal_size,MPIU_SCALAR,&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,0,comm);CHKERRQ(ierr);
194             if (vec_to) {
195               if (imode == ADD_VALUES) {
196                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
197                   array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
198                 }
199               } else {
200                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
201                   array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
202                 }
203               }
204             }
205           } else {
206             if (vec_from) {
207               if (imode == ADD_VALUES) {
208                 MPI_Comm vec_from_comm;
209                 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr);
210                 SETERRQ2(vec_from_comm,PETSC_ERR_SUP,"Unsupported insert mode ADD_VALUES for SCATTER_REVERSE in %s for case %d\n",__FUNCT__,pcbddc->coarse_problem_type);
211               }
212               for (i=0;i<pcbddc->replicated_primal_size;i++) {
213                 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]];
214               }
215             }
216             ierr = MPI_Scatterv(&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,&array_to[0],pcbddc->local_primal_size,MPIU_SCALAR,0,comm);CHKERRQ(ierr);
217           }
218           break;
219         case REPLICATED_BDDC:
220           if (smode == SCATTER_FORWARD) {
221             ierr = MPI_Allgatherv(&array_from[0],pcbddc->local_primal_size,MPIU_SCALAR,&pcbddc->replicated_local_primal_values[0],pcbddc->local_primal_sizes,pcbddc->local_primal_displacements,MPIU_SCALAR,comm);CHKERRQ(ierr);
222             if (imode == ADD_VALUES) {
223               for (i=0;i<pcbddc->replicated_primal_size;i++) {
224                 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
225               }
226             } else {
227               for (i=0;i<pcbddc->replicated_primal_size;i++) {
228                 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
229               }
230             }
231           } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */
232             if (imode == ADD_VALUES) {
233               for (i=0;i<pcbddc->local_primal_size;i++) {
234                 array_to[i]+=array_from[pcbddc->local_primal_indices[i]];
235               }
236             } else {
237               for (i=0;i<pcbddc->local_primal_size;i++) {
238                 array_to[i]=array_from[pcbddc->local_primal_indices[i]];
239               }
240             }
241           }
242           break;
243         case MULTILEVEL_BDDC:
244           break;
245         case PARALLEL_BDDC:
246           break;
247       }
248       if (vec_from) {
249         ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr);
250       }
251       if (vec_to) {
252         ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr);
253       }
254       break;
255   }
256   PetscFunctionReturn(0);
257 }
258 
259 /* uncomment for testing purposes */
260 /* #define PETSC_MISSING_LAPACK_GESVD 1 */
261 #undef __FUNCT__
262 #define __FUNCT__ "PCBDDCConstraintsSetUp"
263 PetscErrorCode PCBDDCConstraintsSetUp(PC pc)
264 {
265   PetscErrorCode    ierr;
266   PC_IS*            pcis = (PC_IS*)(pc->data);
267   PC_BDDC*          pcbddc = (PC_BDDC*)pc->data;
268   Mat_IS*           matis = (Mat_IS*)pc->pmat->data;
269   /* constraint and (optionally) change of basis matrix implemented as SeqAIJ */
270   MatType           impMatType=MATSEQAIJ;
271   /* one and zero */
272   PetscScalar       one=1.0,zero=0.0;
273   /* space to store constraints and their local indices */
274   PetscScalar       *temp_quadrature_constraint;
275   PetscInt          *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B;
276   /* iterators */
277   PetscInt          i,j,k,total_counts,temp_start_ptr;
278   /* stuff to store connected components stored in pcbddc->mat_graph */
279   IS                ISForVertices,*ISForFaces,*ISForEdges,*used_IS;
280   PetscInt          n_ISForFaces,n_ISForEdges;
281   PetscBool         get_faces,get_edges,get_vertices;
282   /* near null space stuff */
283   MatNullSpace      nearnullsp;
284   const Vec         *nearnullvecs;
285   Vec               *localnearnullsp;
286   PetscBool         nnsp_has_cnst;
287   PetscInt          nnsp_size;
288   PetscScalar       *array;
289   /* BLAS integers */
290   PetscBLASInt      Bs,Bt,lwork,lierr,Bone=1;
291   /* LAPACK working arrays for SVD or POD */
292   PetscScalar       *work;
293   PetscReal         *singular_vals;
294 #if defined(PETSC_USE_COMPLEX)
295   PetscReal         *rwork;
296 #endif
297 #if defined(PETSC_MISSING_LAPACK_GESVD)
298   PetscScalar       *temp_basis,*correlation_mat;
299 #endif
300   /* change of basis */
301   PetscInt          *aux_primal_numbering,*aux_primal_minloc,*global_indices;
302   PetscBool         boolforchange,*change_basis,*touched;
303   /* auxiliary stuff */
304   PetscInt          *nnz,*is_indices,*local_to_B;
305   /* some quantities */
306   PetscInt          n_vertices,total_primal_vertices;
307   PetscInt          size_of_constraint,max_size_of_constraint,max_constraints,temp_constraints;
308 
309 
310   PetscFunctionBegin;
311   /* Get index sets for faces, edges and vertices from graph */
312   get_faces = PETSC_TRUE;
313   get_edges = PETSC_TRUE;
314   get_vertices = PETSC_TRUE;
315   if (pcbddc->vertices_flag) {
316     get_faces = PETSC_FALSE;
317     get_edges = PETSC_FALSE;
318   }
319   if (pcbddc->constraints_flag) {
320     get_vertices = PETSC_FALSE;
321   }
322   if (pcbddc->faces_flag) {
323     get_edges = PETSC_FALSE;
324   }
325   if (pcbddc->edges_flag) {
326     get_faces = PETSC_FALSE;
327   }
328   /* default */
329   if (!get_faces && !get_edges && !get_vertices) {
330     get_vertices = PETSC_TRUE;
331   }
332   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
333   /* print some info */
334   if (pcbddc->dbg_flag) {
335     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
336     i = 0;
337     if (ISForVertices) {
338       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
339     }
340     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
341     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
342     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
343     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
344   }
345   /* check if near null space is attached to global mat */
346   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
347   if (nearnullsp) {
348     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
349   } else { /* if near null space is not provided BDDC uses constants by default */
350     nnsp_size = 0;
351     nnsp_has_cnst = PETSC_TRUE;
352   }
353   /* get max number of constraints on a single cc */
354   max_constraints = nnsp_size;
355   if (nnsp_has_cnst) max_constraints++;
356 
357   /*
358        Evaluate maximum storage size needed by the procedure
359        - temp_indices will contain start index of each constraint stored as follows
360        - temp_indices_to_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts
361        - temp_indices_to_constraint_B[temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in boundary numbering) on which the constraint acts
362        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
363                                                                                                                                                          */
364   total_counts = n_ISForFaces+n_ISForEdges;
365   total_counts *= max_constraints;
366   n_vertices = 0;
367   if (ISForVertices) {
368     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
369   }
370   total_counts += n_vertices;
371   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
372   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
373   total_counts = 0;
374   max_size_of_constraint = 0;
375   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
376     if (i<n_ISForEdges) {
377       used_IS = &ISForEdges[i];
378     } else {
379       used_IS = &ISForFaces[i-n_ISForEdges];
380     }
381     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
382     total_counts += j;
383     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
384   }
385   total_counts *= max_constraints;
386   total_counts += n_vertices;
387   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
388   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
389   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
390   /* local to boundary numbering */
391   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
392   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
393   for (i=0;i<pcis->n;i++) local_to_B[i]=-1;
394   for (i=0;i<pcis->n_B;i++) local_to_B[is_indices[i]]=i;
395   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
396   /* get local part of global near null space vectors */
397   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
398   for (k=0;k<nnsp_size;k++) {
399     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
400     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
401     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
402   }
403 
404   /* First we issue queries to allocate optimal workspace for LAPACKgesvd (or LAPACKsyev if SVD is missing) */
405   if (!pcbddc->use_nnsp_true) {
406     PetscScalar temp_work;
407 #if defined(PETSC_MISSING_LAPACK_GESVD)
408     /* Proper Orthogonal Decomposition (POD) using the snapshot method */
409     ierr = PetscMalloc(max_constraints*max_constraints*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
410     ierr = PetscMalloc(max_constraints*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
411     ierr = PetscMalloc(max_size_of_constraint*max_constraints*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
412 #if defined(PETSC_USE_COMPLEX)
413     ierr = PetscMalloc(3*max_constraints*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
414 #endif
415     /* now we evaluate the optimal workspace using query with lwork=-1 */
416     ierr = PetscBLASIntCast(max_constraints,&Bt);CHKERRQ(ierr);
417     lwork = -1;
418     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
419 #if !defined(PETSC_USE_COMPLEX)
420     PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr));
421 #else
422     PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr));
423 #endif
424     ierr = PetscFPTrapPop();CHKERRQ(ierr);
425     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEV Lapack routine %d",(int)lierr);
426 #else /* on missing GESVD */
427     /* SVD */
428     PetscInt max_n,min_n;
429     max_n = max_size_of_constraint;
430     min_n = max_constraints;
431     if (max_size_of_constraint < max_constraints) {
432       min_n = max_size_of_constraint;
433       max_n = max_constraints;
434     }
435     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
436 #if defined(PETSC_USE_COMPLEX)
437     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
438 #endif
439     /* now we evaluate the optimal workspace using query with lwork=-1 */
440     lwork = -1;
441     ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr);
442     ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr);
443     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
444 #if !defined(PETSC_USE_COMPLEX)
445     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,work,&Bt,work,&Bt,&temp_work,&lwork,&lierr));
446 #else
447     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,work,&Bt,work,&Bt,&temp_work,&lwork,rwork,&lierr));
448 #endif
449     ierr = PetscFPTrapPop();CHKERRQ(ierr);
450     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GESVD Lapack routine %d",(int)lierr);
451 #endif /* on missing GESVD */
452     /* Allocate optimal workspace */
453     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
454     ierr = PetscMalloc((PetscInt)lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
455   }
456   /* Now we can loop on constraining sets */
457   total_counts = 0;
458   temp_indices[0] = 0;
459   /* vertices */
460   if (ISForVertices) {
461     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
462     if (nnsp_has_cnst) { /* consider all vertices */
463       for (i=0;i<n_vertices;i++) {
464         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
465         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
466         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
467         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
468         change_basis[total_counts]=PETSC_FALSE;
469         total_counts++;
470       }
471     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
472       PetscBool used_vertex;
473       for (i=0;i<n_vertices;i++) {
474         used_vertex = PETSC_FALSE;
475         k = 0;
476         while (!used_vertex && k<nnsp_size) {
477           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
478           if (PetscAbsScalar(array[is_indices[i]])>0.0) {
479             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
480             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
481             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
482             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
483             change_basis[total_counts]=PETSC_FALSE;
484             total_counts++;
485             used_vertex = PETSC_TRUE;
486           }
487           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
488           k++;
489         }
490       }
491     }
492     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
493     n_vertices = total_counts;
494   }
495 
496   /* edges and faces */
497   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
498     if (i<n_ISForEdges) {
499       used_IS = &ISForEdges[i];
500       boolforchange = pcbddc->use_change_of_basis; /* change or not the basis on the edge */
501     } else {
502       used_IS = &ISForFaces[i-n_ISForEdges];
503       boolforchange = (PetscBool)(pcbddc->use_change_of_basis && pcbddc->use_change_on_faces); /* change or not the basis on the face */
504     }
505     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
506     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
507     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
508     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
509     /* change of basis should not be performed on local periodic nodes */
510     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) boolforchange = PETSC_FALSE;
511     if (nnsp_has_cnst) {
512       PetscScalar quad_value;
513       temp_constraints++;
514       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
515       for (j=0;j<size_of_constraint;j++) {
516         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
517         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
518         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
519       }
520       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
521       change_basis[total_counts]=boolforchange;
522       total_counts++;
523     }
524     for (k=0;k<nnsp_size;k++) {
525       PetscReal real_value;
526       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
527       for (j=0;j<size_of_constraint;j++) {
528         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
529         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
530         temp_quadrature_constraint[temp_indices[total_counts]+j]=array[is_indices[j]];
531       }
532       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
533       /* check if array is null on the connected component */
534       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
535       PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone));
536       if (real_value > 0.0) { /* keep indices and values */
537         temp_constraints++;
538         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
539         change_basis[total_counts]=boolforchange;
540         total_counts++;
541       }
542     }
543     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
544     /* perform SVD on the constraints if use_nnsp_true has not be requested by the user */
545     if (!pcbddc->use_nnsp_true) {
546       PetscReal tol = 1.0e-8; /* tolerance for retaining eigenmodes */
547 
548       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
549       ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr);
550 #if defined(PETSC_MISSING_LAPACK_GESVD)
551       /* SVD: Y = U*S*V^H                -> U (eigenvectors of Y*Y^H) = Y*V*(S)^\dag
552          POD: Y^H*Y = V*D*V^H, D = S^H*S -> U = Y*V*D^(-1/2)
553          -> When PETSC_USE_COMPLEX and PETSC_MISSING_LAPACK_GESVD are defined
554             the constraints basis will differ (by a complex factor with absolute value equal to 1)
555             from that computed using LAPACKgesvd
556          -> This is due to a different computation of eigenvectors in LAPACKheev
557          -> The quality of the POD-computed basis will be the same */
558       ierr = PetscMemzero(correlation_mat,temp_constraints*temp_constraints*sizeof(PetscScalar));CHKERRQ(ierr);
559       /* Store upper triangular part of correlation matrix */
560       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
561       for (j=0;j<temp_constraints;j++) {
562         for (k=0;k<j+1;k++) {
563           PetscStackCallBLAS("BLASdot",correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone));
564         }
565       }
566 #if !defined(PETSC_USE_COMPLEX)
567       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr));
568 #else
569       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,rwork,&lierr));
570 #endif
571       ierr = PetscFPTrapPop();CHKERRQ(ierr);
572       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEV Lapack routine %d",(int)lierr);
573       /* retain eigenvalues greater than tol: note that LAPACKsyev gives eigs in ascending order */
574       j=0;
575       while (j < temp_constraints && singular_vals[j] < tol) j++;
576       total_counts=total_counts-j;
577       if (j<temp_constraints) {
578         PetscInt ii;
579         for (k=j;k<temp_constraints;k++) singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
580         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
581         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,correlation_mat,&Bt,&zero,temp_basis,&Bs));
582         ierr = PetscFPTrapPop();CHKERRQ(ierr);
583         /* scale and copy POD basis into used quadrature memory */
584         for (k=0;k<temp_constraints-j;k++) {
585           for (ii=0;ii<size_of_constraint;ii++) {
586             temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[temp_constraints-1-k]*temp_basis[(temp_constraints-1-k)*size_of_constraint+ii];
587           }
588         }
589       }
590 #else  /* on missing GESVD */
591       PetscInt min_n = temp_constraints;
592       if (min_n > size_of_constraint) min_n = size_of_constraint;
593       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
594 #if !defined(PETSC_USE_COMPLEX)
595       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,work,&Bt,work,&Bt,work,&lwork,&lierr));
596 #else
597       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,work,&Bt,work,&Bt,work,&lwork,rwork,&lierr));
598 #endif
599       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESVD Lapack routine %d",(int)lierr);
600       ierr = PetscFPTrapPop();CHKERRQ(ierr);
601       /* retain eigenvalues greater than tol: note that LAPACKgesvd gives eigs in descending order */
602       j = 0;
603       while (j < min_n && singular_vals[min_n-j-1] < tol) j++;
604       total_counts = total_counts-temp_constraints+min_n-j;
605 #endif /* on missing GESVD */
606     }
607   }
608   /* free index sets of faces, edges and vertices */
609   for (i=0;i<n_ISForFaces;i++) {
610     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
611   }
612   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
613   for (i=0;i<n_ISForEdges;i++) {
614     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
615   }
616   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
617   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
618 
619   /* free workspace */
620   if (!pcbddc->use_nnsp_true) {
621     ierr = PetscFree(work);CHKERRQ(ierr);
622 #if defined(PETSC_USE_COMPLEX)
623     ierr = PetscFree(rwork);CHKERRQ(ierr);
624 #endif
625     ierr = PetscFree(singular_vals);CHKERRQ(ierr);
626 #if defined(PETSC_MISSING_LAPACK_GESVD)
627     ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
628     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
629 #endif
630   }
631   for (k=0;k<nnsp_size;k++) {
632     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
633   }
634   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
635 
636   /* set quantities in pcbddc data structure */
637   /* n_vertices defines the number of subdomain corners in the primal space */
638   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
639   pcbddc->local_primal_size = total_counts;
640   pcbddc->n_vertices = n_vertices;
641   pcbddc->n_constraints = pcbddc->local_primal_size-pcbddc->n_vertices;
642 
643   /* Create constraint matrix */
644   /* The constraint matrix is used to compute the l2g map of primal dofs */
645   /* so we need to set it up properly either with or without change of basis */
646   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
647   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
648   ierr = MatSetSizes(pcbddc->ConstraintMatrix,pcbddc->local_primal_size,pcis->n,pcbddc->local_primal_size,pcis->n);CHKERRQ(ierr);
649   /* array to compute a local numbering of constraints : vertices first then constraints */
650   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
651   /* array to select the proper local node (of minimum index with respect to global ordering) when changing the basis */
652   /* note: it should not be needed since IS for faces and edges are already sorted by global ordering when analyzing the graph but... just in case */
653   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_minloc);CHKERRQ(ierr);
654   /* auxiliary stuff for basis change */
655   ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&global_indices);CHKERRQ(ierr);
656   ierr = PetscMalloc(pcis->n_B*sizeof(PetscBool),&touched);CHKERRQ(ierr);
657   ierr = PetscMemzero(touched,pcis->n_B*sizeof(PetscBool));CHKERRQ(ierr);
658 
659   /* find primal_dofs: subdomain corners plus dofs selected as primal after change of basis */
660   total_primal_vertices=0;
661   for (i=0;i<pcbddc->local_primal_size;i++) {
662     size_of_constraint=temp_indices[i+1]-temp_indices[i];
663     if (size_of_constraint == 1) {
664       touched[temp_indices_to_constraint_B[temp_indices[i]]]=PETSC_TRUE;
665       aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]];
666       aux_primal_minloc[total_primal_vertices]=0;
667       total_primal_vertices++;
668     } else if (change_basis[i]) { /* Same procedure used in PCBDDCGetPrimalConstraintsLocalIdx */
669       PetscInt min_loc,min_index;
670       ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],global_indices);CHKERRQ(ierr);
671       /* find first untouched local node */
672       k = 0;
673       while (touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) k++;
674       min_index = global_indices[k];
675       min_loc = k;
676       /* search the minimum among global nodes already untouched on the cc */
677       for (k=1;k<size_of_constraint;k++) {
678         /* there can be more than one constraint on a single connected component */
679         if (min_index > global_indices[k] && !touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) {
680           min_index = global_indices[k];
681           min_loc = k;
682         }
683       }
684       touched[temp_indices_to_constraint_B[temp_indices[i]+min_loc]] = PETSC_TRUE;
685       aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]+min_loc];
686       aux_primal_minloc[total_primal_vertices]=min_loc;
687       total_primal_vertices++;
688     }
689   }
690   /* free workspace */
691   ierr = PetscFree(global_indices);CHKERRQ(ierr);
692   ierr = PetscFree(touched);CHKERRQ(ierr);
693   /* permute indices in order to have a sorted set of vertices */
694   ierr = PetscSortInt(total_primal_vertices,aux_primal_numbering);
695 
696   /* nonzero structure of constraint matrix */
697   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
698   for (i=0;i<total_primal_vertices;i++) nnz[i]=1;
699   j=total_primal_vertices;
700   for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
701     if (!change_basis[i]) {
702       nnz[j]=temp_indices[i+1]-temp_indices[i];
703       j++;
704     }
705   }
706   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
707   ierr = PetscFree(nnz);CHKERRQ(ierr);
708   /* set values in constraint matrix */
709   for (i=0;i<total_primal_vertices;i++) {
710     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,aux_primal_numbering[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
711   }
712   total_counts = total_primal_vertices;
713   for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
714     if (!change_basis[i]) {
715       size_of_constraint=temp_indices[i+1]-temp_indices[i];
716       ierr = MatSetValues(pcbddc->ConstraintMatrix,1,&total_counts,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],&temp_quadrature_constraint[temp_indices[i]],INSERT_VALUES);CHKERRQ(ierr);
717       total_counts++;
718     }
719   }
720   /* assembling */
721   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
722   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
723   /*
724   ierr = MatView(pcbddc->ConstraintMatrix,(PetscViewer)0);CHKERRQ(ierr);
725   */
726   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
727   if (pcbddc->use_change_of_basis) {
728     PetscBool qr_needed = PETSC_FALSE;
729     /* change of basis acts on local interfaces -> dimension is n_B x n_B */
730     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
731     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
732     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
733     /* work arrays */
734     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
735     for (i=0;i<pcis->n_B;i++) nnz[i]=1;
736     /* nonzeros per row */
737     for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
738       if (change_basis[i]) {
739         qr_needed = PETSC_TRUE;
740         size_of_constraint = temp_indices[i+1]-temp_indices[i];
741         for (j=0;j<size_of_constraint;j++) nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
742       }
743     }
744     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
745     ierr = PetscFree(nnz);CHKERRQ(ierr);
746     /* Set initial identity in the matrix */
747     for (i=0;i<pcis->n_B;i++) {
748       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
749     }
750 
751     /* Now we loop on the constraints which need a change of basis */
752     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
753     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (see Sect 6.2.1) */
754     /* Change of basis matrix T computed via QR decomposition of constraints */
755     if (qr_needed) {
756       /* dual and primal dofs on a single cc */
757       PetscInt     dual_dofs,primal_dofs;
758       /* iterator on aux_primal_minloc (ordered as read from nearnullspace: vertices, edges and then constraints) */
759       PetscInt     primal_counter;
760       /* working stuff for GEQRF */
761       PetscScalar  *qr_basis,*qr_tau,*qr_work,lqr_work_t;
762       PetscBLASInt lqr_work;
763       /* working stuff for UNGQR */
764       PetscScalar  *gqr_work,lgqr_work_t;
765       PetscBLASInt lgqr_work;
766       /* working stuff for TRTRS */
767       PetscScalar  *trs_rhs;
768       /* pointers for values insertion into change of basis matrix */
769       PetscInt     *start_rows,*start_cols;
770       PetscScalar  *start_vals;
771       /* working stuff for values insertion */
772       PetscBool    *is_primal;
773 
774       /* space to store Q */
775       ierr = PetscMalloc((max_size_of_constraint)*(max_size_of_constraint)*sizeof(PetscScalar),&qr_basis);CHKERRQ(ierr);
776       /* first we issue queries for optimal work */
777       ierr = PetscBLASIntCast(max_size_of_constraint,&Bs);CHKERRQ(ierr);
778       ierr = PetscBLASIntCast(max_constraints,&Bt);CHKERRQ(ierr);
779       lqr_work = -1;
780       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Bs,&Bt,qr_basis,&Bs,qr_tau,&lqr_work_t,&lqr_work,&lierr));
781       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GEQRF Lapack routine %d",(int)lierr);
782       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lqr_work_t),&lqr_work);CHKERRQ(ierr);
783       ierr = PetscMalloc((PetscInt)PetscRealPart(lqr_work_t)*sizeof(*qr_work),&qr_work);CHKERRQ(ierr);
784       lgqr_work = -1;
785       if (Bt>Bs) Bt=Bs; /* adjust Bt just for computing optimal work */
786       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Bs,&Bs,&Bt,qr_basis,&Bs,qr_tau,&lgqr_work_t,&lgqr_work,&lierr));
787       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to UNGQR Lapack routine %d",(int)lierr);
788       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lgqr_work_t),&lgqr_work);CHKERRQ(ierr);
789       ierr = PetscMalloc((PetscInt)PetscRealPart(lgqr_work_t)*sizeof(*gqr_work),&gqr_work);CHKERRQ(ierr);
790       /* array to store scaling factors for reflectors */
791       ierr = PetscMalloc(max_constraints*sizeof(*qr_tau),&qr_tau);CHKERRQ(ierr);
792       /* array to store rhs and solution of triangular solver */
793       ierr = PetscMalloc(max_constraints*max_constraints*sizeof(*trs_rhs),&trs_rhs);CHKERRQ(ierr);
794       /* array to store whether a node is primal or not */
795       ierr = PetscMalloc(pcis->n_B*sizeof(*is_primal),&is_primal);CHKERRQ(ierr);
796       ierr = PetscMemzero(is_primal,pcis->n_B*sizeof(*is_primal));CHKERRQ(ierr);
797       for (i=0;i<total_primal_vertices;i++) is_primal[local_to_B[aux_primal_numbering[i]]] = PETSC_TRUE;
798 
799       /* allocating workspace for check */
800       if (pcbddc->dbg_flag) {
801         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
802         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Checking change of basis computation for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
803         ierr = PetscMalloc(max_size_of_constraint*(max_constraints+max_size_of_constraint)*sizeof(*work),&work);CHKERRQ(ierr);
804       }
805 
806       /* loop on constraints and see whether or not they need a change of basis */
807       /* -> using implicit ordering contained in temp_indices data */
808       total_counts = pcbddc->n_vertices;
809       primal_counter = total_counts;
810       while (total_counts<pcbddc->local_primal_size) {
811         primal_dofs = 1;
812         if (change_basis[total_counts]) {
813           /* get all constraints with same support: if more then one constraint is present on the cc then surely indices are stored contiguosly */
814           while (total_counts+primal_dofs < pcbddc->local_primal_size && temp_indices_to_constraint_B[temp_indices[total_counts]] == temp_indices_to_constraint_B[temp_indices[total_counts+primal_dofs]]) {
815             primal_dofs++;
816           }
817           /* get constraint info */
818           size_of_constraint = temp_indices[total_counts+1]-temp_indices[total_counts];
819           dual_dofs = size_of_constraint-primal_dofs;
820           /* get BLAS dims */
821           ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
822           ierr = PetscBLASIntCast(primal_dofs,&Bt);CHKERRQ(ierr);
823 
824           /* copy quadrature constraints for change of basis check */
825           if (pcbddc->dbg_flag) {
826             ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Constraint %d to %d need a change of basis (size %d)\n",total_counts,total_counts+primal_dofs,size_of_constraint);CHKERRQ(ierr);
827             ierr = PetscMemcpy(work,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
828           }
829 
830           /* copy temporary constraints into larger work vector (in order to store all columns of Q) */
831           ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
832 
833           /* compute QR decomposition of constraints */
834           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
835           PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Bs,&Bt,qr_basis,&Bs,qr_tau,qr_work,&lqr_work,&lierr));
836           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GEQRF Lapack routine %d",(int)lierr);
837           ierr = PetscFPTrapPop();CHKERRQ(ierr);
838 
839           /* explictly compute R^-T */
840           ierr = PetscMemzero(trs_rhs,primal_dofs*primal_dofs*sizeof(*trs_rhs));CHKERRQ(ierr);
841           for (j=0;j<primal_dofs;j++) trs_rhs[j*(primal_dofs+1)] = 1.0;
842           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
843           PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&Bt,&Bt,qr_basis,&Bs,trs_rhs,&Bt,&lierr));
844           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TRTRS Lapack routine %d",(int)lierr);
845           ierr = PetscFPTrapPop();CHKERRQ(ierr);
846 
847           /* explcitly compute all columns of Q (Q = [Q1 | Q2] ) overwriting QR factorization in qr_basis */
848           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
849           PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Bs,&Bs,&Bt,qr_basis,&Bs,qr_tau,gqr_work,&lgqr_work,&lierr));
850           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in UNGQR Lapack routine %d",(int)lierr);
851           ierr = PetscFPTrapPop();CHKERRQ(ierr);
852 
853           /* first primal_dofs columns of Q need to be re-scaled in order to be unitary w.r.t constraints
854              i.e. C_{pxn}*Q_{nxn} should be equal to [I_pxp | 0_pxd] (see check below)
855              where n=size_of_constraint, p=primal_dofs, d=dual_dofs (n=p+d), I and 0 identity and null matrix resp. */
856           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
857           PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,qr_basis,&Bs,trs_rhs,&Bt,&zero,&temp_quadrature_constraint[temp_indices[total_counts]],&Bs));
858           ierr = PetscFPTrapPop();CHKERRQ(ierr);
859           ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
860 
861           /* insert values in change of basis matrix respecting global ordering of new primal dofs */
862           start_rows = &temp_indices_to_constraint_B[temp_indices[total_counts]];
863           /* insert cols for primal dofs */
864           for (j=0;j<primal_dofs;j++) {
865             start_vals = &qr_basis[j*size_of_constraint];
866             start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+aux_primal_minloc[primal_counter+j]];
867             ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr);
868           }
869           /* insert cols for dual dofs */
870           for (j=0,k=0;j<dual_dofs;k++) {
871             if (!is_primal[temp_indices_to_constraint_B[temp_indices[total_counts]+k]]) {
872               start_vals = &qr_basis[(primal_dofs+j)*size_of_constraint];
873               start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+k];
874               ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr);
875               j++;
876             }
877           }
878 
879           /* check change of basis */
880           if (pcbddc->dbg_flag) {
881             PetscInt   ii,jj;
882             PetscBool valid_qr=PETSC_TRUE;
883             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
884             PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&Bt,&Bs,&Bs,&one,work,&Bs,qr_basis,&Bs,&zero,&work[size_of_constraint*primal_dofs],&Bt));
885             ierr = PetscFPTrapPop();CHKERRQ(ierr);
886             for (jj=0;jj<size_of_constraint;jj++) {
887               for (ii=0;ii<primal_dofs;ii++) {
888                 if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) valid_qr = PETSC_FALSE;
889                 if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) valid_qr = PETSC_FALSE;
890               }
891             }
892             if (!valid_qr) {
893               ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> wrong change of basis!\n",PetscGlobalRank);CHKERRQ(ierr);
894               for (jj=0;jj<size_of_constraint;jj++) {
895                 for (ii=0;ii<primal_dofs;ii++) {
896                   if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) {
897                     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\tQr basis function %d is not orthogonal to constraint %d (%1.14e)!\n",jj,ii,PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]));
898                   }
899                   if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) {
900                     PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\tQr basis function %d is not unitary w.r.t constraint %d (%1.14e)!\n",jj,ii,PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]));
901                   }
902                 }
903               }
904             } else {
905               ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> right change of basis!\n",PetscGlobalRank);CHKERRQ(ierr);
906             }
907           }
908           /* increment primal counter */
909           primal_counter += primal_dofs;
910         } else {
911           if (pcbddc->dbg_flag) {
912             ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Constraint %d does not need a change of basis (size %d)\n",total_counts,temp_indices[total_counts+1]-temp_indices[total_counts]);CHKERRQ(ierr);
913           }
914         }
915         /* increment constraint counter total_counts */
916         total_counts += primal_dofs;
917       }
918       if (pcbddc->dbg_flag) {
919         ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
920         ierr = PetscFree(work);CHKERRQ(ierr);
921       }
922       /* free workspace */
923       ierr = PetscFree(trs_rhs);CHKERRQ(ierr);
924       ierr = PetscFree(qr_tau);CHKERRQ(ierr);
925       ierr = PetscFree(qr_work);CHKERRQ(ierr);
926       ierr = PetscFree(gqr_work);CHKERRQ(ierr);
927       ierr = PetscFree(is_primal);CHKERRQ(ierr);
928       ierr = PetscFree(qr_basis);CHKERRQ(ierr);
929     }
930     /* assembling */
931     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
932     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
933     /*
934     ierr = MatView(pcbddc->ChangeOfBasisMatrix,(PetscViewer)0);CHKERRQ(ierr);
935     */
936   }
937   /* free workspace no longer needed */
938   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
939   ierr = PetscFree(aux_primal_minloc);CHKERRQ(ierr);
940   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
941   ierr = PetscFree(change_basis);CHKERRQ(ierr);
942   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
943   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
944   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
945   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
946   PetscFunctionReturn(0);
947 }
948 
949 #undef __FUNCT__
950 #define __FUNCT__ "PCBDDCAnalyzeInterface"
951 PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
952 {
953   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
954   PC_IS       *pcis = (PC_IS*)pc->data;
955   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
956   PetscInt    bs,ierr,i,vertex_size;
957   PetscViewer viewer=pcbddc->dbg_viewer;
958 
959   PetscFunctionBegin;
960   /* Init local Graph struct */
961   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
962 
963   /* Check validity of the csr graph passed in by the user */
964   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
965     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
966   }
967   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
968   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
969     Mat mat_adj;
970     const PetscInt *xadj,*adjncy;
971     PetscBool flg_row=PETSC_TRUE;
972 
973     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
974     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
975     if (!flg_row) {
976       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
977     }
978     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
979     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
980     if (!flg_row) {
981       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
982     }
983     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
984   }
985 
986   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
987   vertex_size = 1;
988   if (!pcbddc->n_ISForDofs) {
989     IS *custom_ISForDofs;
990 
991     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
992     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
993     for (i=0;i<bs;i++) {
994       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
995     }
996     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
997     /* remove my references to IS objects */
998     for (i=0;i<bs;i++) {
999       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
1000     }
1001     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
1002   } else { /* mat block size as vertex size (used for elasticity) */
1003     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
1004   }
1005 
1006   /* Setup of Graph */
1007   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
1008 
1009   /* Graph's connected components analysis */
1010   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
1011 
1012   /* print some info to stdout */
1013   if (pcbddc->dbg_flag) {
1014     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
1015   }
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 #undef __FUNCT__
1020 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
1021 PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
1022 {
1023   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1024   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
1025   PetscErrorCode ierr;
1026 
1027   PetscFunctionBegin;
1028   n = 0;
1029   vertices = 0;
1030   if (pcbddc->ConstraintMatrix) {
1031     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
1032     for (i=0;i<local_primal_size;i++) {
1033       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1034       if (size_of_constraint == 1) n++;
1035       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1036     }
1037     ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
1038     n = 0;
1039     for (i=0;i<local_primal_size;i++) {
1040       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1041       if (size_of_constraint == 1) {
1042         vertices[n++]=row_cmat_indices[0];
1043       }
1044       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1045     }
1046   }
1047   *n_vertices = n;
1048   *vertices_idx = vertices;
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 #undef __FUNCT__
1053 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
1054 PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
1055 {
1056   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1057   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
1058   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
1059   PetscBool      *touched;
1060   PetscErrorCode ierr;
1061 
1062   PetscFunctionBegin;
1063   n = 0;
1064   constraints_index = 0;
1065   if (pcbddc->ConstraintMatrix) {
1066     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
1067     max_size_of_constraint = 0;
1068     for (i=0;i<local_primal_size;i++) {
1069       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1070       if (size_of_constraint > 1) {
1071         n++;
1072       }
1073       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1074       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1075     }
1076     ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1077     ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1078     ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1079     ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1080     n = 0;
1081     for (i=0;i<local_primal_size;i++) {
1082       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1083       if (size_of_constraint > 1) {
1084         ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
1085         min_index = row_cmat_global_indices[0];
1086         min_loc = 0;
1087         for (j=1;j<size_of_constraint;j++) {
1088           /* there can be more than one constraint on a single connected component */
1089           if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1090             min_index = row_cmat_global_indices[j];
1091             min_loc = j;
1092           }
1093         }
1094         touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1095         constraints_index[n++] = row_cmat_indices[min_loc];
1096       }
1097       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1098     }
1099   }
1100   ierr = PetscFree(touched);CHKERRQ(ierr);
1101   ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1102   *n_constraints = n;
1103   *constraints_idx = constraints_index;
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 /* the next two functions has been adapted from pcis.c */
1108 #undef __FUNCT__
1109 #define __FUNCT__ "PCBDDCApplySchur"
1110 PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1111 {
1112   PetscErrorCode ierr;
1113   PC_IS          *pcis = (PC_IS*)(pc->data);
1114 
1115   PetscFunctionBegin;
1116   if (!vec2_B) { vec2_B = v; }
1117   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1118   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1119   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1120   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1121   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 #undef __FUNCT__
1126 #define __FUNCT__ "PCBDDCApplySchurTranspose"
1127 PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1128 {
1129   PetscErrorCode ierr;
1130   PC_IS          *pcis = (PC_IS*)(pc->data);
1131 
1132   PetscFunctionBegin;
1133   if (!vec2_B) { vec2_B = v; }
1134   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1135   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1136   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1137   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1138   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "PCBDDCSubsetNumbering"
1144 PetscErrorCode PCBDDCSubsetNumbering(MPI_Comm comm,ISLocalToGlobalMapping l2gmap, PetscInt n_local_dofs, PetscInt local_dofs[], PetscInt local_dofs_mult[], PetscInt* n_global_subset, PetscInt* global_numbering_subset[])
1145 {
1146   Vec            local_vec,global_vec;
1147   IS             seqis,paris;
1148   VecScatter     scatter_ctx;
1149   PetscScalar    *array;
1150   PetscInt       *temp_global_dofs;
1151   PetscScalar    globalsum;
1152   PetscInt       i,j,s;
1153   PetscInt       nlocals,first_index,old_index,max_local;
1154   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1155   PetscMPIInt    *dof_sizes,*dof_displs;
1156   PetscBool      first_found;
1157   PetscErrorCode ierr;
1158 
1159   PetscFunctionBegin;
1160   /* mpi buffers */
1161   MPI_Comm_size(comm,&size_prec_comm);
1162   MPI_Comm_rank(comm,&rank_prec_comm);
1163   j = ( !rank_prec_comm ? size_prec_comm : 0);
1164   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1165   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1166   /* get maximum size of subset */
1167   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1168   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1169   max_local = 0;
1170   if (n_local_dofs) {
1171     max_local = temp_global_dofs[0];
1172     for (i=1;i<n_local_dofs;i++) {
1173       if (max_local < temp_global_dofs[i] ) {
1174         max_local = temp_global_dofs[i];
1175       }
1176     }
1177   }
1178   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1179   max_global++;
1180   max_local = 0;
1181   if (n_local_dofs) {
1182     max_local = local_dofs[0];
1183     for (i=1;i<n_local_dofs;i++) {
1184       if (max_local < local_dofs[i] ) {
1185         max_local = local_dofs[i];
1186       }
1187     }
1188   }
1189   max_local++;
1190   /* allocate workspace */
1191   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1192   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1193   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1194   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1195   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1196   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1197   /* create scatter */
1198   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1199   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1200   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1201   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1202   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1203   /* init array */
1204   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1205   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1206   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1207   if (local_dofs_mult) {
1208     for (i=0;i<n_local_dofs;i++) {
1209       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1210     }
1211   } else {
1212     for (i=0;i<n_local_dofs;i++) {
1213       array[local_dofs[i]]=1.0;
1214     }
1215   }
1216   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1217   /* scatter into global vec and get total number of global dofs */
1218   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1219   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1220   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
1221   *n_global_subset = (PetscInt)PetscRealPart(globalsum);
1222   /* Fill global_vec with cumulative function for global numbering */
1223   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1224   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1225   nlocals = 0;
1226   first_index = -1;
1227   first_found = PETSC_FALSE;
1228   for (i=0;i<s;i++) {
1229     if (!first_found && PetscRealPart(array[i]) > 0.0) {
1230       first_found = PETSC_TRUE;
1231       first_index = i;
1232     }
1233     nlocals += (PetscInt)PetscRealPart(array[i]);
1234   }
1235   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1236   if (!rank_prec_comm) {
1237     dof_displs[0]=0;
1238     for (i=1;i<size_prec_comm;i++) {
1239       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1240     }
1241   }
1242   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1243   if (first_found) {
1244     array[first_index] += (PetscScalar)nlocals;
1245     old_index = first_index;
1246     for (i=first_index+1;i<s;i++) {
1247       if (PetscRealPart(array[i]) > 0.0) {
1248         array[i] += array[old_index];
1249         old_index = i;
1250       }
1251     }
1252   }
1253   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1254   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1255   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1256   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1257   /* get global ordering of local dofs */
1258   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1259   if (local_dofs_mult) {
1260     for (i=0;i<n_local_dofs;i++) {
1261       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i];
1262     }
1263   } else {
1264     for (i=0;i<n_local_dofs;i++) {
1265       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1;
1266     }
1267   }
1268   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1269   /* free workspace */
1270   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1271   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1272   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1273   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1274   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1275   /* return pointer to global ordering of local dofs */
1276   *global_numbering_subset = temp_global_dofs;
1277   PetscFunctionReturn(0);
1278 }
1279