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