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