xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision 3f08241adfa0b18368962bb20d212635c956b227)
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   PetscBLASInt      Blas_N,Blas_M,Blas_K;
292   PetscBLASInt      Blas_LDA,Blas_LDB,Blas_LDC;
293   /* LAPACK working arrays for SVD or POD */
294   PetscBool         skip_lapack;
295   PetscScalar       *work;
296   PetscReal         *singular_vals;
297 #if defined(PETSC_USE_COMPLEX)
298   PetscReal         *rwork;
299 #endif
300 #if defined(PETSC_MISSING_LAPACK_GESVD)
301   PetscBLASInt      Bone_2=1;
302   PetscScalar       *temp_basis,*correlation_mat;
303 #else
304   PetscBLASInt      dummy_int_1=1,dummy_int_2=1;
305   PetscScalar       dummy_scalar_1=0.0,dummy_scalar_2=0.0;
306 #endif
307   /* change of basis */
308   PetscInt          *aux_primal_numbering,*aux_primal_minloc,*global_indices;
309   PetscBool         boolforchange,*change_basis,*touched;
310   /* auxiliary stuff */
311   PetscInt          *nnz,*is_indices,*local_to_B;
312   /* some quantities */
313   PetscInt          n_vertices,total_primal_vertices;
314   PetscInt          size_of_constraint,max_size_of_constraint,max_constraints,temp_constraints;
315 
316 
317   PetscFunctionBegin;
318   /* Get index sets for faces, edges and vertices from graph */
319   get_faces = PETSC_TRUE;
320   get_edges = PETSC_TRUE;
321   get_vertices = PETSC_TRUE;
322   if (pcbddc->vertices_flag) {
323     get_faces = PETSC_FALSE;
324     get_edges = PETSC_FALSE;
325   }
326   if (pcbddc->constraints_flag) {
327     get_vertices = PETSC_FALSE;
328   }
329   if (pcbddc->faces_flag) {
330     get_edges = PETSC_FALSE;
331   }
332   if (pcbddc->edges_flag) {
333     get_faces = PETSC_FALSE;
334   }
335   /* default */
336   if (!get_faces && !get_edges && !get_vertices) {
337     get_vertices = PETSC_TRUE;
338   }
339   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
340   /* print some info */
341   if (pcbddc->dbg_flag) {
342     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
343     i = 0;
344     if (ISForVertices) {
345       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
346     }
347     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
348     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
349     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
350     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
351   }
352   /* check if near null space is attached to global mat */
353   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
354   if (nearnullsp) {
355     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
356   } else { /* if near null space is not provided BDDC uses constants by default */
357     nnsp_size = 0;
358     nnsp_has_cnst = PETSC_TRUE;
359   }
360   /* get max number of constraints on a single cc */
361   max_constraints = nnsp_size;
362   if (nnsp_has_cnst) max_constraints++;
363 
364   /*
365        Evaluate maximum storage size needed by the procedure
366        - temp_indices will contain start index of each constraint stored as follows
367        - temp_indices_to_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts
368        - 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
369        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
370                                                                                                                                                          */
371   total_counts = n_ISForFaces+n_ISForEdges;
372   total_counts *= max_constraints;
373   n_vertices = 0;
374   if (ISForVertices) {
375     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
376   }
377   total_counts += n_vertices;
378   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
379   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
380   total_counts = 0;
381   max_size_of_constraint = 0;
382   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
383     if (i<n_ISForEdges) {
384       used_IS = &ISForEdges[i];
385     } else {
386       used_IS = &ISForFaces[i-n_ISForEdges];
387     }
388     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
389     total_counts += j;
390     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
391   }
392   total_counts *= max_constraints;
393   total_counts += n_vertices;
394   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
395   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
396   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
397   /* local to boundary numbering */
398   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
399   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
400   for (i=0;i<pcis->n;i++) local_to_B[i]=-1;
401   for (i=0;i<pcis->n_B;i++) local_to_B[is_indices[i]]=i;
402   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
403   /* get local part of global near null space vectors */
404   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
405   for (k=0;k<nnsp_size;k++) {
406     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
407     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
408     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
409   }
410 
411   /* whether or not to skip lapack calls */
412   skip_lapack = PETSC_TRUE;
413   if (n_ISForFaces+n_ISForEdges) skip_lapack = PETSC_FALSE;
414 
415   /* First we issue queries to allocate optimal workspace for LAPACKgesvd (or LAPACKsyev if SVD is missing) */
416   if (!pcbddc->use_nnsp_true && !skip_lapack) {
417     PetscScalar temp_work;
418 #if defined(PETSC_MISSING_LAPACK_GESVD)
419     /* Proper Orthogonal Decomposition (POD) using the snapshot method */
420     ierr = PetscMalloc(max_constraints*max_constraints*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
421     ierr = PetscMalloc(max_constraints*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
422     ierr = PetscMalloc(max_size_of_constraint*max_constraints*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
423 #if defined(PETSC_USE_COMPLEX)
424     ierr = PetscMalloc(3*max_constraints*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
425 #endif
426     /* now we evaluate the optimal workspace using query with lwork=-1 */
427     ierr = PetscBLASIntCast(max_constraints,&Blas_N);CHKERRQ(ierr);
428     ierr = PetscBLASIntCast(max_constraints,&Blas_LDA);CHKERRQ(ierr);
429     lwork = -1;
430     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
431 #if !defined(PETSC_USE_COMPLEX)
432     PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,&temp_work,&lwork,&lierr));
433 #else
434     PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,&temp_work,&lwork,rwork,&lierr));
435 #endif
436     ierr = PetscFPTrapPop();CHKERRQ(ierr);
437     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEV Lapack routine %d",(int)lierr);
438 #else /* on missing GESVD */
439     /* SVD */
440     PetscInt max_n,min_n;
441     max_n = max_size_of_constraint;
442     min_n = max_constraints;
443     if (max_size_of_constraint < max_constraints) {
444       min_n = max_size_of_constraint;
445       max_n = max_constraints;
446     }
447     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
448 #if defined(PETSC_USE_COMPLEX)
449     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
450 #endif
451     /* now we evaluate the optimal workspace using query with lwork=-1 */
452     lwork = -1;
453     ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr);
454     ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr);
455     ierr = PetscBLASIntCast(max_n,&Blas_LDA);CHKERRQ(ierr);
456     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
457 #if !defined(PETSC_USE_COMPLEX)
458     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Blas_LDA,singular_vals,&dummy_scalar_1,&dummy_int_1,&dummy_scalar_2,&dummy_int_2,&temp_work,&lwork,&lierr));
459 #else
460     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Blas_LDA,singular_vals,&dummy_scalar_1,&dummy_int_1,&dummy_scalar_2,&dummy_int_2,&temp_work,&lwork,rwork,&lierr));
461 #endif
462     ierr = PetscFPTrapPop();CHKERRQ(ierr);
463     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GESVD Lapack routine %d",(int)lierr);
464 #endif /* on missing GESVD */
465     /* Allocate optimal workspace */
466     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
467     ierr = PetscMalloc((PetscInt)lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
468   }
469   /* Now we can loop on constraining sets */
470   total_counts = 0;
471   temp_indices[0] = 0;
472   /* vertices */
473   if (ISForVertices) {
474     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
475     if (nnsp_has_cnst) { /* consider all vertices */
476       for (i=0;i<n_vertices;i++) {
477         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
478         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
479         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
480         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
481         change_basis[total_counts]=PETSC_FALSE;
482         total_counts++;
483       }
484     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
485       PetscBool used_vertex;
486       for (i=0;i<n_vertices;i++) {
487         used_vertex = PETSC_FALSE;
488         k = 0;
489         while (!used_vertex && k<nnsp_size) {
490           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
491           if (PetscAbsScalar(array[is_indices[i]])>0.0) {
492             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
493             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
494             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
495             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
496             change_basis[total_counts]=PETSC_FALSE;
497             total_counts++;
498             used_vertex = PETSC_TRUE;
499           }
500           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
501           k++;
502         }
503       }
504     }
505     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
506     n_vertices = total_counts;
507   }
508 
509   /* edges and faces */
510   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
511     if (i<n_ISForEdges) {
512       used_IS = &ISForEdges[i];
513       boolforchange = pcbddc->use_change_of_basis; /* change or not the basis on the edge */
514     } else {
515       used_IS = &ISForFaces[i-n_ISForEdges];
516       boolforchange = (PetscBool)(pcbddc->use_change_of_basis && pcbddc->use_change_on_faces); /* change or not the basis on the face */
517     }
518     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
519     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
520     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
521     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
522     /* change of basis should not be performed on local periodic nodes */
523     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) boolforchange = PETSC_FALSE;
524     if (nnsp_has_cnst) {
525       PetscScalar quad_value;
526       temp_constraints++;
527       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
528       for (j=0;j<size_of_constraint;j++) {
529         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
530         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
531         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
532       }
533       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
534       change_basis[total_counts]=boolforchange;
535       total_counts++;
536     }
537     for (k=0;k<nnsp_size;k++) {
538       PetscReal real_value;
539       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
540       for (j=0;j<size_of_constraint;j++) {
541         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
542         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
543         temp_quadrature_constraint[temp_indices[total_counts]+j]=array[is_indices[j]];
544       }
545       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
546       /* check if array is null on the connected component */
547       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
548       PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone));
549       if (real_value > 0.0) { /* keep indices and values */
550         temp_constraints++;
551         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
552         change_basis[total_counts]=boolforchange;
553         total_counts++;
554       }
555     }
556     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
557     /* perform SVD on the constraints if use_nnsp_true has not be requested by the user */
558     if (!pcbddc->use_nnsp_true) {
559       PetscReal tol = 1.0e-8; /* tolerance for retaining eigenmodes */
560 
561       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
562       ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr);
563 #if defined(PETSC_MISSING_LAPACK_GESVD)
564       /* SVD: Y = U*S*V^H                -> U (eigenvectors of Y*Y^H) = Y*V*(S)^\dag
565          POD: Y^H*Y = V*D*V^H, D = S^H*S -> U = Y*V*D^(-1/2)
566          -> When PETSC_USE_COMPLEX and PETSC_MISSING_LAPACK_GESVD are defined
567             the constraints basis will differ (by a complex factor with absolute value equal to 1)
568             from that computed using LAPACKgesvd
569          -> This is due to a different computation of eigenvectors in LAPACKheev
570          -> The quality of the POD-computed basis will be the same */
571       ierr = PetscBLASIntCast(temp_constraints,&Blas_N);CHKERRQ(ierr);
572       ierr = PetscBLASIntCast(temp_constraints,&Blas_LDA);CHKERRQ(ierr);
573       ierr = PetscMemzero(correlation_mat,temp_constraints*temp_constraints*sizeof(PetscScalar));CHKERRQ(ierr);
574       /* Store upper triangular part of correlation matrix */
575       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
576       for (j=0;j<temp_constraints;j++) {
577         for (k=0;k<j+1;k++) {
578           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_2));
579         }
580       }
581 #if !defined(PETSC_USE_COMPLEX)
582       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,work,&lwork,&lierr));
583 #else
584       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,work,&lwork,rwork,&lierr));
585 #endif
586       ierr = PetscFPTrapPop();CHKERRQ(ierr);
587       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEV Lapack routine %d",(int)lierr);
588       /* retain eigenvalues greater than tol: note that LAPACKsyev gives eigs in ascending order */
589       j=0;
590       while (j < temp_constraints && singular_vals[j] < tol) j++;
591       total_counts=total_counts-j;
592       ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
593       ierr = PetscBLASIntCast(temp_constraints,&Blas_N);CHKERRQ(ierr);
594       ierr = PetscBLASIntCast(temp_constraints,&Blas_K);CHKERRQ(ierr);
595       ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
596       ierr = PetscBLASIntCast(temp_constraints,&Blas_LDB);CHKERRQ(ierr);
597       ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDC);CHKERRQ(ierr);
598       if (j<temp_constraints) {
599         PetscInt ii;
600         for (k=j;k<temp_constraints;k++) singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
601         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
602         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&Blas_M,&Blas_N,&Blas_K,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Blas_LDA,correlation_mat,&Blas_LDB,&zero,temp_basis,&Blas_LDC));
603         ierr = PetscFPTrapPop();CHKERRQ(ierr);
604         /* scale and copy POD basis into used quadrature memory */
605         for (k=0;k<temp_constraints-j;k++) {
606           for (ii=0;ii<size_of_constraint;ii++) {
607             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];
608           }
609         }
610       }
611 #else  /* on missing GESVD */
612       ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
613       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
614 #if !defined(PETSC_USE_COMPLEX)
615       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Blas_LDA,singular_vals,&dummy_scalar_1,&dummy_int_1,&dummy_scalar_2,&dummy_int_2,work,&lwork,&lierr));
616 #else
617       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Blas_LDA,singular_vals,&dummy_scalar_1,&dummy_int_1,&dummy_scalar_2,&dummy_int_2,work,&lwork,rwork,&lierr));
618 #endif
619       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESVD Lapack routine %d",(int)lierr);
620       ierr = PetscFPTrapPop();CHKERRQ(ierr);
621       /* retain eigenvalues greater than tol: note that LAPACKgesvd gives eigs in descending order */
622       PetscInt min_n = temp_constraints;
623       if (min_n > size_of_constraint) min_n = size_of_constraint;
624       j = 0;
625       while (j < min_n && singular_vals[min_n-j-1] < tol) j++;
626       total_counts = total_counts-temp_constraints+min_n-j;
627 #endif /* on missing GESVD */
628     }
629   }
630   /* free index sets of faces, edges and vertices */
631   for (i=0;i<n_ISForFaces;i++) {
632     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
633   }
634   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
635   for (i=0;i<n_ISForEdges;i++) {
636     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
637   }
638   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
639   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
640 
641   /* free workspace */
642   if (!pcbddc->use_nnsp_true && !skip_lapack) {
643     ierr = PetscFree(work);CHKERRQ(ierr);
644 #if defined(PETSC_USE_COMPLEX)
645     ierr = PetscFree(rwork);CHKERRQ(ierr);
646 #endif
647     ierr = PetscFree(singular_vals);CHKERRQ(ierr);
648 #if defined(PETSC_MISSING_LAPACK_GESVD)
649     ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
650     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
651 #endif
652   }
653   for (k=0;k<nnsp_size;k++) {
654     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
655   }
656   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
657 
658   /* set quantities in pcbddc data structure */
659   /* n_vertices defines the number of subdomain corners in the primal space */
660   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
661   pcbddc->local_primal_size = total_counts;
662   pcbddc->n_vertices = n_vertices;
663   pcbddc->n_constraints = pcbddc->local_primal_size-pcbddc->n_vertices;
664 
665   /* Create constraint matrix */
666   /* The constraint matrix is used to compute the l2g map of primal dofs */
667   /* so we need to set it up properly either with or without change of basis */
668   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
669   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
670   ierr = MatSetSizes(pcbddc->ConstraintMatrix,pcbddc->local_primal_size,pcis->n,pcbddc->local_primal_size,pcis->n);CHKERRQ(ierr);
671   /* array to compute a local numbering of constraints : vertices first then constraints */
672   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
673   /* array to select the proper local node (of minimum index with respect to global ordering) when changing the basis */
674   /* 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 */
675   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_minloc);CHKERRQ(ierr);
676   /* auxiliary stuff for basis change */
677   ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&global_indices);CHKERRQ(ierr);
678   ierr = PetscMalloc(pcis->n_B*sizeof(PetscBool),&touched);CHKERRQ(ierr);
679   ierr = PetscMemzero(touched,pcis->n_B*sizeof(PetscBool));CHKERRQ(ierr);
680 
681   /* find primal_dofs: subdomain corners plus dofs selected as primal after change of basis */
682   total_primal_vertices=0;
683   for (i=0;i<pcbddc->local_primal_size;i++) {
684     size_of_constraint=temp_indices[i+1]-temp_indices[i];
685     if (size_of_constraint == 1) {
686       touched[temp_indices_to_constraint_B[temp_indices[i]]]=PETSC_TRUE;
687       aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]];
688       aux_primal_minloc[total_primal_vertices]=0;
689       total_primal_vertices++;
690     } else if (change_basis[i]) { /* Same procedure used in PCBDDCGetPrimalConstraintsLocalIdx */
691       PetscInt min_loc,min_index;
692       ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],global_indices);CHKERRQ(ierr);
693       /* find first untouched local node */
694       k = 0;
695       while (touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) k++;
696       min_index = global_indices[k];
697       min_loc = k;
698       /* search the minimum among global nodes already untouched on the cc */
699       for (k=1;k<size_of_constraint;k++) {
700         /* there can be more than one constraint on a single connected component */
701         if (min_index > global_indices[k] && !touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) {
702           min_index = global_indices[k];
703           min_loc = k;
704         }
705       }
706       touched[temp_indices_to_constraint_B[temp_indices[i]+min_loc]] = PETSC_TRUE;
707       aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]+min_loc];
708       aux_primal_minloc[total_primal_vertices]=min_loc;
709       total_primal_vertices++;
710     }
711   }
712   /* free workspace */
713   ierr = PetscFree(global_indices);CHKERRQ(ierr);
714   ierr = PetscFree(touched);CHKERRQ(ierr);
715   /* permute indices in order to have a sorted set of vertices */
716   ierr = PetscSortInt(total_primal_vertices,aux_primal_numbering);
717 
718   /* nonzero structure of constraint matrix */
719   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
720   for (i=0;i<total_primal_vertices;i++) nnz[i]=1;
721   j=total_primal_vertices;
722   for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
723     if (!change_basis[i]) {
724       nnz[j]=temp_indices[i+1]-temp_indices[i];
725       j++;
726     }
727   }
728   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
729   ierr = PetscFree(nnz);CHKERRQ(ierr);
730   /* set values in constraint matrix */
731   for (i=0;i<total_primal_vertices;i++) {
732     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,aux_primal_numbering[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
733   }
734   total_counts = total_primal_vertices;
735   for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
736     if (!change_basis[i]) {
737       size_of_constraint=temp_indices[i+1]-temp_indices[i];
738       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);
739       total_counts++;
740     }
741   }
742   /* assembling */
743   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
744   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
745   /*
746   ierr = MatView(pcbddc->ConstraintMatrix,(PetscViewer)0);CHKERRQ(ierr);
747   */
748   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
749   if (pcbddc->use_change_of_basis) {
750     PetscBool qr_needed = PETSC_FALSE;
751     /* change of basis acts on local interfaces -> dimension is n_B x n_B */
752     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
753     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
754     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
755     /* work arrays */
756     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
757     for (i=0;i<pcis->n_B;i++) nnz[i]=1;
758     /* nonzeros per row */
759     for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
760       if (change_basis[i]) {
761         qr_needed = PETSC_TRUE;
762         size_of_constraint = temp_indices[i+1]-temp_indices[i];
763         for (j=0;j<size_of_constraint;j++) nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
764       }
765     }
766     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
767     ierr = PetscFree(nnz);CHKERRQ(ierr);
768     /* Set initial identity in the matrix */
769     for (i=0;i<pcis->n_B;i++) {
770       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
771     }
772 
773     /* Now we loop on the constraints which need a change of basis */
774     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
775     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (see Sect 6.2.1) */
776     /* Change of basis matrix T computed via QR decomposition of constraints */
777     if (qr_needed) {
778       /* dual and primal dofs on a single cc */
779       PetscInt     dual_dofs,primal_dofs;
780       /* iterator on aux_primal_minloc (ordered as read from nearnullspace: vertices, edges and then constraints) */
781       PetscInt     primal_counter;
782       /* working stuff for GEQRF */
783       PetscScalar  *qr_basis,*qr_tau,*qr_work,lqr_work_t;
784       PetscBLASInt lqr_work;
785       /* working stuff for UNGQR */
786       PetscScalar  *gqr_work,lgqr_work_t;
787       PetscBLASInt lgqr_work;
788       /* working stuff for TRTRS */
789       PetscScalar  *trs_rhs;
790       PetscBLASInt Blas_NRHS;
791       /* pointers for values insertion into change of basis matrix */
792       PetscInt     *start_rows,*start_cols;
793       PetscScalar  *start_vals;
794       /* working stuff for values insertion */
795       PetscBool    *is_primal;
796 
797       /* space to store Q */
798       ierr = PetscMalloc((max_size_of_constraint)*(max_size_of_constraint)*sizeof(PetscScalar),&qr_basis);CHKERRQ(ierr);
799       /* first we issue queries for optimal work */
800       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_M);CHKERRQ(ierr);
801       ierr = PetscBLASIntCast(max_constraints,&Blas_N);CHKERRQ(ierr);
802       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
803       lqr_work = -1;
804       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Blas_M,&Blas_N,qr_basis,&Blas_LDA,qr_tau,&lqr_work_t,&lqr_work,&lierr));
805       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GEQRF Lapack routine %d",(int)lierr);
806       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lqr_work_t),&lqr_work);CHKERRQ(ierr);
807       ierr = PetscMalloc((PetscInt)PetscRealPart(lqr_work_t)*sizeof(*qr_work),&qr_work);CHKERRQ(ierr);
808       lgqr_work = -1;
809       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_M);CHKERRQ(ierr);
810       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_N);CHKERRQ(ierr);
811       ierr = PetscBLASIntCast(max_constraints,&Blas_K);CHKERRQ(ierr);
812       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
813       if (Blas_K>Blas_M) Blas_K=Blas_M; /* adjust just for computing optimal work */
814       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Blas_M,&Blas_N,&Blas_K,qr_basis,&Blas_LDA,qr_tau,&lgqr_work_t,&lgqr_work,&lierr));
815       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to UNGQR Lapack routine %d",(int)lierr);
816       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lgqr_work_t),&lgqr_work);CHKERRQ(ierr);
817       ierr = PetscMalloc((PetscInt)PetscRealPart(lgqr_work_t)*sizeof(*gqr_work),&gqr_work);CHKERRQ(ierr);
818       /* array to store scaling factors for reflectors */
819       ierr = PetscMalloc(max_constraints*sizeof(*qr_tau),&qr_tau);CHKERRQ(ierr);
820       /* array to store rhs and solution of triangular solver */
821       ierr = PetscMalloc(max_constraints*max_constraints*sizeof(*trs_rhs),&trs_rhs);CHKERRQ(ierr);
822       /* array to store whether a node is primal or not */
823       ierr = PetscMalloc(pcis->n_B*sizeof(*is_primal),&is_primal);CHKERRQ(ierr);
824       ierr = PetscMemzero(is_primal,pcis->n_B*sizeof(*is_primal));CHKERRQ(ierr);
825       for (i=0;i<total_primal_vertices;i++) is_primal[local_to_B[aux_primal_numbering[i]]] = PETSC_TRUE;
826 
827       /* allocating workspace for check */
828       if (pcbddc->dbg_flag) {
829         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
830         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Checking change of basis computation for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
831         ierr = PetscMalloc(max_size_of_constraint*(max_constraints+max_size_of_constraint)*sizeof(*work),&work);CHKERRQ(ierr);
832       }
833 
834       /* loop on constraints and see whether or not they need a change of basis */
835       /* -> using implicit ordering contained in temp_indices data */
836       total_counts = pcbddc->n_vertices;
837       primal_counter = total_counts;
838       while (total_counts<pcbddc->local_primal_size) {
839         primal_dofs = 1;
840         if (change_basis[total_counts]) {
841           /* get all constraints with same support: if more then one constraint is present on the cc then surely indices are stored contiguosly */
842           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]]) {
843             primal_dofs++;
844           }
845           /* get constraint info */
846           size_of_constraint = temp_indices[total_counts+1]-temp_indices[total_counts];
847           dual_dofs = size_of_constraint-primal_dofs;
848 
849           /* copy quadrature constraints for change of basis check */
850           if (pcbddc->dbg_flag) {
851             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);
852             ierr = PetscMemcpy(work,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
853           }
854 
855           /* copy temporary constraints into larger work vector (in order to store all columns of Q) */
856           ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
857 
858           /* compute QR decomposition of constraints */
859           ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
860           ierr = PetscBLASIntCast(primal_dofs,&Blas_N);CHKERRQ(ierr);
861           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
862           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
863           PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Blas_M,&Blas_N,qr_basis,&Blas_LDA,qr_tau,qr_work,&lqr_work,&lierr));
864           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GEQRF Lapack routine %d",(int)lierr);
865           ierr = PetscFPTrapPop();CHKERRQ(ierr);
866 
867           /* explictly compute R^-T */
868           ierr = PetscMemzero(trs_rhs,primal_dofs*primal_dofs*sizeof(*trs_rhs));CHKERRQ(ierr);
869           for (j=0;j<primal_dofs;j++) trs_rhs[j*(primal_dofs+1)] = 1.0;
870           ierr = PetscBLASIntCast(primal_dofs,&Blas_N);CHKERRQ(ierr);
871           ierr = PetscBLASIntCast(primal_dofs,&Blas_NRHS);CHKERRQ(ierr);
872           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
873           ierr = PetscBLASIntCast(primal_dofs,&Blas_LDB);CHKERRQ(ierr);
874           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
875           PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&Blas_N,&Blas_NRHS,qr_basis,&Blas_LDA,trs_rhs,&Blas_LDB,&lierr));
876           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TRTRS Lapack routine %d",(int)lierr);
877           ierr = PetscFPTrapPop();CHKERRQ(ierr);
878 
879           /* explcitly compute all columns of Q (Q = [Q1 | Q2] ) overwriting QR factorization in qr_basis */
880           ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
881           ierr = PetscBLASIntCast(size_of_constraint,&Blas_N);CHKERRQ(ierr);
882           ierr = PetscBLASIntCast(primal_dofs,&Blas_K);CHKERRQ(ierr);
883           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
884           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
885           PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Blas_M,&Blas_N,&Blas_K,qr_basis,&Blas_LDA,qr_tau,gqr_work,&lgqr_work,&lierr));
886           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in UNGQR Lapack routine %d",(int)lierr);
887           ierr = PetscFPTrapPop();CHKERRQ(ierr);
888 
889           /* first primal_dofs columns of Q need to be re-scaled in order to be unitary w.r.t constraints
890              i.e. C_{pxn}*Q_{nxn} should be equal to [I_pxp | 0_pxd] (see check below)
891              where n=size_of_constraint, p=primal_dofs, d=dual_dofs (n=p+d), I and 0 identity and null matrix resp. */
892           ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
893           ierr = PetscBLASIntCast(primal_dofs,&Blas_N);CHKERRQ(ierr);
894           ierr = PetscBLASIntCast(primal_dofs,&Blas_K);CHKERRQ(ierr);
895           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
896           ierr = PetscBLASIntCast(primal_dofs,&Blas_LDB);CHKERRQ(ierr);
897           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDC);CHKERRQ(ierr);
898           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
899           PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&Blas_M,&Blas_N,&Blas_K,&one,qr_basis,&Blas_LDA,trs_rhs,&Blas_LDB,&zero,&temp_quadrature_constraint[temp_indices[total_counts]],&Blas_LDC));
900           ierr = PetscFPTrapPop();CHKERRQ(ierr);
901           ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
902 
903           /* insert values in change of basis matrix respecting global ordering of new primal dofs */
904           start_rows = &temp_indices_to_constraint_B[temp_indices[total_counts]];
905           /* insert cols for primal dofs */
906           for (j=0;j<primal_dofs;j++) {
907             start_vals = &qr_basis[j*size_of_constraint];
908             start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+aux_primal_minloc[primal_counter+j]];
909             ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr);
910           }
911           /* insert cols for dual dofs */
912           for (j=0,k=0;j<dual_dofs;k++) {
913             if (!is_primal[temp_indices_to_constraint_B[temp_indices[total_counts]+k]]) {
914               start_vals = &qr_basis[(primal_dofs+j)*size_of_constraint];
915               start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+k];
916               ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr);
917               j++;
918             }
919           }
920 
921           /* check change of basis */
922           if (pcbddc->dbg_flag) {
923             PetscInt   ii,jj;
924             PetscBool valid_qr=PETSC_TRUE;
925             ierr = PetscBLASIntCast(primal_dofs,&Blas_M);CHKERRQ(ierr);
926             ierr = PetscBLASIntCast(size_of_constraint,&Blas_N);CHKERRQ(ierr);
927             ierr = PetscBLASIntCast(size_of_constraint,&Blas_K);CHKERRQ(ierr);
928             ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
929             ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDB);CHKERRQ(ierr);
930             ierr = PetscBLASIntCast(primal_dofs,&Blas_LDC);CHKERRQ(ierr);
931             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
932             PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&Blas_M,&Blas_N,&Blas_K,&one,work,&Blas_LDA,qr_basis,&Blas_LDB,&zero,&work[size_of_constraint*primal_dofs],&Blas_LDC));
933             ierr = PetscFPTrapPop();CHKERRQ(ierr);
934             for (jj=0;jj<size_of_constraint;jj++) {
935               for (ii=0;ii<primal_dofs;ii++) {
936                 if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) valid_qr = PETSC_FALSE;
937                 if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) valid_qr = PETSC_FALSE;
938               }
939             }
940             if (!valid_qr) {
941               ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> wrong change of basis!\n",PetscGlobalRank);CHKERRQ(ierr);
942               for (jj=0;jj<size_of_constraint;jj++) {
943                 for (ii=0;ii<primal_dofs;ii++) {
944                   if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) {
945                     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]));
946                   }
947                   if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) {
948                     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]));
949                   }
950                 }
951               }
952             } else {
953               ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> right change of basis!\n",PetscGlobalRank);CHKERRQ(ierr);
954             }
955           }
956           /* increment primal counter */
957           primal_counter += primal_dofs;
958         } else {
959           if (pcbddc->dbg_flag) {
960             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);
961           }
962         }
963         /* increment constraint counter total_counts */
964         total_counts += primal_dofs;
965       }
966       if (pcbddc->dbg_flag) {
967         ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
968         ierr = PetscFree(work);CHKERRQ(ierr);
969       }
970       /* free workspace */
971       ierr = PetscFree(trs_rhs);CHKERRQ(ierr);
972       ierr = PetscFree(qr_tau);CHKERRQ(ierr);
973       ierr = PetscFree(qr_work);CHKERRQ(ierr);
974       ierr = PetscFree(gqr_work);CHKERRQ(ierr);
975       ierr = PetscFree(is_primal);CHKERRQ(ierr);
976       ierr = PetscFree(qr_basis);CHKERRQ(ierr);
977     }
978     /* assembling */
979     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
980     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
981     /*
982     ierr = MatView(pcbddc->ChangeOfBasisMatrix,(PetscViewer)0);CHKERRQ(ierr);
983     */
984   }
985   /* free workspace no longer needed */
986   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
987   ierr = PetscFree(aux_primal_minloc);CHKERRQ(ierr);
988   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
989   ierr = PetscFree(change_basis);CHKERRQ(ierr);
990   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
991   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
992   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
993   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
994   PetscFunctionReturn(0);
995 }
996 
997 #undef __FUNCT__
998 #define __FUNCT__ "PCBDDCAnalyzeInterface"
999 PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
1000 {
1001   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
1002   PC_IS       *pcis = (PC_IS*)pc->data;
1003   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
1004   PetscInt    bs,ierr,i,vertex_size;
1005   PetscViewer viewer=pcbddc->dbg_viewer;
1006 
1007   PetscFunctionBegin;
1008   /* Init local Graph struct */
1009   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
1010 
1011   /* Check validity of the csr graph passed in by the user */
1012   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
1013     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
1014   }
1015   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
1016   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
1017     Mat mat_adj;
1018     const PetscInt *xadj,*adjncy;
1019     PetscBool flg_row=PETSC_TRUE;
1020 
1021     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
1022     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1023     if (!flg_row) {
1024       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
1025     }
1026     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
1027     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1028     if (!flg_row) {
1029       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
1030     }
1031     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
1032   }
1033 
1034   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
1035   vertex_size = 1;
1036   if (!pcbddc->n_ISForDofs) {
1037     IS *custom_ISForDofs;
1038 
1039     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
1040     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
1041     for (i=0;i<bs;i++) {
1042       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
1043     }
1044     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
1045     /* remove my references to IS objects */
1046     for (i=0;i<bs;i++) {
1047       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
1048     }
1049     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
1050   } else { /* mat block size as vertex size (used for elasticity) */
1051     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
1052   }
1053 
1054   /* Setup of Graph */
1055   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
1056 
1057   /* Graph's connected components analysis */
1058   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
1059 
1060   /* print some info to stdout */
1061   if (pcbddc->dbg_flag) {
1062     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
1063   }
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 #undef __FUNCT__
1068 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
1069 PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
1070 {
1071   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1072   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
1073   PetscErrorCode ierr;
1074 
1075   PetscFunctionBegin;
1076   n = 0;
1077   vertices = 0;
1078   if (pcbddc->ConstraintMatrix) {
1079     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
1080     for (i=0;i<local_primal_size;i++) {
1081       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1082       if (size_of_constraint == 1) n++;
1083       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1084     }
1085     ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
1086     n = 0;
1087     for (i=0;i<local_primal_size;i++) {
1088       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1089       if (size_of_constraint == 1) {
1090         vertices[n++]=row_cmat_indices[0];
1091       }
1092       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1093     }
1094   }
1095   *n_vertices = n;
1096   *vertices_idx = vertices;
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #undef __FUNCT__
1101 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
1102 PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
1103 {
1104   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1105   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
1106   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
1107   PetscBool      *touched;
1108   PetscErrorCode ierr;
1109 
1110   PetscFunctionBegin;
1111   n = 0;
1112   constraints_index = 0;
1113   if (pcbddc->ConstraintMatrix) {
1114     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
1115     max_size_of_constraint = 0;
1116     for (i=0;i<local_primal_size;i++) {
1117       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1118       if (size_of_constraint > 1) {
1119         n++;
1120       }
1121       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1122       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1123     }
1124     ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1125     ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1126     ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1127     ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1128     n = 0;
1129     for (i=0;i<local_primal_size;i++) {
1130       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1131       if (size_of_constraint > 1) {
1132         ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
1133         min_index = row_cmat_global_indices[0];
1134         min_loc = 0;
1135         for (j=1;j<size_of_constraint;j++) {
1136           /* there can be more than one constraint on a single connected component */
1137           if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1138             min_index = row_cmat_global_indices[j];
1139             min_loc = j;
1140           }
1141         }
1142         touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1143         constraints_index[n++] = row_cmat_indices[min_loc];
1144       }
1145       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1146     }
1147   }
1148   ierr = PetscFree(touched);CHKERRQ(ierr);
1149   ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1150   *n_constraints = n;
1151   *constraints_idx = constraints_index;
1152   PetscFunctionReturn(0);
1153 }
1154 
1155 /* the next two functions has been adapted from pcis.c */
1156 #undef __FUNCT__
1157 #define __FUNCT__ "PCBDDCApplySchur"
1158 PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1159 {
1160   PetscErrorCode ierr;
1161   PC_IS          *pcis = (PC_IS*)(pc->data);
1162 
1163   PetscFunctionBegin;
1164   if (!vec2_B) { vec2_B = v; }
1165   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1166   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1167   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1168   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1169   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 #undef __FUNCT__
1174 #define __FUNCT__ "PCBDDCApplySchurTranspose"
1175 PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1176 {
1177   PetscErrorCode ierr;
1178   PC_IS          *pcis = (PC_IS*)(pc->data);
1179 
1180   PetscFunctionBegin;
1181   if (!vec2_B) { vec2_B = v; }
1182   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1183   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1184   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1185   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1186   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1187   PetscFunctionReturn(0);
1188 }
1189 
1190 #undef __FUNCT__
1191 #define __FUNCT__ "PCBDDCSubsetNumbering"
1192 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[])
1193 {
1194   Vec            local_vec,global_vec;
1195   IS             seqis,paris;
1196   VecScatter     scatter_ctx;
1197   PetscScalar    *array;
1198   PetscInt       *temp_global_dofs;
1199   PetscScalar    globalsum;
1200   PetscInt       i,j,s;
1201   PetscInt       nlocals,first_index,old_index,max_local;
1202   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1203   PetscMPIInt    *dof_sizes,*dof_displs;
1204   PetscBool      first_found;
1205   PetscErrorCode ierr;
1206 
1207   PetscFunctionBegin;
1208   /* mpi buffers */
1209   MPI_Comm_size(comm,&size_prec_comm);
1210   MPI_Comm_rank(comm,&rank_prec_comm);
1211   j = ( !rank_prec_comm ? size_prec_comm : 0);
1212   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1213   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1214   /* get maximum size of subset */
1215   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1216   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1217   max_local = 0;
1218   if (n_local_dofs) {
1219     max_local = temp_global_dofs[0];
1220     for (i=1;i<n_local_dofs;i++) {
1221       if (max_local < temp_global_dofs[i] ) {
1222         max_local = temp_global_dofs[i];
1223       }
1224     }
1225   }
1226   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1227   max_global++;
1228   max_local = 0;
1229   if (n_local_dofs) {
1230     max_local = local_dofs[0];
1231     for (i=1;i<n_local_dofs;i++) {
1232       if (max_local < local_dofs[i] ) {
1233         max_local = local_dofs[i];
1234       }
1235     }
1236   }
1237   max_local++;
1238   /* allocate workspace */
1239   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1240   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1241   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1242   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1243   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1244   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1245   /* create scatter */
1246   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1247   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1248   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1249   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1250   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1251   /* init array */
1252   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1253   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1254   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1255   if (local_dofs_mult) {
1256     for (i=0;i<n_local_dofs;i++) {
1257       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1258     }
1259   } else {
1260     for (i=0;i<n_local_dofs;i++) {
1261       array[local_dofs[i]]=1.0;
1262     }
1263   }
1264   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1265   /* scatter into global vec and get total number of global dofs */
1266   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1267   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1268   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
1269   *n_global_subset = (PetscInt)PetscRealPart(globalsum);
1270   /* Fill global_vec with cumulative function for global numbering */
1271   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1272   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1273   nlocals = 0;
1274   first_index = -1;
1275   first_found = PETSC_FALSE;
1276   for (i=0;i<s;i++) {
1277     if (!first_found && PetscRealPart(array[i]) > 0.0) {
1278       first_found = PETSC_TRUE;
1279       first_index = i;
1280     }
1281     nlocals += (PetscInt)PetscRealPart(array[i]);
1282   }
1283   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1284   if (!rank_prec_comm) {
1285     dof_displs[0]=0;
1286     for (i=1;i<size_prec_comm;i++) {
1287       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1288     }
1289   }
1290   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1291   if (first_found) {
1292     array[first_index] += (PetscScalar)nlocals;
1293     old_index = first_index;
1294     for (i=first_index+1;i<s;i++) {
1295       if (PetscRealPart(array[i]) > 0.0) {
1296         array[i] += array[old_index];
1297         old_index = i;
1298       }
1299     }
1300   }
1301   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1302   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1303   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1304   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1305   /* get global ordering of local dofs */
1306   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1307   if (local_dofs_mult) {
1308     for (i=0;i<n_local_dofs;i++) {
1309       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i];
1310     }
1311   } else {
1312     for (i=0;i<n_local_dofs;i++) {
1313       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1;
1314     }
1315   }
1316   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1317   /* free workspace */
1318   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1319   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1320   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1321   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1322   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1323   /* return pointer to global ordering of local dofs */
1324   *global_numbering_subset = temp_global_dofs;
1325   PetscFunctionReturn(0);
1326 }
1327