xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision 62a6ff1d89f223ad7991b4be5dee3de016785c7b)
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->coarse_vec);CHKERRQ(ierr);
49   ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr);
50   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
51   ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr);
52   ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr);
53   ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr);
54   ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr);
55   ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr);
56   ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr);
57   ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr);
58   ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr);
59   ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr);
60   ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr);
61   ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr);
62   ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr);
63   ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr);
64   ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr);
65   ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr);
66   ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr);
67   ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr);
68   ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr);
69   ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr);
70   ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "PCBDDCCreateWorkVectors"
76 PetscErrorCode PCBDDCCreateWorkVectors(PC pc)
77 {
78   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
79   PC_IS          *pcis = (PC_IS*)pc->data;
80   VecType        impVecType;
81   PetscInt       n_vertices,n_constraints,local_primal_size,n_R;
82   PetscErrorCode ierr;
83 
84   PetscFunctionBegin;
85   ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,NULL);CHKERRQ(ierr);
86   ierr = PCBDDCGetPrimalConstraintsLocalIdx(pc,&n_constraints,NULL);CHKERRQ(ierr);
87   local_primal_size = n_constraints+n_vertices;
88   n_R = pcis->n-n_vertices;
89   /* local work vectors */
90   ierr = VecGetType(pcis->vec1_N,&impVecType);CHKERRQ(ierr);
91   ierr = VecDuplicate(pcis->vec1_D,&pcbddc->vec4_D);CHKERRQ(ierr);
92   ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_R);CHKERRQ(ierr);
93   ierr = VecSetSizes(pcbddc->vec1_R,PETSC_DECIDE,n_R);CHKERRQ(ierr);
94   ierr = VecSetType(pcbddc->vec1_R,impVecType);CHKERRQ(ierr);
95   ierr = VecDuplicate(pcbddc->vec1_R,&pcbddc->vec2_R);CHKERRQ(ierr);
96   ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_P);CHKERRQ(ierr);
97   ierr = VecSetSizes(pcbddc->vec1_P,PETSC_DECIDE,local_primal_size);CHKERRQ(ierr);
98   ierr = VecSetType(pcbddc->vec1_P,impVecType);CHKERRQ(ierr);
99   if (n_constraints) {
100     ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_C);CHKERRQ(ierr);
101     ierr = VecSetSizes(pcbddc->vec1_C,PETSC_DECIDE,n_constraints);CHKERRQ(ierr);
102     ierr = VecSetType(pcbddc->vec1_C,impVecType);CHKERRQ(ierr);
103   }
104   PetscFunctionReturn(0);
105 }
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "PCBDDCSetUpLocalMatrices"
109 PetscErrorCode PCBDDCSetUpLocalMatrices(PC pc)
110 {
111   PC_IS*            pcis = (PC_IS*)(pc->data);
112   PC_BDDC*          pcbddc = (PC_BDDC*)pc->data;
113   Mat_IS*           matis = (Mat_IS*)pc->pmat->data;
114   /* manage repeated solves */
115   MatReuse          reuse;
116   MatStructure      matstruct;
117   PetscErrorCode    ierr;
118 
119   PetscFunctionBegin;
120   /* get mat flags */
121   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
122   reuse = MAT_INITIAL_MATRIX;
123   if (pc->setupcalled) {
124     /* when matstruct is SAME_PRECONDITIONER, we shouldn't be here */
125     if (matstruct == SAME_PRECONDITIONER) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen");
126     if (matstruct == SAME_NONZERO_PATTERN) {
127       reuse = MAT_REUSE_MATRIX;
128     } else {
129       reuse = MAT_INITIAL_MATRIX;
130     }
131   }
132   if (reuse == MAT_INITIAL_MATRIX) {
133     ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
134     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
135     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
136     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
137     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
138   }
139 
140   /* transform local matrices if needed */
141   if (pcbddc->use_change_of_basis) {
142     Mat         change_mat_all;
143     PetscScalar *row_cmat_values;
144     PetscInt    *row_cmat_indices;
145     PetscInt    *nnz,*is_indices,*temp_indices;
146     PetscInt    i,j,k,n_D,n_B;
147 
148     /* Get Non-overlapping dimensions */
149     n_B = pcis->n_B;
150     n_D = pcis->n-n_B;
151 
152     /* compute nonzero structure of change of basis on all local nodes */
153     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
154     ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
155     for (i=0;i<n_D;i++) nnz[is_indices[i]] = 1;
156     ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
157     ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
158     k=1;
159     for (i=0;i<n_B;i++) {
160       ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr);
161       nnz[is_indices[i]]=j;
162       if (k < j) k = j;
163       ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr);
164     }
165     ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
166     /* assemble change of basis matrix on the whole set of local dofs */
167     ierr = PetscMalloc(k*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
168     ierr = MatCreate(PETSC_COMM_SELF,&change_mat_all);CHKERRQ(ierr);
169     ierr = MatSetSizes(change_mat_all,pcis->n,pcis->n,pcis->n,pcis->n);CHKERRQ(ierr);
170     ierr = MatSetType(change_mat_all,MATSEQAIJ);CHKERRQ(ierr);
171     ierr = MatSeqAIJSetPreallocation(change_mat_all,0,nnz);CHKERRQ(ierr);
172     ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
173     for (i=0;i<n_D;i++) {
174       ierr = MatSetValue(change_mat_all,is_indices[i],is_indices[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
175     }
176     ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
177     ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
178     for (i=0;i<n_B;i++) {
179       ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr);
180       for (k=0; k<j; k++) temp_indices[k]=is_indices[row_cmat_indices[k]];
181       ierr = MatSetValues(change_mat_all,1,&is_indices[i],j,temp_indices,row_cmat_values,INSERT_VALUES);CHKERRQ(ierr);
182       ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr);
183     }
184     ierr = MatAssemblyBegin(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
185     ierr = MatAssemblyEnd(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
186     /* TODO: HOW TO WORK WITH BAIJ? PtAP not provided */
187     ierr = MatGetBlockSize(matis->A,&i);CHKERRQ(ierr);
188     if (i==1) {
189       ierr = MatPtAP(matis->A,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
190     } else {
191       Mat work_mat;
192       ierr = MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&work_mat);CHKERRQ(ierr);
193       ierr = MatPtAP(work_mat,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
194       ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
195     }
196     ierr = MatDestroy(&change_mat_all);CHKERRQ(ierr);
197     ierr = PetscFree(nnz);CHKERRQ(ierr);
198     ierr = PetscFree(temp_indices);CHKERRQ(ierr);
199   } else {
200     /* without change of basis, the local matrix is unchanged */
201     if (!pcbddc->local_mat) {
202       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
203       pcbddc->local_mat = matis->A;
204     }
205   }
206 
207   /* get submatrices */
208   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
209   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
210   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
211   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "PCBDDCSetUpLocalScatters"
217 PetscErrorCode PCBDDCSetUpLocalScatters(PC pc,IS* is_R_local_n)
218 {
219   PC_IS*         pcis = (PC_IS*)(pc->data);
220   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
221   IS             is_R_local,is_aux1,is_aux2;
222   PetscInt       *vertices,*aux_array1,*aux_array2,*is_indices,*idx_R_local;
223   PetscInt       n_vertices,n_constraints,i,j,n_R,n_D,n_B;
224   PetscBool      *array_bool;
225   PetscErrorCode ierr;
226 
227   PetscFunctionBegin;
228   /* Set Non-overlapping dimensions */
229   n_B = pcis->n_B; n_D = pcis->n - n_B;
230   /* get vertex indices from constraint matrix */
231   ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,&vertices);CHKERRQ(ierr);
232   /* Set number of constraints */
233   n_constraints = pcbddc->local_primal_size-n_vertices;
234   /* Dohrmann's notation: dofs splitted in R (Remaining: all dofs but the vertices) and V (Vertices) */
235   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&array_bool);CHKERRQ(ierr);
236   for (i=0;i<pcis->n;i++) array_bool[i] = PETSC_TRUE;
237   for (i=0;i<n_vertices;i++) array_bool[vertices[i]] = PETSC_FALSE;
238   ierr = PetscMalloc((pcis->n-n_vertices)*sizeof(PetscInt),&idx_R_local);CHKERRQ(ierr);
239   for (i=0, n_R=0; i<pcis->n; i++) {
240     if (array_bool[i]) {
241       idx_R_local[n_R] = i;
242       n_R++;
243     }
244   }
245   ierr = PetscFree(vertices);CHKERRQ(ierr);
246   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_R,idx_R_local,PETSC_OWN_POINTER,&is_R_local);CHKERRQ(ierr);
247 
248   /* print some info if requested */
249   if (pcbddc->dbg_flag) {
250     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
251     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
252     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d local dimensions\n",PetscGlobalRank);CHKERRQ(ierr);
253     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"local_size = %d, dirichlet_size = %d, boundary_size = %d\n",pcis->n,n_D,n_B);CHKERRQ(ierr);
254     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"r_size = %d, v_size = %d, constraints = %d, local_primal_size = %d\n",n_R,n_vertices,n_constraints,pcbddc->local_primal_size);CHKERRQ(ierr);
255     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"pcbddc->n_vertices = %d, pcbddc->n_constraints = %d\n",pcbddc->n_vertices,pcbddc->n_constraints);CHKERRQ(ierr);
256     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
257   }
258 
259   /* VecScatters pcbddc->R_to_B and (optionally) pcbddc->R_to_D */
260   ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr);
261   ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array2);CHKERRQ(ierr);
262   ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
263   for (i=0; i<n_D; i++) array_bool[is_indices[i]] = PETSC_FALSE;
264   ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
265   for (i=0, j=0; i<n_R; i++) {
266     if (array_bool[idx_R_local[i]]) {
267       aux_array1[j] = i;
268       j++;
269     }
270   }
271   ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr);
272   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
273   for (i=0, j=0; i<n_B; i++) {
274     if (array_bool[is_indices[i]]) {
275       aux_array2[j] = i; j++;
276     }
277   }
278   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
279   ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array2,PETSC_OWN_POINTER,&is_aux2);CHKERRQ(ierr);
280   ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_B,is_aux2,&pcbddc->R_to_B);CHKERRQ(ierr);
281   ierr = ISDestroy(&is_aux1);CHKERRQ(ierr);
282   ierr = ISDestroy(&is_aux2);CHKERRQ(ierr);
283 
284   if (pcbddc->inexact_prec_type || pcbddc->dbg_flag ) {
285     ierr = PetscMalloc(n_D*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr);
286     for (i=0, j=0; i<n_R; i++) {
287       if (!array_bool[idx_R_local[i]]) {
288         aux_array1[j] = i;
289         j++;
290       }
291     }
292     ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr);
293     ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_D,(IS)0,&pcbddc->R_to_D);CHKERRQ(ierr);
294     ierr = ISDestroy(&is_aux1);CHKERRQ(ierr);
295   }
296   ierr = PetscFree(array_bool);CHKERRQ(ierr);
297   *is_R_local_n = is_R_local;
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNCT__
302 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
303 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool use)
304 {
305   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
306 
307   PetscFunctionBegin;
308   pcbddc->use_exact_dirichlet=use;
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "PCBDDCSetUpLocalSolvers"
314 PetscErrorCode PCBDDCSetUpLocalSolvers(PC pc, IS is_I_local, IS is_R_local)
315 {
316   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
317   PC_IS          *pcis = (PC_IS*)pc->data;
318   PC             pc_temp;
319   Mat            A_RR;
320   Vec            vec1,vec2,vec3;
321   MatStructure   matstruct;
322   PetscScalar    m_one = -1.0;
323   PetscReal      value;
324   PetscInt       n_D,n_R,use_exact,use_exact_reduced;
325   PetscErrorCode ierr;
326 
327   PetscFunctionBegin;
328   /* Creating PC contexts for local Dirichlet and Neumann problems */
329   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
330 
331   /* DIRICHLET PROBLEM */
332   /* Matrix for Dirichlet problem is pcis->A_II */
333   ierr = ISGetSize(is_I_local,&n_D);CHKERRQ(ierr);
334   if (!pcbddc->ksp_D) { /* create object if not yet build */
335     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_D);CHKERRQ(ierr);
336     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
337     /* default */
338     ierr = KSPSetType(pcbddc->ksp_D,KSPPREONLY);CHKERRQ(ierr);
339     ierr = KSPSetOptionsPrefix(pcbddc->ksp_D,"dirichlet_");CHKERRQ(ierr);
340     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
341     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
342     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
343   }
344   ierr = KSPSetOperators(pcbddc->ksp_D,pcis->A_II,pcis->A_II,matstruct);CHKERRQ(ierr);
345   /* Allow user's customization */
346   ierr = KSPSetFromOptions(pcbddc->ksp_D);CHKERRQ(ierr);
347   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
348   if (!n_D) {
349     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
350     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
351   }
352   /* Set Up KSP for Dirichlet problem of BDDC */
353   ierr = KSPSetUp(pcbddc->ksp_D);CHKERRQ(ierr);
354   /* set ksp_D into pcis data */
355   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
356   ierr = PetscObjectReference((PetscObject)pcbddc->ksp_D);CHKERRQ(ierr);
357   pcis->ksp_D = pcbddc->ksp_D;
358 
359   /* NEUMANN PROBLEM */
360   /* Matrix for Neumann problem is A_RR -> we need to create it */
361   ierr = ISGetSize(is_R_local,&n_R);CHKERRQ(ierr);
362   ierr = MatGetSubMatrix(pcbddc->local_mat,is_R_local,is_R_local,MAT_INITIAL_MATRIX,&A_RR);CHKERRQ(ierr);
363   if (!pcbddc->ksp_R) { /* create object if not yet build */
364     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_R);CHKERRQ(ierr);
365     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_R,(PetscObject)pc,1);CHKERRQ(ierr);
366     /* default */
367     ierr = KSPSetType(pcbddc->ksp_R,KSPPREONLY);CHKERRQ(ierr);
368     ierr = KSPSetOptionsPrefix(pcbddc->ksp_R,"neumann_");CHKERRQ(ierr);
369     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
370     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
371     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
372   }
373   ierr = KSPSetOperators(pcbddc->ksp_R,A_RR,A_RR,matstruct);CHKERRQ(ierr);
374   /* Allow user's customization */
375   ierr = KSPSetFromOptions(pcbddc->ksp_R);CHKERRQ(ierr);
376   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
377   if (!n_R) {
378     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
379     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
380   }
381   /* Set Up KSP for Neumann problem of BDDC */
382   ierr = KSPSetUp(pcbddc->ksp_R);CHKERRQ(ierr);
383 
384   /* check Dirichlet and Neumann solvers and adapt them if a nullspace correction is needed */
385 
386   /* Dirichlet */
387   ierr = MatGetVecs(pcis->A_II,&vec1,&vec2);CHKERRQ(ierr);
388   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
389   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
390   ierr = MatMult(pcis->A_II,vec1,vec2);CHKERRQ(ierr);
391   ierr = KSPSolve(pcbddc->ksp_D,vec2,vec3);CHKERRQ(ierr);
392   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
393   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
394   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
395   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
396   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
397   /* need to be adapted? */
398   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
399   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
400   ierr = PCBDDCSetUseExactDirichlet(pc,(PetscBool)use_exact_reduced);CHKERRQ(ierr);
401   /* print info */
402   if (pcbddc->dbg_flag) {
403     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
404     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
405     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Checking solution of Dirichlet and Neumann problems\n");CHKERRQ(ierr);
406     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
407     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
408   }
409   if (n_D && pcbddc->NullSpace && !use_exact_reduced && !pcbddc->inexact_prec_type) {
410     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_I_local);CHKERRQ(ierr);
411   }
412 
413   /* Neumann */
414   ierr = MatGetVecs(A_RR,&vec1,&vec2);CHKERRQ(ierr);
415   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
416   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
417   ierr = MatMult(A_RR,vec1,vec2);CHKERRQ(ierr);
418   ierr = KSPSolve(pcbddc->ksp_R,vec2,vec3);CHKERRQ(ierr);
419   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
420   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
421   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
422   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
423   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
424   /* need to be adapted? */
425   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
426   if (PetscAbsReal(value) > 1.e-4) use_exact = 0;
427   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
428   /* print info */
429   if (pcbddc->dbg_flag) {
430     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for  Neumann  solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
431     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
432   }
433   if (n_R && pcbddc->NullSpace && !use_exact_reduced) { /* is it the right logic? */
434     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_R_local);CHKERRQ(ierr);
435   }
436 
437   /* free Neumann problem's matrix */
438   ierr = MatDestroy(&A_RR);CHKERRQ(ierr);
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "PCBDDCSolveSaddlePoint"
444 static PetscErrorCode  PCBDDCSolveSaddlePoint(PC pc)
445 {
446   PetscErrorCode ierr;
447   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
448 
449   PetscFunctionBegin;
450   ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr);
451   if (pcbddc->local_auxmat1) {
452     ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr);
453     ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr);
454   }
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner"
460 PetscErrorCode  PCBDDCApplyInterfacePreconditioner(PC pc)
461 {
462   PetscErrorCode ierr;
463   PC_BDDC*        pcbddc = (PC_BDDC*)(pc->data);
464   PC_IS*            pcis = (PC_IS*)  (pc->data);
465   const PetscScalar zero = 0.0;
466 
467   PetscFunctionBegin;
468   /* Application of PHI^T (or PSI^T)  */
469   if (pcbddc->coarse_psi_B) {
470     ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
471     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
472   } else {
473     ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
474     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
475   }
476   /* Scatter data of coarse_rhs */
477   if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); }
478   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
479 
480   /* Local solution on R nodes */
481   ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr);
482   ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
483   ierr = VecScatterEnd  (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
484   if (pcbddc->inexact_prec_type) {
485     ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
486     ierr = VecScatterEnd  (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
487   }
488   ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr);
489   ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
490   ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
491   ierr = VecScatterEnd  (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
492   if (pcbddc->inexact_prec_type) {
493     ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
494     ierr = VecScatterEnd  (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
495   }
496 
497   /* Coarse solution */
498   ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
499   if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */
500     ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr);
501   }
502   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
503   ierr = PCBDDCScatterCoarseDataEnd  (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
504 
505   /* Sum contributions from two levels */
506   ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
507   if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "PCBDDCScatterCoarseDataBegin"
513 PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
514 {
515   PetscErrorCode ierr;
516   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
517 
518   PetscFunctionBegin;
519   switch (pcbddc->coarse_communications_type) {
520     case SCATTERS_BDDC:
521       ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
522       break;
523     case GATHERS_BDDC:
524       break;
525   }
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "PCBDDCScatterCoarseDataEnd"
531 PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
532 {
533   PetscErrorCode ierr;
534   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
535   PetscScalar*   array_to;
536   PetscScalar*   array_from;
537   MPI_Comm       comm;
538   PetscInt       i;
539 
540   PetscFunctionBegin;
541   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
542   switch (pcbddc->coarse_communications_type) {
543     case SCATTERS_BDDC:
544       ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
545       break;
546     case GATHERS_BDDC:
547       if (vec_from) {
548         ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr);
549       }
550       if (vec_to) {
551         ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr);
552       }
553       switch(pcbddc->coarse_problem_type){
554         case SEQUENTIAL_BDDC:
555           if (smode == SCATTER_FORWARD) {
556             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);
557             if (vec_to) {
558               if (imode == ADD_VALUES) {
559                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
560                   array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
561                 }
562               } else {
563                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
564                   array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
565                 }
566               }
567             }
568           } else {
569             if (vec_from) {
570               if (imode == ADD_VALUES) {
571                 MPI_Comm vec_from_comm;
572                 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr);
573                 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);
574               }
575               for (i=0;i<pcbddc->replicated_primal_size;i++) {
576                 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]];
577               }
578             }
579             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);
580           }
581           break;
582         case REPLICATED_BDDC:
583           if (smode == SCATTER_FORWARD) {
584             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);
585             if (imode == ADD_VALUES) {
586               for (i=0;i<pcbddc->replicated_primal_size;i++) {
587                 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
588               }
589             } else {
590               for (i=0;i<pcbddc->replicated_primal_size;i++) {
591                 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
592               }
593             }
594           } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */
595             if (imode == ADD_VALUES) {
596               for (i=0;i<pcbddc->local_primal_size;i++) {
597                 array_to[i]+=array_from[pcbddc->local_primal_indices[i]];
598               }
599             } else {
600               for (i=0;i<pcbddc->local_primal_size;i++) {
601                 array_to[i]=array_from[pcbddc->local_primal_indices[i]];
602               }
603             }
604           }
605           break;
606         case MULTILEVEL_BDDC:
607           break;
608         case PARALLEL_BDDC:
609           break;
610       }
611       if (vec_from) {
612         ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr);
613       }
614       if (vec_to) {
615         ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr);
616       }
617       break;
618   }
619   PetscFunctionReturn(0);
620 }
621 
622 /* uncomment for testing purposes */
623 /* #define PETSC_MISSING_LAPACK_GESVD 1 */
624 #undef __FUNCT__
625 #define __FUNCT__ "PCBDDCConstraintsSetUp"
626 PetscErrorCode PCBDDCConstraintsSetUp(PC pc)
627 {
628   PetscErrorCode    ierr;
629   PC_IS*            pcis = (PC_IS*)(pc->data);
630   PC_BDDC*          pcbddc = (PC_BDDC*)pc->data;
631   Mat_IS*           matis = (Mat_IS*)pc->pmat->data;
632   /* constraint and (optionally) change of basis matrix implemented as SeqAIJ */
633   MatType           impMatType=MATSEQAIJ;
634   /* one and zero */
635   PetscScalar       one=1.0,zero=0.0;
636   /* space to store constraints and their local indices */
637   PetscScalar       *temp_quadrature_constraint;
638   PetscInt          *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B;
639   /* iterators */
640   PetscInt          i,j,k,total_counts,temp_start_ptr;
641   /* stuff to store connected components stored in pcbddc->mat_graph */
642   IS                ISForVertices,*ISForFaces,*ISForEdges,*used_IS;
643   PetscInt          n_ISForFaces,n_ISForEdges;
644   PetscBool         get_faces,get_edges,get_vertices;
645   /* near null space stuff */
646   MatNullSpace      nearnullsp;
647   const Vec         *nearnullvecs;
648   Vec               *localnearnullsp;
649   PetscBool         nnsp_has_cnst;
650   PetscInt          nnsp_size;
651   PetscScalar       *array;
652   /* BLAS integers */
653   PetscBLASInt      lwork,lierr;
654   PetscBLASInt      Blas_N,Blas_M,Blas_K,Blas_one=1;
655   PetscBLASInt      Blas_LDA,Blas_LDB,Blas_LDC;
656   /* LAPACK working arrays for SVD or POD */
657   PetscBool         skip_lapack;
658   PetscScalar       *work;
659   PetscReal         *singular_vals;
660 #if defined(PETSC_USE_COMPLEX)
661   PetscReal         *rwork;
662 #endif
663 #if defined(PETSC_MISSING_LAPACK_GESVD)
664   PetscBLASInt      Blas_one_2=1;
665   PetscScalar       *temp_basis,*correlation_mat;
666 #else
667   PetscBLASInt      dummy_int_1=1,dummy_int_2=1;
668   PetscScalar       dummy_scalar_1=0.0,dummy_scalar_2=0.0;
669 #endif
670   /* change of basis */
671   PetscInt          *aux_primal_numbering,*aux_primal_minloc,*global_indices;
672   PetscBool         boolforchange,*change_basis,*touched;
673   /* auxiliary stuff */
674   PetscInt          *nnz,*is_indices,*local_to_B;
675   /* some quantities */
676   PetscInt          n_vertices,total_primal_vertices;
677   PetscInt          size_of_constraint,max_size_of_constraint,max_constraints,temp_constraints;
678 
679 
680   PetscFunctionBegin;
681   /* Get index sets for faces, edges and vertices from graph */
682   get_faces = PETSC_TRUE;
683   get_edges = PETSC_TRUE;
684   get_vertices = PETSC_TRUE;
685   if (pcbddc->vertices_flag) {
686     get_faces = PETSC_FALSE;
687     get_edges = PETSC_FALSE;
688   }
689   if (pcbddc->constraints_flag) {
690     get_vertices = PETSC_FALSE;
691   }
692   if (pcbddc->faces_flag) {
693     get_edges = PETSC_FALSE;
694   }
695   if (pcbddc->edges_flag) {
696     get_faces = PETSC_FALSE;
697   }
698   /* default */
699   if (!get_faces && !get_edges && !get_vertices) {
700     get_vertices = PETSC_TRUE;
701   }
702   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
703   /* print some info */
704   if (pcbddc->dbg_flag) {
705     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
706     i = 0;
707     if (ISForVertices) {
708       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
709     }
710     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
711     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
712     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
713     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
714   }
715   /* check if near null space is attached to global mat */
716   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
717   if (nearnullsp) {
718     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
719   } else { /* if near null space is not provided BDDC uses constants by default */
720     nnsp_size = 0;
721     nnsp_has_cnst = PETSC_TRUE;
722   }
723   /* get max number of constraints on a single cc */
724   max_constraints = nnsp_size;
725   if (nnsp_has_cnst) max_constraints++;
726 
727   /*
728        Evaluate maximum storage size needed by the procedure
729        - temp_indices will contain start index of each constraint stored as follows
730        - temp_indices_to_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts
731        - 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
732        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
733                                                                                                                                                          */
734   total_counts = n_ISForFaces+n_ISForEdges;
735   total_counts *= max_constraints;
736   n_vertices = 0;
737   if (ISForVertices) {
738     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
739   }
740   total_counts += n_vertices;
741   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
742   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
743   total_counts = 0;
744   max_size_of_constraint = 0;
745   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
746     if (i<n_ISForEdges) {
747       used_IS = &ISForEdges[i];
748     } else {
749       used_IS = &ISForFaces[i-n_ISForEdges];
750     }
751     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
752     total_counts += j;
753     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
754   }
755   total_counts *= max_constraints;
756   total_counts += n_vertices;
757   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
758   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
759   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
760   /* local to boundary numbering */
761   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
762   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
763   for (i=0;i<pcis->n;i++) local_to_B[i]=-1;
764   for (i=0;i<pcis->n_B;i++) local_to_B[is_indices[i]]=i;
765   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
766   /* get local part of global near null space vectors */
767   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
768   for (k=0;k<nnsp_size;k++) {
769     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
770     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
771     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
772   }
773 
774   /* whether or not to skip lapack calls */
775   skip_lapack = PETSC_TRUE;
776   if (n_ISForFaces+n_ISForEdges) skip_lapack = PETSC_FALSE;
777 
778   /* First we issue queries to allocate optimal workspace for LAPACKgesvd (or LAPACKsyev if SVD is missing) */
779   if (!pcbddc->use_nnsp_true && !skip_lapack) {
780     PetscScalar temp_work;
781 #if defined(PETSC_MISSING_LAPACK_GESVD)
782     /* Proper Orthogonal Decomposition (POD) using the snapshot method */
783     ierr = PetscMalloc(max_constraints*max_constraints*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
784     ierr = PetscMalloc(max_constraints*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
785     ierr = PetscMalloc(max_size_of_constraint*max_constraints*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
786 #if defined(PETSC_USE_COMPLEX)
787     ierr = PetscMalloc(3*max_constraints*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
788 #endif
789     /* now we evaluate the optimal workspace using query with lwork=-1 */
790     ierr = PetscBLASIntCast(max_constraints,&Blas_N);CHKERRQ(ierr);
791     ierr = PetscBLASIntCast(max_constraints,&Blas_LDA);CHKERRQ(ierr);
792     lwork = -1;
793     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
794 #if !defined(PETSC_USE_COMPLEX)
795     PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,&temp_work,&lwork,&lierr));
796 #else
797     PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,&temp_work,&lwork,rwork,&lierr));
798 #endif
799     ierr = PetscFPTrapPop();CHKERRQ(ierr);
800     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEV Lapack routine %d",(int)lierr);
801 #else /* on missing GESVD */
802     /* SVD */
803     PetscInt max_n,min_n;
804     max_n = max_size_of_constraint;
805     min_n = max_constraints;
806     if (max_size_of_constraint < max_constraints) {
807       min_n = max_size_of_constraint;
808       max_n = max_constraints;
809     }
810     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
811 #if defined(PETSC_USE_COMPLEX)
812     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
813 #endif
814     /* now we evaluate the optimal workspace using query with lwork=-1 */
815     lwork = -1;
816     ierr = PetscBLASIntCast(max_n,&Blas_M);CHKERRQ(ierr);
817     ierr = PetscBLASIntCast(min_n,&Blas_N);CHKERRQ(ierr);
818     ierr = PetscBLASIntCast(max_n,&Blas_LDA);CHKERRQ(ierr);
819     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
820 #if !defined(PETSC_USE_COMPLEX)
821     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Blas_M,&Blas_N,&temp_quadrature_constraint[0],&Blas_LDA,singular_vals,&dummy_scalar_1,&dummy_int_1,&dummy_scalar_2,&dummy_int_2,&temp_work,&lwork,&lierr));
822 #else
823     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Blas_M,&Blas_N,&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));
824 #endif
825     ierr = PetscFPTrapPop();CHKERRQ(ierr);
826     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GESVD Lapack routine %d",(int)lierr);
827 #endif /* on missing GESVD */
828     /* Allocate optimal workspace */
829     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
830     ierr = PetscMalloc((PetscInt)lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
831   }
832   /* Now we can loop on constraining sets */
833   total_counts = 0;
834   temp_indices[0] = 0;
835   /* vertices */
836   if (ISForVertices) {
837     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
838     if (nnsp_has_cnst) { /* consider all vertices */
839       for (i=0;i<n_vertices;i++) {
840         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
841         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
842         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
843         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
844         change_basis[total_counts]=PETSC_FALSE;
845         total_counts++;
846       }
847     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
848       PetscBool used_vertex;
849       for (i=0;i<n_vertices;i++) {
850         used_vertex = PETSC_FALSE;
851         k = 0;
852         while (!used_vertex && k<nnsp_size) {
853           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
854           if (PetscAbsScalar(array[is_indices[i]])>0.0) {
855             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
856             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
857             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
858             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
859             change_basis[total_counts]=PETSC_FALSE;
860             total_counts++;
861             used_vertex = PETSC_TRUE;
862           }
863           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
864           k++;
865         }
866       }
867     }
868     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
869     n_vertices = total_counts;
870   }
871 
872   /* edges and faces */
873   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
874     if (i<n_ISForEdges) {
875       used_IS = &ISForEdges[i];
876       boolforchange = pcbddc->use_change_of_basis; /* change or not the basis on the edge */
877     } else {
878       used_IS = &ISForFaces[i-n_ISForEdges];
879       boolforchange = (PetscBool)(pcbddc->use_change_of_basis && pcbddc->use_change_on_faces); /* change or not the basis on the face */
880     }
881     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
882     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
883     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
884     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
885     /* change of basis should not be performed on local periodic nodes */
886     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) boolforchange = PETSC_FALSE;
887     if (nnsp_has_cnst) {
888       PetscScalar quad_value;
889       temp_constraints++;
890       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
891       for (j=0;j<size_of_constraint;j++) {
892         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
893         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
894         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
895       }
896       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
897       change_basis[total_counts]=boolforchange;
898       total_counts++;
899     }
900     for (k=0;k<nnsp_size;k++) {
901       PetscReal real_value;
902       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
903       for (j=0;j<size_of_constraint;j++) {
904         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
905         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
906         temp_quadrature_constraint[temp_indices[total_counts]+j]=array[is_indices[j]];
907       }
908       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array);CHKERRQ(ierr);
909       /* check if array is null on the connected component */
910       ierr = PetscBLASIntCast(size_of_constraint,&Blas_N);CHKERRQ(ierr);
911       PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Blas_N,&temp_quadrature_constraint[temp_indices[total_counts]],&Blas_one));
912       if (real_value > 0.0) { /* keep indices and values */
913         temp_constraints++;
914         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
915         change_basis[total_counts]=boolforchange;
916         total_counts++;
917       }
918     }
919     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
920     /* perform SVD on the constraints if use_nnsp_true has not be requested by the user */
921     if (!pcbddc->use_nnsp_true) {
922       PetscReal tol = 1.0e-8; /* tolerance for retaining eigenmodes */
923 
924 #if defined(PETSC_MISSING_LAPACK_GESVD)
925       /* SVD: Y = U*S*V^H                -> U (eigenvectors of Y*Y^H) = Y*V*(S)^\dag
926          POD: Y^H*Y = V*D*V^H, D = S^H*S -> U = Y*V*D^(-1/2)
927          -> When PETSC_USE_COMPLEX and PETSC_MISSING_LAPACK_GESVD are defined
928             the constraints basis will differ (by a complex factor with absolute value equal to 1)
929             from that computed using LAPACKgesvd
930          -> This is due to a different computation of eigenvectors in LAPACKheev
931          -> The quality of the POD-computed basis will be the same */
932       ierr = PetscMemzero(correlation_mat,temp_constraints*temp_constraints*sizeof(PetscScalar));CHKERRQ(ierr);
933       /* Store upper triangular part of correlation matrix */
934       ierr = PetscBLASIntCast(size_of_constraint,&Blas_N);CHKERRQ(ierr);
935       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
936       for (j=0;j<temp_constraints;j++) {
937         for (k=0;k<j+1;k++) {
938           PetscStackCallBLAS("BLASdot",correlation_mat[j*temp_constraints+k]=BLASdot_(&Blas_N,&temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Blas_one,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Blas_one_2));
939         }
940       }
941       /* compute eigenvalues and eigenvectors of correlation matrix */
942       ierr = PetscBLASIntCast(temp_constraints,&Blas_N);CHKERRQ(ierr);
943       ierr = PetscBLASIntCast(temp_constraints,&Blas_LDA);CHKERRQ(ierr);
944 #if !defined(PETSC_USE_COMPLEX)
945       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,work,&lwork,&lierr));
946 #else
947       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&Blas_N,correlation_mat,&Blas_LDA,singular_vals,work,&lwork,rwork,&lierr));
948 #endif
949       ierr = PetscFPTrapPop();CHKERRQ(ierr);
950       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEV Lapack routine %d",(int)lierr);
951       /* retain eigenvalues greater than tol: note that LAPACKsyev gives eigs in ascending order */
952       j=0;
953       while (j < temp_constraints && singular_vals[j] < tol) j++;
954       total_counts=total_counts-j;
955       /* scale and copy POD basis into used quadrature memory */
956       ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
957       ierr = PetscBLASIntCast(temp_constraints,&Blas_N);CHKERRQ(ierr);
958       ierr = PetscBLASIntCast(temp_constraints,&Blas_K);CHKERRQ(ierr);
959       ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
960       ierr = PetscBLASIntCast(temp_constraints,&Blas_LDB);CHKERRQ(ierr);
961       ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDC);CHKERRQ(ierr);
962       if (j<temp_constraints) {
963         PetscInt ii;
964         for (k=j;k<temp_constraints;k++) singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
965         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
966         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));
967         ierr = PetscFPTrapPop();CHKERRQ(ierr);
968         for (k=0;k<temp_constraints-j;k++) {
969           for (ii=0;ii<size_of_constraint;ii++) {
970             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];
971           }
972         }
973       }
974 #else  /* on missing GESVD */
975       ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
976       ierr = PetscBLASIntCast(temp_constraints,&Blas_N);CHKERRQ(ierr);
977       ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
978       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
979 #if !defined(PETSC_USE_COMPLEX)
980       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Blas_M,&Blas_N,&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));
981 #else
982       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Blas_M,&Blas_N,&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));
983 #endif
984       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESVD Lapack routine %d",(int)lierr);
985       ierr = PetscFPTrapPop();CHKERRQ(ierr);
986       /* retain eigenvalues greater than tol: note that LAPACKgesvd gives eigs in descending order */
987       k = temp_constraints;
988       if (k > size_of_constraint) k = size_of_constraint;
989       j = 0;
990       while (j < k && singular_vals[k-j-1] < tol) j++;
991       total_counts = total_counts-temp_constraints+k-j;
992 #endif /* on missing GESVD */
993     }
994   }
995   /* free index sets of faces, edges and vertices */
996   for (i=0;i<n_ISForFaces;i++) {
997     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
998   }
999   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
1000   for (i=0;i<n_ISForEdges;i++) {
1001     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
1002   }
1003   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
1004   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
1005 
1006   /* free workspace */
1007   if (!pcbddc->use_nnsp_true && !skip_lapack) {
1008     ierr = PetscFree(work);CHKERRQ(ierr);
1009 #if defined(PETSC_USE_COMPLEX)
1010     ierr = PetscFree(rwork);CHKERRQ(ierr);
1011 #endif
1012     ierr = PetscFree(singular_vals);CHKERRQ(ierr);
1013 #if defined(PETSC_MISSING_LAPACK_GESVD)
1014     ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
1015     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
1016 #endif
1017   }
1018   for (k=0;k<nnsp_size;k++) {
1019     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
1020   }
1021   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
1022 
1023   /* set quantities in pcbddc data structure */
1024   /* n_vertices defines the number of subdomain corners in the primal space */
1025   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
1026   pcbddc->local_primal_size = total_counts;
1027   pcbddc->n_vertices = n_vertices;
1028   pcbddc->n_constraints = pcbddc->local_primal_size-pcbddc->n_vertices;
1029 
1030   /* Create constraint matrix */
1031   /* The constraint matrix is used to compute the l2g map of primal dofs */
1032   /* so we need to set it up properly either with or without change of basis */
1033   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
1034   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
1035   ierr = MatSetSizes(pcbddc->ConstraintMatrix,pcbddc->local_primal_size,pcis->n,pcbddc->local_primal_size,pcis->n);CHKERRQ(ierr);
1036   /* array to compute a local numbering of constraints : vertices first then constraints */
1037   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
1038   /* array to select the proper local node (of minimum index with respect to global ordering) when changing the basis */
1039   /* 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 */
1040   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&aux_primal_minloc);CHKERRQ(ierr);
1041   /* auxiliary stuff for basis change */
1042   ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&global_indices);CHKERRQ(ierr);
1043   ierr = PetscMalloc(pcis->n_B*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1044   ierr = PetscMemzero(touched,pcis->n_B*sizeof(PetscBool));CHKERRQ(ierr);
1045 
1046   /* find primal_dofs: subdomain corners plus dofs selected as primal after change of basis */
1047   total_primal_vertices=0;
1048   for (i=0;i<pcbddc->local_primal_size;i++) {
1049     size_of_constraint=temp_indices[i+1]-temp_indices[i];
1050     if (size_of_constraint == 1) {
1051       touched[temp_indices_to_constraint_B[temp_indices[i]]]=PETSC_TRUE;
1052       aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]];
1053       aux_primal_minloc[total_primal_vertices]=0;
1054       total_primal_vertices++;
1055     } else if (change_basis[i]) { /* Same procedure used in PCBDDCGetPrimalConstraintsLocalIdx */
1056       PetscInt min_loc,min_index;
1057       ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,&temp_indices_to_constraint[temp_indices[i]],global_indices);CHKERRQ(ierr);
1058       /* find first untouched local node */
1059       k = 0;
1060       while (touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) k++;
1061       min_index = global_indices[k];
1062       min_loc = k;
1063       /* search the minimum among global nodes already untouched on the cc */
1064       for (k=1;k<size_of_constraint;k++) {
1065         /* there can be more than one constraint on a single connected component */
1066         if (min_index > global_indices[k] && !touched[temp_indices_to_constraint_B[temp_indices[i]+k]]) {
1067           min_index = global_indices[k];
1068           min_loc = k;
1069         }
1070       }
1071       touched[temp_indices_to_constraint_B[temp_indices[i]+min_loc]] = PETSC_TRUE;
1072       aux_primal_numbering[total_primal_vertices]=temp_indices_to_constraint[temp_indices[i]+min_loc];
1073       aux_primal_minloc[total_primal_vertices]=min_loc;
1074       total_primal_vertices++;
1075     }
1076   }
1077   /* free workspace */
1078   ierr = PetscFree(global_indices);CHKERRQ(ierr);
1079   ierr = PetscFree(touched);CHKERRQ(ierr);
1080   /* permute indices in order to have a sorted set of vertices */
1081   ierr = PetscSortInt(total_primal_vertices,aux_primal_numbering);
1082 
1083   /* nonzero structure of constraint matrix */
1084   ierr = PetscMalloc(pcbddc->local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
1085   for (i=0;i<total_primal_vertices;i++) nnz[i]=1;
1086   j=total_primal_vertices;
1087   for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
1088     if (!change_basis[i]) {
1089       nnz[j]=temp_indices[i+1]-temp_indices[i];
1090       j++;
1091     }
1092   }
1093   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
1094   ierr = PetscFree(nnz);CHKERRQ(ierr);
1095   /* set values in constraint matrix */
1096   for (i=0;i<total_primal_vertices;i++) {
1097     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,aux_primal_numbering[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
1098   }
1099   total_counts = total_primal_vertices;
1100   for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
1101     if (!change_basis[i]) {
1102       size_of_constraint=temp_indices[i+1]-temp_indices[i];
1103       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);
1104       total_counts++;
1105     }
1106   }
1107   /* assembling */
1108   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1109   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1110   /*
1111   ierr = MatView(pcbddc->ConstraintMatrix,(PetscViewer)0);CHKERRQ(ierr);
1112   */
1113   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
1114   if (pcbddc->use_change_of_basis) {
1115     PetscBool qr_needed = PETSC_FALSE;
1116     /* change of basis acts on local interfaces -> dimension is n_B x n_B */
1117     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
1118     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
1119     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
1120     /* work arrays */
1121     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
1122     for (i=0;i<pcis->n_B;i++) nnz[i]=1;
1123     /* nonzeros per row */
1124     for (i=pcbddc->n_vertices;i<pcbddc->local_primal_size;i++) {
1125       if (change_basis[i]) {
1126         qr_needed = PETSC_TRUE;
1127         size_of_constraint = temp_indices[i+1]-temp_indices[i];
1128         for (j=0;j<size_of_constraint;j++) nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
1129       }
1130     }
1131     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
1132     ierr = PetscFree(nnz);CHKERRQ(ierr);
1133     /* Set initial identity in the matrix */
1134     for (i=0;i<pcis->n_B;i++) {
1135       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
1136     }
1137 
1138     /* Now we loop on the constraints which need a change of basis */
1139     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
1140     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (see Sect 6.2.1) */
1141     /* Change of basis matrix T computed via QR decomposition of constraints */
1142     if (qr_needed) {
1143       /* dual and primal dofs on a single cc */
1144       PetscInt     dual_dofs,primal_dofs;
1145       /* iterator on aux_primal_minloc (ordered as read from nearnullspace: vertices, edges and then constraints) */
1146       PetscInt     primal_counter;
1147       /* working stuff for GEQRF */
1148       PetscScalar  *qr_basis,*qr_tau,*qr_work,lqr_work_t;
1149       PetscBLASInt lqr_work;
1150       /* working stuff for UNGQR */
1151       PetscScalar  *gqr_work,lgqr_work_t;
1152       PetscBLASInt lgqr_work;
1153       /* working stuff for TRTRS */
1154       PetscScalar  *trs_rhs;
1155       PetscBLASInt Blas_NRHS;
1156       /* pointers for values insertion into change of basis matrix */
1157       PetscInt     *start_rows,*start_cols;
1158       PetscScalar  *start_vals;
1159       /* working stuff for values insertion */
1160       PetscBool    *is_primal;
1161 
1162       /* space to store Q */
1163       ierr = PetscMalloc((max_size_of_constraint)*(max_size_of_constraint)*sizeof(PetscScalar),&qr_basis);CHKERRQ(ierr);
1164       /* first we issue queries for optimal work */
1165       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_M);CHKERRQ(ierr);
1166       ierr = PetscBLASIntCast(max_constraints,&Blas_N);CHKERRQ(ierr);
1167       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1168       lqr_work = -1;
1169       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Blas_M,&Blas_N,qr_basis,&Blas_LDA,qr_tau,&lqr_work_t,&lqr_work,&lierr));
1170       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GEQRF Lapack routine %d",(int)lierr);
1171       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lqr_work_t),&lqr_work);CHKERRQ(ierr);
1172       ierr = PetscMalloc((PetscInt)PetscRealPart(lqr_work_t)*sizeof(*qr_work),&qr_work);CHKERRQ(ierr);
1173       lgqr_work = -1;
1174       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_M);CHKERRQ(ierr);
1175       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_N);CHKERRQ(ierr);
1176       ierr = PetscBLASIntCast(max_constraints,&Blas_K);CHKERRQ(ierr);
1177       ierr = PetscBLASIntCast(max_size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1178       if (Blas_K>Blas_M) Blas_K=Blas_M; /* adjust just for computing optimal work */
1179       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Blas_M,&Blas_N,&Blas_K,qr_basis,&Blas_LDA,qr_tau,&lgqr_work_t,&lgqr_work,&lierr));
1180       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to UNGQR Lapack routine %d",(int)lierr);
1181       ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lgqr_work_t),&lgqr_work);CHKERRQ(ierr);
1182       ierr = PetscMalloc((PetscInt)PetscRealPart(lgqr_work_t)*sizeof(*gqr_work),&gqr_work);CHKERRQ(ierr);
1183       /* array to store scaling factors for reflectors */
1184       ierr = PetscMalloc(max_constraints*sizeof(*qr_tau),&qr_tau);CHKERRQ(ierr);
1185       /* array to store rhs and solution of triangular solver */
1186       ierr = PetscMalloc(max_constraints*max_constraints*sizeof(*trs_rhs),&trs_rhs);CHKERRQ(ierr);
1187       /* array to store whether a node is primal or not */
1188       ierr = PetscMalloc(pcis->n_B*sizeof(*is_primal),&is_primal);CHKERRQ(ierr);
1189       ierr = PetscMemzero(is_primal,pcis->n_B*sizeof(*is_primal));CHKERRQ(ierr);
1190       for (i=0;i<total_primal_vertices;i++) is_primal[local_to_B[aux_primal_numbering[i]]] = PETSC_TRUE;
1191 
1192       /* allocating workspace for check */
1193       if (pcbddc->dbg_flag) {
1194         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
1195         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Checking change of basis computation for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
1196         ierr = PetscMalloc(max_size_of_constraint*(max_constraints+max_size_of_constraint)*sizeof(*work),&work);CHKERRQ(ierr);
1197       }
1198 
1199       /* loop on constraints and see whether or not they need a change of basis */
1200       /* -> using implicit ordering contained in temp_indices data */
1201       total_counts = pcbddc->n_vertices;
1202       primal_counter = total_counts;
1203       while (total_counts<pcbddc->local_primal_size) {
1204         primal_dofs = 1;
1205         if (change_basis[total_counts]) {
1206           /* get all constraints with same support: if more then one constraint is present on the cc then surely indices are stored contiguosly */
1207           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]]) {
1208             primal_dofs++;
1209           }
1210           /* get constraint info */
1211           size_of_constraint = temp_indices[total_counts+1]-temp_indices[total_counts];
1212           dual_dofs = size_of_constraint-primal_dofs;
1213 
1214           /* copy quadrature constraints for change of basis check */
1215           if (pcbddc->dbg_flag) {
1216             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);
1217             ierr = PetscMemcpy(work,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
1218           }
1219 
1220           /* copy temporary constraints into larger work vector (in order to store all columns of Q) */
1221           ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
1222 
1223           /* compute QR decomposition of constraints */
1224           ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
1225           ierr = PetscBLASIntCast(primal_dofs,&Blas_N);CHKERRQ(ierr);
1226           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1227           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1228           PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Blas_M,&Blas_N,qr_basis,&Blas_LDA,qr_tau,qr_work,&lqr_work,&lierr));
1229           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GEQRF Lapack routine %d",(int)lierr);
1230           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1231 
1232           /* explictly compute R^-T */
1233           ierr = PetscMemzero(trs_rhs,primal_dofs*primal_dofs*sizeof(*trs_rhs));CHKERRQ(ierr);
1234           for (j=0;j<primal_dofs;j++) trs_rhs[j*(primal_dofs+1)] = 1.0;
1235           ierr = PetscBLASIntCast(primal_dofs,&Blas_N);CHKERRQ(ierr);
1236           ierr = PetscBLASIntCast(primal_dofs,&Blas_NRHS);CHKERRQ(ierr);
1237           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1238           ierr = PetscBLASIntCast(primal_dofs,&Blas_LDB);CHKERRQ(ierr);
1239           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1240           PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U","T","N",&Blas_N,&Blas_NRHS,qr_basis,&Blas_LDA,trs_rhs,&Blas_LDB,&lierr));
1241           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TRTRS Lapack routine %d",(int)lierr);
1242           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1243 
1244           /* explcitly compute all columns of Q (Q = [Q1 | Q2] ) overwriting QR factorization in qr_basis */
1245           ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
1246           ierr = PetscBLASIntCast(size_of_constraint,&Blas_N);CHKERRQ(ierr);
1247           ierr = PetscBLASIntCast(primal_dofs,&Blas_K);CHKERRQ(ierr);
1248           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1249           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1250           PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Blas_M,&Blas_N,&Blas_K,qr_basis,&Blas_LDA,qr_tau,gqr_work,&lgqr_work,&lierr));
1251           if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in UNGQR Lapack routine %d",(int)lierr);
1252           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1253 
1254           /* first primal_dofs columns of Q need to be re-scaled in order to be unitary w.r.t constraints
1255              i.e. C_{pxn}*Q_{nxn} should be equal to [I_pxp | 0_pxd] (see check below)
1256              where n=size_of_constraint, p=primal_dofs, d=dual_dofs (n=p+d), I and 0 identity and null matrix resp. */
1257           ierr = PetscBLASIntCast(size_of_constraint,&Blas_M);CHKERRQ(ierr);
1258           ierr = PetscBLASIntCast(primal_dofs,&Blas_N);CHKERRQ(ierr);
1259           ierr = PetscBLASIntCast(primal_dofs,&Blas_K);CHKERRQ(ierr);
1260           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1261           ierr = PetscBLASIntCast(primal_dofs,&Blas_LDB);CHKERRQ(ierr);
1262           ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDC);CHKERRQ(ierr);
1263           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1264           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));
1265           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1266           ierr = PetscMemcpy(qr_basis,&temp_quadrature_constraint[temp_indices[total_counts]],size_of_constraint*primal_dofs*sizeof(PetscScalar));CHKERRQ(ierr);
1267 
1268           /* insert values in change of basis matrix respecting global ordering of new primal dofs */
1269           start_rows = &temp_indices_to_constraint_B[temp_indices[total_counts]];
1270           /* insert cols for primal dofs */
1271           for (j=0;j<primal_dofs;j++) {
1272             start_vals = &qr_basis[j*size_of_constraint];
1273             start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+aux_primal_minloc[primal_counter+j]];
1274             ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr);
1275           }
1276           /* insert cols for dual dofs */
1277           for (j=0,k=0;j<dual_dofs;k++) {
1278             if (!is_primal[temp_indices_to_constraint_B[temp_indices[total_counts]+k]]) {
1279               start_vals = &qr_basis[(primal_dofs+j)*size_of_constraint];
1280               start_cols = &temp_indices_to_constraint_B[temp_indices[total_counts]+k];
1281               ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,start_rows,1,start_cols,start_vals,INSERT_VALUES);CHKERRQ(ierr);
1282               j++;
1283             }
1284           }
1285 
1286           /* check change of basis */
1287           if (pcbddc->dbg_flag) {
1288             PetscInt   ii,jj;
1289             PetscBool valid_qr=PETSC_TRUE;
1290             ierr = PetscBLASIntCast(primal_dofs,&Blas_M);CHKERRQ(ierr);
1291             ierr = PetscBLASIntCast(size_of_constraint,&Blas_N);CHKERRQ(ierr);
1292             ierr = PetscBLASIntCast(size_of_constraint,&Blas_K);CHKERRQ(ierr);
1293             ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDA);CHKERRQ(ierr);
1294             ierr = PetscBLASIntCast(size_of_constraint,&Blas_LDB);CHKERRQ(ierr);
1295             ierr = PetscBLASIntCast(primal_dofs,&Blas_LDC);CHKERRQ(ierr);
1296             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1297             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));
1298             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1299             for (jj=0;jj<size_of_constraint;jj++) {
1300               for (ii=0;ii<primal_dofs;ii++) {
1301                 if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) valid_qr = PETSC_FALSE;
1302                 if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) valid_qr = PETSC_FALSE;
1303               }
1304             }
1305             if (!valid_qr) {
1306               ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> wrong change of basis!\n",PetscGlobalRank);CHKERRQ(ierr);
1307               for (jj=0;jj<size_of_constraint;jj++) {
1308                 for (ii=0;ii<primal_dofs;ii++) {
1309                   if (ii != jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]) > 1.e-12) {
1310                     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]));
1311                   }
1312                   if (ii == jj && PetscAbsScalar(work[size_of_constraint*primal_dofs+jj*primal_dofs+ii]-1.0) > 1.e-12) {
1313                     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]));
1314                   }
1315                 }
1316               }
1317             } else {
1318               ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"\t-> right change of basis!\n",PetscGlobalRank);CHKERRQ(ierr);
1319             }
1320           }
1321           /* increment primal counter */
1322           primal_counter += primal_dofs;
1323         } else {
1324           if (pcbddc->dbg_flag) {
1325             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);
1326           }
1327         }
1328         /* increment constraint counter total_counts */
1329         total_counts += primal_dofs;
1330       }
1331       if (pcbddc->dbg_flag) {
1332         ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
1333         ierr = PetscFree(work);CHKERRQ(ierr);
1334       }
1335       /* free workspace */
1336       ierr = PetscFree(trs_rhs);CHKERRQ(ierr);
1337       ierr = PetscFree(qr_tau);CHKERRQ(ierr);
1338       ierr = PetscFree(qr_work);CHKERRQ(ierr);
1339       ierr = PetscFree(gqr_work);CHKERRQ(ierr);
1340       ierr = PetscFree(is_primal);CHKERRQ(ierr);
1341       ierr = PetscFree(qr_basis);CHKERRQ(ierr);
1342     }
1343     /* assembling */
1344     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1345     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1346     /*
1347     ierr = MatView(pcbddc->ChangeOfBasisMatrix,(PetscViewer)0);CHKERRQ(ierr);
1348     */
1349   }
1350   /* free workspace */
1351   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
1352   ierr = PetscFree(aux_primal_minloc);CHKERRQ(ierr);
1353   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
1354   ierr = PetscFree(change_basis);CHKERRQ(ierr);
1355   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
1356   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
1357   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
1358   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
1359   PetscFunctionReturn(0);
1360 }
1361 
1362 #undef __FUNCT__
1363 #define __FUNCT__ "PCBDDCAnalyzeInterface"
1364 PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
1365 {
1366   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
1367   PC_IS       *pcis = (PC_IS*)pc->data;
1368   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
1369   PetscInt    bs,ierr,i,vertex_size;
1370   PetscViewer viewer=pcbddc->dbg_viewer;
1371 
1372   PetscFunctionBegin;
1373   /* Init local Graph struct */
1374   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
1375 
1376   /* Check validity of the csr graph passed in by the user */
1377   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
1378     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
1379   }
1380   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
1381   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
1382     Mat mat_adj;
1383     const PetscInt *xadj,*adjncy;
1384     PetscBool flg_row=PETSC_TRUE;
1385 
1386     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
1387     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1388     if (!flg_row) {
1389       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
1390     }
1391     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
1392     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1393     if (!flg_row) {
1394       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
1395     }
1396     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
1397   }
1398 
1399   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
1400   vertex_size = 1;
1401   if (!pcbddc->n_ISForDofs) {
1402     IS *custom_ISForDofs;
1403 
1404     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
1405     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
1406     for (i=0;i<bs;i++) {
1407       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
1408     }
1409     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
1410     /* remove my references to IS objects */
1411     for (i=0;i<bs;i++) {
1412       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
1413     }
1414     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
1415   } else { /* mat block size as vertex size (used for elasticity) */
1416     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
1417   }
1418 
1419   /* Setup of Graph */
1420   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
1421 
1422   /* Graph's connected components analysis */
1423   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
1424 
1425   /* print some info to stdout */
1426   if (pcbddc->dbg_flag) {
1427     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
1428   }
1429   PetscFunctionReturn(0);
1430 }
1431 
1432 #undef __FUNCT__
1433 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
1434 PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
1435 {
1436   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1437   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
1438   PetscErrorCode ierr;
1439 
1440   PetscFunctionBegin;
1441   n = 0;
1442   vertices = 0;
1443   if (pcbddc->ConstraintMatrix) {
1444     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
1445     for (i=0;i<local_primal_size;i++) {
1446       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1447       if (size_of_constraint == 1) n++;
1448       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1449     }
1450     if (vertices_idx) {
1451       ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
1452       n = 0;
1453       for (i=0;i<local_primal_size;i++) {
1454         ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1455         if (size_of_constraint == 1) {
1456           vertices[n++]=row_cmat_indices[0];
1457         }
1458         ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1459       }
1460     }
1461   }
1462   *n_vertices = n;
1463   if (vertices_idx) *vertices_idx = vertices;
1464   PetscFunctionReturn(0);
1465 }
1466 
1467 #undef __FUNCT__
1468 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
1469 PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
1470 {
1471   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1472   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
1473   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
1474   PetscBool      *touched;
1475   PetscErrorCode ierr;
1476 
1477   PetscFunctionBegin;
1478   n = 0;
1479   constraints_index = 0;
1480   if (pcbddc->ConstraintMatrix) {
1481     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
1482     max_size_of_constraint = 0;
1483     for (i=0;i<local_primal_size;i++) {
1484       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1485       if (size_of_constraint > 1) {
1486         n++;
1487       }
1488       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1489       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1490     }
1491     if (constraints_idx) {
1492       ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1493       ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1494       ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1495       ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1496       n = 0;
1497       for (i=0;i<local_primal_size;i++) {
1498         ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1499         if (size_of_constraint > 1) {
1500           ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
1501           /* find first untouched local node */
1502           j = 0;
1503           while(touched[row_cmat_indices[j]]) j++;
1504           min_index = row_cmat_global_indices[j];
1505           min_loc = j;
1506           /* search the minimum among nodes not yet touched on the connected component
1507              since there can be more than one constraint on a single cc */
1508           for (j=1;j<size_of_constraint;j++) {
1509             if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1510               min_index = row_cmat_global_indices[j];
1511               min_loc = j;
1512             }
1513           }
1514           touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1515           constraints_index[n++] = row_cmat_indices[min_loc];
1516         }
1517         ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1518       }
1519       ierr = PetscFree(touched);CHKERRQ(ierr);
1520       ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1521     }
1522   }
1523   *n_constraints = n;
1524   if (constraints_idx) *constraints_idx = constraints_index;
1525   PetscFunctionReturn(0);
1526 }
1527 
1528 /* the next two functions has been adapted from pcis.c */
1529 #undef __FUNCT__
1530 #define __FUNCT__ "PCBDDCApplySchur"
1531 PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1532 {
1533   PetscErrorCode ierr;
1534   PC_IS          *pcis = (PC_IS*)(pc->data);
1535 
1536   PetscFunctionBegin;
1537   if (!vec2_B) { vec2_B = v; }
1538   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1539   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1540   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1541   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1542   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1543   PetscFunctionReturn(0);
1544 }
1545 
1546 #undef __FUNCT__
1547 #define __FUNCT__ "PCBDDCApplySchurTranspose"
1548 PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1549 {
1550   PetscErrorCode ierr;
1551   PC_IS          *pcis = (PC_IS*)(pc->data);
1552 
1553   PetscFunctionBegin;
1554   if (!vec2_B) { vec2_B = v; }
1555   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1556   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1557   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1558   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1559   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1560   PetscFunctionReturn(0);
1561 }
1562 
1563 #undef __FUNCT__
1564 #define __FUNCT__ "PCBDDCSubsetNumbering"
1565 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[])
1566 {
1567   Vec            local_vec,global_vec;
1568   IS             seqis,paris;
1569   VecScatter     scatter_ctx;
1570   PetscScalar    *array;
1571   PetscInt       *temp_global_dofs;
1572   PetscScalar    globalsum;
1573   PetscInt       i,j,s;
1574   PetscInt       nlocals,first_index,old_index,max_local;
1575   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1576   PetscMPIInt    *dof_sizes,*dof_displs;
1577   PetscBool      first_found;
1578   PetscErrorCode ierr;
1579 
1580   PetscFunctionBegin;
1581   /* mpi buffers */
1582   MPI_Comm_size(comm,&size_prec_comm);
1583   MPI_Comm_rank(comm,&rank_prec_comm);
1584   j = ( !rank_prec_comm ? size_prec_comm : 0);
1585   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1586   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1587   /* get maximum size of subset */
1588   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1589   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1590   max_local = 0;
1591   if (n_local_dofs) {
1592     max_local = temp_global_dofs[0];
1593     for (i=1;i<n_local_dofs;i++) {
1594       if (max_local < temp_global_dofs[i] ) {
1595         max_local = temp_global_dofs[i];
1596       }
1597     }
1598   }
1599   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1600   max_global++;
1601   max_local = 0;
1602   if (n_local_dofs) {
1603     max_local = local_dofs[0];
1604     for (i=1;i<n_local_dofs;i++) {
1605       if (max_local < local_dofs[i] ) {
1606         max_local = local_dofs[i];
1607       }
1608     }
1609   }
1610   max_local++;
1611   /* allocate workspace */
1612   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1613   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1614   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1615   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1616   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1617   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1618   /* create scatter */
1619   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1620   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1621   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1622   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1623   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1624   /* init array */
1625   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1626   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1627   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1628   if (local_dofs_mult) {
1629     for (i=0;i<n_local_dofs;i++) {
1630       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1631     }
1632   } else {
1633     for (i=0;i<n_local_dofs;i++) {
1634       array[local_dofs[i]]=1.0;
1635     }
1636   }
1637   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1638   /* scatter into global vec and get total number of global dofs */
1639   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1640   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1641   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
1642   *n_global_subset = (PetscInt)PetscRealPart(globalsum);
1643   /* Fill global_vec with cumulative function for global numbering */
1644   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1645   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1646   nlocals = 0;
1647   first_index = -1;
1648   first_found = PETSC_FALSE;
1649   for (i=0;i<s;i++) {
1650     if (!first_found && PetscRealPart(array[i]) > 0.0) {
1651       first_found = PETSC_TRUE;
1652       first_index = i;
1653     }
1654     nlocals += (PetscInt)PetscRealPart(array[i]);
1655   }
1656   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1657   if (!rank_prec_comm) {
1658     dof_displs[0]=0;
1659     for (i=1;i<size_prec_comm;i++) {
1660       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1661     }
1662   }
1663   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1664   if (first_found) {
1665     array[first_index] += (PetscScalar)nlocals;
1666     old_index = first_index;
1667     for (i=first_index+1;i<s;i++) {
1668       if (PetscRealPart(array[i]) > 0.0) {
1669         array[i] += array[old_index];
1670         old_index = i;
1671       }
1672     }
1673   }
1674   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1675   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1676   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1677   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1678   /* get global ordering of local dofs */
1679   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1680   if (local_dofs_mult) {
1681     for (i=0;i<n_local_dofs;i++) {
1682       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i];
1683     }
1684   } else {
1685     for (i=0;i<n_local_dofs;i++) {
1686       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1;
1687     }
1688   }
1689   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1690   /* free workspace */
1691   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1692   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1693   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1694   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1695   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1696   /* return pointer to global ordering of local dofs */
1697   *global_numbering_subset = temp_global_dofs;
1698   PetscFunctionReturn(0);
1699 }
1700 
1701 #undef __FUNCT__
1702 #define __FUNCT__ "PCBDDCOrthonormalizeVecs"
1703 PetscErrorCode PCBDDCOrthonormalizeVecs(PetscInt n, Vec vecs[])
1704 {
1705   PetscInt       i,j;
1706   PetscScalar    *alphas;
1707   PetscErrorCode ierr;
1708 
1709   PetscFunctionBegin;
1710   /* this implements stabilized Gram-Schmidt */
1711   ierr = PetscMalloc(n*sizeof(PetscScalar),&alphas);CHKERRQ(ierr);
1712   for (i=0;i<n;i++) {
1713     ierr = VecNormalize(vecs[i],NULL);CHKERRQ(ierr);
1714     if (i<n) { ierr = VecMDot(vecs[i],n-i-1,&vecs[i+1],&alphas[i+1]);CHKERRQ(ierr); }
1715     for (j=i+1;j<n;j++) { ierr = VecAXPY(vecs[j],PetscConj(-alphas[j]),vecs[i]);CHKERRQ(ierr); }
1716   }
1717   ierr = PetscFree(alphas);CHKERRQ(ierr);
1718   PetscFunctionReturn(0);
1719 }
1720 
1721 
1722