xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision aa0d41d463ae17ff6de1da0e4fb594266144f895)
1 #include "bddc.h"
2 #include "bddcprivate.h"
3 #include <petscblaslapack.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "PCBDDCResetCustomization"
7 PetscErrorCode PCBDDCResetCustomization(PC pc)
8 {
9   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
10   PetscInt       i;
11   PetscErrorCode ierr;
12 
13   PetscFunctionBegin;
14   ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
15   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
16   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
17   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
18   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
19   for (i=0;i<pcbddc->n_ISForDofs;i++) {
20     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
21   }
22   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "PCBDDCResetTopography"
28 PetscErrorCode PCBDDCResetTopography(PC pc)
29 {
30   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
31   PetscErrorCode ierr;
32 
33   PetscFunctionBegin;
34   ierr = MatDestroy(&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
35   ierr = MatDestroy(&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
36   ierr = PCBDDCGraphReset(pcbddc->mat_graph);CHKERRQ(ierr);
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "PCBDDCResetSolvers"
42 PetscErrorCode PCBDDCResetSolvers(PC pc)
43 {
44   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
49   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
50   ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr);
51   ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr);
52   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
53   ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr);
54   ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr);
55   ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr);
56   ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr);
57   ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr);
58   ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr);
59   ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr);
60   ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr);
61   ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr);
62   ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr);
63   ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr);
64   ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr);
65   ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr);
66   ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr);
67   ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr);
68   ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr);
69   ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr);
70   ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr);
71   ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr);
72   ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "PCBDDCSetUpLocalMatrices"
78 PetscErrorCode PCBDDCSetUpLocalMatrices(PC pc)
79 {
80   PC_IS*            pcis = (PC_IS*)(pc->data);
81   PC_BDDC*          pcbddc = (PC_BDDC*)pc->data;
82   Mat_IS*           matis = (Mat_IS*)pc->pmat->data;
83   /* manage repeated solves */
84   MatReuse          reuse;
85   MatStructure      matstruct;
86   PetscErrorCode    ierr;
87 
88   PetscFunctionBegin;
89   /* get mat flags */
90   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
91   reuse = MAT_INITIAL_MATRIX;
92   if (pc->setupcalled) {
93     /* when matstruct is SAME_PRECONDITIONER, we shouldn't be here */
94     if (matstruct == SAME_PRECONDITIONER) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen");
95     if (matstruct == SAME_NONZERO_PATTERN) {
96       reuse = MAT_REUSE_MATRIX;
97     } else {
98       reuse = MAT_INITIAL_MATRIX;
99     }
100   }
101   if (reuse == MAT_INITIAL_MATRIX) {
102     ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
103     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
104     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
105     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
106     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
107   }
108 
109   /* transform local matrices if needed */
110   if (pcbddc->use_change_of_basis) {
111     Mat         change_mat_all;
112     PetscScalar *row_cmat_values;
113     PetscInt    *row_cmat_indices;
114     PetscInt    *nnz,*is_indices,*temp_indices;
115     PetscInt    i,j,k,n_D,n_B;
116 
117     /* Get Non-overlapping dimensions */
118     n_B = pcis->n_B;
119     n_D = pcis->n-n_B;
120 
121     /* compute nonzero structure of change of basis on all local nodes */
122     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
123     ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
124     for (i=0;i<n_D;i++) nnz[is_indices[i]] = 1;
125     ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
126     ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
127     k=1;
128     for (i=0;i<n_B;i++) {
129       ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr);
130       nnz[is_indices[i]]=j;
131       if (k < j) k = j;
132       ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr);
133     }
134     ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
135     /* assemble change of basis matrix on the whole set of local dofs */
136     ierr = PetscMalloc(k*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
137     ierr = MatCreate(PETSC_COMM_SELF,&change_mat_all);CHKERRQ(ierr);
138     ierr = MatSetSizes(change_mat_all,pcis->n,pcis->n,pcis->n,pcis->n);CHKERRQ(ierr);
139     ierr = MatSetType(change_mat_all,MATSEQAIJ);CHKERRQ(ierr);
140     ierr = MatSeqAIJSetPreallocation(change_mat_all,0,nnz);CHKERRQ(ierr);
141     ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
142     for (i=0;i<n_D;i++) {
143       ierr = MatSetValue(change_mat_all,is_indices[i],is_indices[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
144     }
145     ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
146     ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
147     for (i=0;i<n_B;i++) {
148       ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr);
149       for (k=0; k<j; k++) temp_indices[k]=is_indices[row_cmat_indices[k]];
150       ierr = MatSetValues(change_mat_all,1,&is_indices[i],j,temp_indices,row_cmat_values,INSERT_VALUES);CHKERRQ(ierr);
151       ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr);
152     }
153     ierr = MatAssemblyBegin(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154     ierr = MatAssemblyEnd(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155     /* TODO: HOW TO WORK WITH BAIJ? PtAP not provided */
156     ierr = MatGetBlockSize(matis->A,&i);CHKERRQ(ierr);
157     if (i==1) {
158       ierr = MatPtAP(matis->A,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
159     } else {
160       Mat work_mat;
161       ierr = MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&work_mat);CHKERRQ(ierr);
162       ierr = MatPtAP(work_mat,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
163       ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
164     }
165     ierr = MatDestroy(&change_mat_all);CHKERRQ(ierr);
166     ierr = PetscFree(nnz);CHKERRQ(ierr);
167     ierr = PetscFree(temp_indices);CHKERRQ(ierr);
168   } else {
169     /* without change of basis, the local matrix is unchanged */
170     if (!pcbddc->local_mat) {
171       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
172       pcbddc->local_mat = matis->A;
173     }
174   }
175 
176   /* get submatrices */
177   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
178   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
179   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
180   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
186 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool use)
187 {
188   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
189 
190   PetscFunctionBegin;
191   pcbddc->use_exact_dirichlet=use;
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "PCBDDCSetUpLocalSolvers"
197 PetscErrorCode PCBDDCSetUpLocalSolvers(PC pc, IS is_I_local, IS is_R_local)
198 {
199   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
200   PC_IS          *pcis = (PC_IS*)pc->data;
201   PC             pc_temp;
202   Mat            A_RR;
203   Vec            vec1,vec2,vec3;
204   MatStructure   matstruct;
205   PetscScalar    m_one = -1.0;
206   PetscReal      value;
207   PetscInt       n_D,n_R,use_exact,use_exact_reduced;
208   PetscErrorCode ierr;
209 
210   PetscFunctionBegin;
211   /* Creating PC contexts for local Dirichlet and Neumann problems */
212   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
213 
214   /* DIRICHLET PROBLEM */
215   /* Matrix for Dirichlet problem is pcis->A_II */
216   ierr = ISGetSize(is_I_local,&n_D);CHKERRQ(ierr);
217   if (!pcbddc->ksp_D) { /* create object if not yet build */
218     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_D);CHKERRQ(ierr);
219     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
220     /* default */
221     ierr = KSPSetType(pcbddc->ksp_D,KSPPREONLY);CHKERRQ(ierr);
222     ierr = KSPSetOptionsPrefix(pcbddc->ksp_D,"dirichlet_");CHKERRQ(ierr);
223     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
224     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
225     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
226   }
227   ierr = KSPSetOperators(pcbddc->ksp_D,pcis->A_II,pcis->A_II,matstruct);CHKERRQ(ierr);
228   /* Allow user's customization */
229   ierr = KSPSetFromOptions(pcbddc->ksp_D);CHKERRQ(ierr);
230   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
231   if (!n_D) {
232     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
233     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
234   }
235   /* Set Up KSP for Dirichlet problem of BDDC */
236   ierr = KSPSetUp(pcbddc->ksp_D);CHKERRQ(ierr);
237   /* set ksp_D into pcis data */
238   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
239   ierr = PetscObjectReference((PetscObject)pcbddc->ksp_D);CHKERRQ(ierr);
240   pcis->ksp_D = pcbddc->ksp_D;
241 
242   /* NEUMANN PROBLEM */
243   /* Matrix for Neumann problem is A_RR -> we need to create it */
244   ierr = ISGetSize(is_R_local,&n_R);CHKERRQ(ierr);
245   ierr = MatGetSubMatrix(pcbddc->local_mat,is_R_local,is_R_local,MAT_INITIAL_MATRIX,&A_RR);CHKERRQ(ierr);
246   if (!pcbddc->ksp_R) { /* create object if not yet build */
247     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_R);CHKERRQ(ierr);
248     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_R,(PetscObject)pc,1);CHKERRQ(ierr);
249     /* default */
250     ierr = KSPSetType(pcbddc->ksp_R,KSPPREONLY);CHKERRQ(ierr);
251     ierr = KSPSetOptionsPrefix(pcbddc->ksp_R,"neumann_");CHKERRQ(ierr);
252     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
253     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
254     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
255   }
256   ierr = KSPSetOperators(pcbddc->ksp_R,A_RR,A_RR,matstruct);CHKERRQ(ierr);
257   /* Allow user's customization */
258   ierr = KSPSetFromOptions(pcbddc->ksp_R);CHKERRQ(ierr);
259   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
260   if (!n_R) {
261     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
262     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
263   }
264   /* Set Up KSP for Neumann problem of BDDC */
265   ierr = KSPSetUp(pcbddc->ksp_R);CHKERRQ(ierr);
266 
267   /* check Dirichlet and Neumann solvers and adapt them if a nullspace correction is needed */
268 
269   /* Dirichlet */
270   ierr = MatGetVecs(pcis->A_II,&vec1,&vec2);CHKERRQ(ierr);
271   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
272   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
273   ierr = MatMult(pcis->A_II,vec1,vec2);CHKERRQ(ierr);
274   ierr = KSPSolve(pcbddc->ksp_D,vec2,vec3);CHKERRQ(ierr);
275   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
276   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
277   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
278   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
279   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
280   /* need to be adapted? */
281   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
282   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
283   ierr = PCBDDCSetUseExactDirichlet(pc,(PetscBool)use_exact_reduced);CHKERRQ(ierr);
284   /* print info */
285   if (pcbddc->dbg_flag) {
286     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
287     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
288     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Checking solution of Dirichlet and Neumann problems\n");CHKERRQ(ierr);
289     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
290     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
291   }
292   if (n_D && pcbddc->NullSpace && !use_exact_reduced && !pcbddc->inexact_prec_type) {
293     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_I_local);CHKERRQ(ierr);
294   }
295 
296   /* Neumann */
297   ierr = MatGetVecs(A_RR,&vec1,&vec2);CHKERRQ(ierr);
298   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
299   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
300   ierr = MatMult(A_RR,vec1,vec2);CHKERRQ(ierr);
301   ierr = KSPSolve(pcbddc->ksp_R,vec2,vec3);CHKERRQ(ierr);
302   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
303   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
304   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
305   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
306   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
307   /* need to be adapted? */
308   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
309   if (PetscAbsReal(value) > 1.e-4) use_exact = 0;
310   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
311   /* print info */
312   if (pcbddc->dbg_flag) {
313     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for  Neumann  solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
314     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
315   }
316   if (n_R && pcbddc->NullSpace && !use_exact_reduced) { /* is it the right logic? */
317     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_R_local);CHKERRQ(ierr);
318   }
319 
320   /* free Neumann problem's matrix */
321   ierr = MatDestroy(&A_RR);CHKERRQ(ierr);
322   PetscFunctionReturn(0);
323 }
324 
325 #undef __FUNCT__
326 #define __FUNCT__ "PCBDDCSolveSaddlePoint"
327 static PetscErrorCode  PCBDDCSolveSaddlePoint(PC pc)
328 {
329   PetscErrorCode ierr;
330   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
331 
332   PetscFunctionBegin;
333   ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr);
334   if (pcbddc->local_auxmat1) {
335     ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr);
336     ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr);
337   }
338   PetscFunctionReturn(0);
339 }
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner"
343 PetscErrorCode  PCBDDCApplyInterfacePreconditioner(PC pc)
344 {
345   PetscErrorCode ierr;
346   PC_BDDC*        pcbddc = (PC_BDDC*)(pc->data);
347   PC_IS*            pcis = (PC_IS*)  (pc->data);
348   const PetscScalar zero = 0.0;
349 
350   PetscFunctionBegin;
351   /* Application of PHI^T (or PSI^T)  */
352   if (pcbddc->coarse_psi_B) {
353     ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
354     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
355   } else {
356     ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
357     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
358   }
359   /* Scatter data of coarse_rhs */
360   if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); }
361   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
362 
363   /* Local solution on R nodes */
364   ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr);
365   ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
366   ierr = VecScatterEnd  (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
367   if (pcbddc->inexact_prec_type) {
368     ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
369     ierr = VecScatterEnd  (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
370   }
371   ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr);
372   ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
373   ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
374   ierr = VecScatterEnd  (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
375   if (pcbddc->inexact_prec_type) {
376     ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
377     ierr = VecScatterEnd  (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
378   }
379 
380   /* Coarse solution */
381   ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
382   if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */
383     ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr);
384   }
385   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
386   ierr = PCBDDCScatterCoarseDataEnd  (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
387 
388   /* Sum contributions from two levels */
389   ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
390   if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
391   PetscFunctionReturn(0);
392 }
393 
394 #undef __FUNCT__
395 #define __FUNCT__ "PCBDDCScatterCoarseDataBegin"
396 PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
397 {
398   PetscErrorCode ierr;
399   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
400 
401   PetscFunctionBegin;
402   switch (pcbddc->coarse_communications_type) {
403     case SCATTERS_BDDC:
404       ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
405       break;
406     case GATHERS_BDDC:
407       break;
408   }
409   PetscFunctionReturn(0);
410 }
411 
412 #undef __FUNCT__
413 #define __FUNCT__ "PCBDDCScatterCoarseDataEnd"
414 PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
415 {
416   PetscErrorCode ierr;
417   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
418   PetscScalar*   array_to;
419   PetscScalar*   array_from;
420   MPI_Comm       comm;
421   PetscInt       i;
422 
423   PetscFunctionBegin;
424   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
425   switch (pcbddc->coarse_communications_type) {
426     case SCATTERS_BDDC:
427       ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
428       break;
429     case GATHERS_BDDC:
430       if (vec_from) {
431         ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr);
432       }
433       if (vec_to) {
434         ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr);
435       }
436       switch(pcbddc->coarse_problem_type){
437         case SEQUENTIAL_BDDC:
438           if (smode == SCATTER_FORWARD) {
439             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);
440             if (vec_to) {
441               if (imode == ADD_VALUES) {
442                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
443                   array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
444                 }
445               } else {
446                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
447                   array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
448                 }
449               }
450             }
451           } else {
452             if (vec_from) {
453               if (imode == ADD_VALUES) {
454                 MPI_Comm vec_from_comm;
455                 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr);
456                 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);
457               }
458               for (i=0;i<pcbddc->replicated_primal_size;i++) {
459                 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]];
460               }
461             }
462             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);
463           }
464           break;
465         case REPLICATED_BDDC:
466           if (smode == SCATTER_FORWARD) {
467             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);
468             if (imode == ADD_VALUES) {
469               for (i=0;i<pcbddc->replicated_primal_size;i++) {
470                 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
471               }
472             } else {
473               for (i=0;i<pcbddc->replicated_primal_size;i++) {
474                 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
475               }
476             }
477           } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */
478             if (imode == ADD_VALUES) {
479               for (i=0;i<pcbddc->local_primal_size;i++) {
480                 array_to[i]+=array_from[pcbddc->local_primal_indices[i]];
481               }
482             } else {
483               for (i=0;i<pcbddc->local_primal_size;i++) {
484                 array_to[i]=array_from[pcbddc->local_primal_indices[i]];
485               }
486             }
487           }
488           break;
489         case MULTILEVEL_BDDC:
490           break;
491         case PARALLEL_BDDC:
492           break;
493       }
494       if (vec_from) {
495         ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr);
496       }
497       if (vec_to) {
498         ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr);
499       }
500       break;
501   }
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "PCBDDCConstraintsSetUp"
507 PetscErrorCode PCBDDCConstraintsSetUp(PC pc)
508 {
509   PetscErrorCode ierr;
510   PC_IS*         pcis = (PC_IS*)(pc->data);
511   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
512   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
513   PetscInt       *nnz,*is_indices;
514   PetscScalar    *temp_quadrature_constraint;
515   PetscInt       *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B;
516   PetscInt       local_primal_size,i,j,k,total_counts,max_size_of_constraint;
517   PetscInt       n_vertices,size_of_constraint;
518   PetscReal      real_value;
519   PetscBool      nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true;
520   PetscInt       nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges;
521   IS             *used_IS,ISForVertices,*ISForFaces,*ISForEdges;
522   MatType        impMatType=MATSEQAIJ;
523   PetscBLASInt   Bs,Bt,lwork,lierr;
524   PetscReal      tol=1.0e-8;
525   MatNullSpace   nearnullsp;
526   const Vec      *nearnullvecs;
527   Vec            *localnearnullsp;
528   PetscScalar    *work,*temp_basis,*array_vector,*correlation_mat;
529   PetscReal      *rwork,*singular_vals;
530   PetscBLASInt   Bone=1,*ipiv;
531   Vec            temp_vec;
532   Mat            temp_mat;
533   KSP            temp_ksp;
534   PC             temp_pc;
535   PetscInt       s,start_constraint,dual_dofs;
536   PetscBool      compute_submatrix,useksp=PETSC_FALSE;
537   PetscInt       *aux_primal_permutation,*aux_primal_numbering;
538   PetscBool      boolforchange,*change_basis;
539 /* some ugly conditional declarations */
540 #if defined(PETSC_MISSING_LAPACK_GESVD)
541   PetscScalar    one=1.0,zero=0.0;
542   PetscInt       ii;
543   PetscScalar    *singular_vectors;
544   PetscBLASInt   *iwork,*ifail;
545   PetscReal      dummy_real,abs_tol;
546   PetscBLASInt   eigs_found;
547 #endif
548   PetscBLASInt   dummy_int;
549   PetscScalar    dummy_scalar;
550   PetscBool      used_vertex,get_faces,get_edges,get_vertices;
551 
552   PetscFunctionBegin;
553   /* Get index sets for faces, edges and vertices from graph */
554   get_faces = PETSC_TRUE;
555   get_edges = PETSC_TRUE;
556   get_vertices = PETSC_TRUE;
557   if (pcbddc->vertices_flag) {
558     get_faces = PETSC_FALSE;
559     get_edges = PETSC_FALSE;
560   }
561   if (pcbddc->constraints_flag) {
562     get_vertices = PETSC_FALSE;
563   }
564   if (pcbddc->faces_flag) {
565     get_edges = PETSC_FALSE;
566   }
567   if (pcbddc->edges_flag) {
568     get_faces = PETSC_FALSE;
569   }
570   /* default */
571   if (!get_faces && !get_edges && !get_vertices) {
572     get_vertices = PETSC_TRUE;
573   }
574   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
575   if (pcbddc->dbg_flag) {
576     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
577     i = 0;
578     if (ISForVertices) {
579       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
580     }
581     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
582     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
583     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
584     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
585   }
586   /* check if near null space is attached to global mat */
587   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
588   if (nearnullsp) {
589     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
590   } else { /* if near null space is not provided it uses constants */
591     nnsp_has_cnst = PETSC_TRUE;
592     use_nnsp_true = PETSC_TRUE;
593   }
594   if (nnsp_has_cnst) {
595     nnsp_addone = 1;
596   }
597   /*
598        Evaluate maximum storage size needed by the procedure
599        - temp_indices will contain start index of each constraint stored as follows
600        - temp_indices_to_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts
601        - 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
602        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
603                                                                                                                                                          */
604   total_counts = n_ISForFaces+n_ISForEdges;
605   total_counts *= (nnsp_addone+nnsp_size);
606   n_vertices = 0;
607   if (ISForVertices) {
608     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
609   }
610   total_counts += n_vertices;
611   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
612   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
613   total_counts = 0;
614   max_size_of_constraint = 0;
615   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
616     if (i<n_ISForEdges) {
617       used_IS = &ISForEdges[i];
618     } else {
619       used_IS = &ISForFaces[i-n_ISForEdges];
620     }
621     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
622     total_counts += j;
623     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
624   }
625   total_counts *= (nnsp_addone+nnsp_size);
626   total_counts += n_vertices;
627   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
628   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
629   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
630   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
631   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
632   for (i=0;i<pcis->n;i++) {
633     local_to_B[i]=-1;
634   }
635   for (i=0;i<pcis->n_B;i++) {
636     local_to_B[is_indices[i]]=i;
637   }
638   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
639 
640   /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */
641   rwork = 0;
642   work = 0;
643   singular_vals = 0;
644   temp_basis = 0;
645   correlation_mat = 0;
646   if (!pcbddc->use_nnsp_true) {
647     PetscScalar temp_work;
648 #if defined(PETSC_MISSING_LAPACK_GESVD)
649     /* POD */
650     PetscInt max_n;
651     max_n = nnsp_addone+nnsp_size;
652     /* using some techniques borrowed from Proper Orthogonal Decomposition */
653     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
654     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr);
655     ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
656     ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
657 #if defined(PETSC_USE_COMPLEX)
658     ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
659 #endif
660     ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr);
661     ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr);
662     /* now we evaluate the optimal workspace using query with lwork=-1 */
663     ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr);
664     lwork=-1;
665     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
666 #if !defined(PETSC_USE_COMPLEX)
667     abs_tol=1.e-8;
668 /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */
669     PetscStackCallBLAS("LAPACKsyevx",LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int,&abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,&temp_work,&lwork,iwork,ifail,&lierr));
670 #else
671 /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */
672 /*  LAPACK call is missing here! TODO */
673     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
674 #endif
675     if ( lierr ) {
676       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr);
677     }
678     ierr = PetscFPTrapPop();CHKERRQ(ierr);
679 #else /* on missing GESVD */
680     /* SVD */
681     PetscInt max_n,min_n;
682     max_n = max_size_of_constraint;
683     min_n = nnsp_addone+nnsp_size;
684     if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) {
685       min_n = max_size_of_constraint;
686       max_n = nnsp_addone+nnsp_size;
687     }
688     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
689 #if defined(PETSC_USE_COMPLEX)
690     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
691 #endif
692     /* now we evaluate the optimal workspace using query with lwork=-1 */
693     lwork=-1;
694     ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr);
695     ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr);
696     dummy_int = Bs;
697     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
698 #if !defined(PETSC_USE_COMPLEX)
699     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,&lierr));
700 #else
701     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[0],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,&temp_work,&lwork,rwork,&lierr));
702 #endif
703     if ( lierr ) {
704       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr);
705     }
706     ierr = PetscFPTrapPop();CHKERRQ(ierr);
707 #endif
708     /* Allocate optimal workspace */
709     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
710     total_counts = (PetscInt)lwork;
711     ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr);
712   }
713   /* get local part of global near null space vectors */
714   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
715   for (k=0;k<nnsp_size;k++) {
716     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
717     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
718     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
719   }
720   /* Now we can loop on constraining sets */
721   total_counts = 0;
722   temp_indices[0] = 0;
723   /* vertices */
724   if (ISForVertices) {
725     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
726     if (nnsp_has_cnst) { /* consider all vertices */
727       for (i=0;i<n_vertices;i++) {
728         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
729         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
730         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
731         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
732         change_basis[total_counts]=PETSC_FALSE;
733         total_counts++;
734       }
735     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
736       for (i=0;i<n_vertices;i++) {
737         used_vertex=PETSC_FALSE;
738         k=0;
739         while (!used_vertex && k<nnsp_size) {
740           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
741           if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) {
742             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
743             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
744             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
745             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
746             change_basis[total_counts]=PETSC_FALSE;
747             total_counts++;
748             used_vertex=PETSC_TRUE;
749           }
750           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
751           k++;
752         }
753       }
754     }
755     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
756     n_vertices = total_counts;
757   }
758   /* edges and faces */
759   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
760     if (i<n_ISForEdges) {
761       used_IS = &ISForEdges[i];
762       boolforchange = pcbddc->use_change_of_basis;
763     } else {
764       used_IS = &ISForFaces[i-n_ISForEdges];
765       boolforchange = pcbddc->use_change_on_faces;
766     }
767     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
768     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
769     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
770     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
771     /* HACK: change of basis should not performed on local periodic nodes */
772     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) {
773       boolforchange = PETSC_FALSE;
774     }
775     if (nnsp_has_cnst) {
776       PetscScalar quad_value;
777       temp_constraints++;
778       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
779       for (j=0;j<size_of_constraint;j++) {
780         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
781         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
782         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
783       }
784       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
785       change_basis[total_counts]=boolforchange;
786       total_counts++;
787     }
788     for (k=0;k<nnsp_size;k++) {
789       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
790       for (j=0;j<size_of_constraint;j++) {
791         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
792         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
793         temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]];
794       }
795       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
796       real_value = 1.0;
797       if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */
798         ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
799         PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone));
800       }
801       if (real_value > 0.0) { /* keep indices and values */
802         temp_constraints++;
803         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
804         change_basis[total_counts]=boolforchange;
805         total_counts++;
806       }
807     }
808     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
809     /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */
810     if (!use_nnsp_true) {
811       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
812       ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr);
813 
814 #if defined(PETSC_MISSING_LAPACK_GESVD)
815       ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr);
816       /* Store upper triangular part of correlation matrix */
817       for (j=0;j<temp_constraints;j++) {
818         for (k=0;k<j+1;k++) {
819           PetscStackCallBLAS("BLASdot",correlation_mat[j*temp_constraints+k]=BLASdot_(&Bs,&temp_quadrature_constraint[temp_indices[temp_start_ptr+j]],&Bone,&temp_quadrature_constraint[temp_indices[temp_start_ptr+k]],&Bone));
820 
821         }
822       }
823       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
824 #if !defined(PETSC_USE_COMPLEX)
825 /*      LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */
826       PetscStackCallBLAS("LAPACKsyevx",LAPACKsyevx_("V","A","U",&Bt,correlation_mat,&Bt,&dummy_real,&dummy_real,&dummy_int,&dummy_int,&abs_tol,&eigs_found,singular_vals,singular_vectors,&Bt,work,&lwork,iwork,ifail,&lierr));
827 #else
828 /*  LAPACK call is missing here! TODO */
829       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
830 #endif
831       if (lierr) {
832         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr);
833       }
834       ierr = PetscFPTrapPop();CHKERRQ(ierr);
835       /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */
836       j=0;
837       while (j < Bt && singular_vals[j] < tol) j++;
838       total_counts=total_counts-j;
839       if (j<temp_constraints) {
840         for (k=j;k<Bt;k++) {
841           singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
842         }
843         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
844         PetscStackCallBLAS("BLASgemm_",BLASgemm_("N","N",&Bs,&Bt,&Bt,&one,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,correlation_mat,&Bt,&zero,temp_basis,&Bs));
845         ierr = PetscFPTrapPop();CHKERRQ(ierr);
846         /* copy POD basis into used quadrature memory */
847         for (k=0;k<Bt-j;k++) {
848           for (ii=0;ii<size_of_constraint;ii++) {
849             temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[Bt-1-k]*temp_basis[(Bt-1-k)*size_of_constraint+ii];
850           }
851         }
852       }
853 
854 #else  /* on missing GESVD */
855       PetscInt min_n = temp_constraints;
856       if (min_n > size_of_constraint) min_n = size_of_constraint;
857       dummy_int = Bs;
858       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
859 #if !defined(PETSC_USE_COMPLEX)
860       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,&lierr));
861 #else
862       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&Bs,&Bt,&temp_quadrature_constraint[temp_indices[temp_start_ptr]],&Bs,singular_vals,&dummy_scalar,&dummy_int,&dummy_scalar,&dummy_int,work,&lwork,rwork,&lierr));
863 #endif
864       if (lierr) {
865         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr);
866       }
867       ierr = PetscFPTrapPop();CHKERRQ(ierr);
868       /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */
869       j = 0;
870       while (j < min_n && singular_vals[min_n-j-1] < tol) j++;
871       total_counts = total_counts-(PetscInt)Bt+(min_n-j);
872 #endif
873     }
874   }
875   /* free index sets of faces, edges and vertices */
876   for (i=0;i<n_ISForFaces;i++) {
877     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
878   }
879   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
880   for (i=0;i<n_ISForEdges;i++) {
881     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
882   }
883   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
884   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
885 
886   /* set quantities in pcbddc data structure */
887   /* n_vertices defines the number of point primal dofs */
888   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
889   local_primal_size = total_counts;
890   pcbddc->n_vertices = n_vertices;
891   pcbddc->n_constraints = total_counts-n_vertices;
892   pcbddc->local_primal_size = local_primal_size;
893 
894   /* Create constraint matrix */
895   /* The constraint matrix is used to compute the l2g map of primal dofs */
896   /* so we need to set it up properly either with or without change of basis */
897   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
898   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
899   ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr);
900   /* compute a local numbering of constraints : vertices first then constraints */
901   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
902   ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
903   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
904   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr);
905   total_counts=0;
906   /* find vertices: subdomain corners plus dofs with basis changed */
907   for (i=0;i<local_primal_size;i++) {
908     size_of_constraint=temp_indices[i+1]-temp_indices[i];
909     if (change_basis[i] || size_of_constraint == 1) {
910       k=0;
911       while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) {
912         k=k+1;
913       }
914       j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1];
915       array_vector[j] = 1.0;
916       aux_primal_numbering[total_counts]=j;
917       aux_primal_permutation[total_counts]=total_counts;
918       total_counts++;
919     }
920   }
921   ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
922   /* permute indices in order to have a sorted set of vertices */
923   ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation);
924   /* nonzero structure */
925   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
926   for (i=0;i<total_counts;i++) {
927     nnz[i]=1;
928   }
929   j=total_counts;
930   for (i=n_vertices;i<local_primal_size;i++) {
931     if (!change_basis[i]) {
932       nnz[j]=temp_indices[i+1]-temp_indices[i];
933       j++;
934     }
935   }
936   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
937   ierr = PetscFree(nnz);CHKERRQ(ierr);
938   /* set values in constraint matrix */
939   for (i=0;i<total_counts;i++) {
940     j = aux_primal_permutation[i];
941     k = aux_primal_numbering[j];
942     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr);
943   }
944   for (i=n_vertices;i<local_primal_size;i++) {
945     if (!change_basis[i]) {
946       size_of_constraint=temp_indices[i+1]-temp_indices[i];
947       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);
948       total_counts++;
949     }
950   }
951   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
952   ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr);
953   /* assembling */
954   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
955   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
956 
957   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
958   if (pcbddc->use_change_of_basis) {
959     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
960     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
961     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
962     /* work arrays */
963     /* we need to reuse these arrays, so we free them */
964     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
965     ierr = PetscFree(work);CHKERRQ(ierr);
966     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
967     ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
968     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr);
969     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr);
970     for (i=0;i<pcis->n_B;i++) {
971       nnz[i]=1;
972     }
973     /* Overestimated nonzeros per row */
974     k=1;
975     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
976       if (change_basis[i]) {
977         size_of_constraint = temp_indices[i+1]-temp_indices[i];
978         if (k < size_of_constraint) {
979           k = size_of_constraint;
980         }
981         for (j=0;j<size_of_constraint;j++) {
982           nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
983         }
984       }
985     }
986     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
987     ierr = PetscFree(nnz);CHKERRQ(ierr);
988     /* Temporary array to store indices */
989     ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr);
990     /* Set initial identity in the matrix */
991     for (i=0;i<pcis->n_B;i++) {
992       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
993     }
994     /* Now we loop on the constraints which need a change of basis */
995     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
996     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */
997     temp_constraints = 0;
998     if (pcbddc->n_vertices < local_primal_size) {
999       temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]];
1000     }
1001     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
1002       if (change_basis[i]) {
1003         compute_submatrix = PETSC_FALSE;
1004         useksp = PETSC_FALSE;
1005         if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) {
1006           temp_constraints++;
1007           if (i == local_primal_size -1 ||  temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) {
1008             compute_submatrix = PETSC_TRUE;
1009           }
1010         }
1011         if (compute_submatrix) {
1012           if (temp_constraints > 1 || pcbddc->use_nnsp_true) {
1013             useksp = PETSC_TRUE;
1014           }
1015           size_of_constraint = temp_indices[i+1]-temp_indices[i];
1016           if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */
1017             ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr);
1018             ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr);
1019             ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr);
1020             ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr);
1021           }
1022           /* First _size_of_constraint-temp_constraints_ columns */
1023           dual_dofs = size_of_constraint-temp_constraints;
1024           start_constraint = i+1-temp_constraints;
1025           for (s=0;s<dual_dofs;s++) {
1026             is_indices[0] = s;
1027             for (j=0;j<temp_constraints;j++) {
1028               for (k=0;k<temp_constraints;k++) {
1029                 temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1];
1030               }
1031               work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s];
1032               is_indices[j+1]=s+j+1;
1033             }
1034             Bt = temp_constraints;
1035             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1036             PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr));
1037             if ( lierr ) {
1038               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr);
1039             }
1040             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1041             j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s];
1042             ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,temp_constraints,&temp_indices_to_constraint_B[temp_indices[start_constraint]+s+1],1,&j,work,INSERT_VALUES);CHKERRQ(ierr);
1043             if (useksp) {
1044               /* temp mat with transposed rows and columns */
1045               ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr);
1046               ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr);
1047             }
1048           }
1049           if (useksp) {
1050             /* last rows of temp_mat */
1051             for (j=0;j<size_of_constraint;j++) {
1052               is_indices[j] = j;
1053             }
1054             for (s=0;s<temp_constraints;s++) {
1055               k = s + dual_dofs;
1056               ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr);
1057             }
1058             ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1059             ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1060             ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr);
1061             ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr);
1062             ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1063             ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr);
1064             ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr);
1065             ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr);
1066             ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr);
1067             for (s=0;s<temp_constraints;s++) {
1068               ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr);
1069               ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr);
1070               ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr);
1071               ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr);
1072               ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr);
1073               ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr);
1074               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
1075               /* last columns of change of basis matrix associated to new primal dofs */
1076               ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,array_vector,INSERT_VALUES);CHKERRQ(ierr);
1077               ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr);
1078             }
1079             ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
1080             ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr);
1081             ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1082           } else {
1083             /* last columns of change of basis matrix associated to new primal dofs */
1084             for (s=0;s<temp_constraints;s++) {
1085               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
1086               ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,size_of_constraint,&temp_indices_to_constraint_B[temp_indices[start_constraint+s]],1,&j,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr);
1087             }
1088           }
1089           /* prepare for the next cycle */
1090           temp_constraints = 0;
1091           if (i != local_primal_size -1 ) {
1092             temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]];
1093           }
1094         }
1095       }
1096     }
1097     /* assembling */
1098     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1099     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1100     ierr = PetscFree(ipiv);CHKERRQ(ierr);
1101     ierr = PetscFree(is_indices);CHKERRQ(ierr);
1102   }
1103   /* free workspace no longer needed */
1104   ierr = PetscFree(rwork);CHKERRQ(ierr);
1105   ierr = PetscFree(work);CHKERRQ(ierr);
1106   ierr = PetscFree(temp_basis);CHKERRQ(ierr);
1107   ierr = PetscFree(singular_vals);CHKERRQ(ierr);
1108   ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
1109   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
1110   ierr = PetscFree(change_basis);CHKERRQ(ierr);
1111   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
1112   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
1113   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
1114   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
1115 #if defined(PETSC_MISSING_LAPACK_GESVD)
1116   ierr = PetscFree(iwork);CHKERRQ(ierr);
1117   ierr = PetscFree(ifail);CHKERRQ(ierr);
1118   ierr = PetscFree(singular_vectors);CHKERRQ(ierr);
1119 #endif
1120   for (k=0;k<nnsp_size;k++) {
1121     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
1122   }
1123   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNCT__
1128 #define __FUNCT__ "PCBDDCAnalyzeInterface"
1129 PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
1130 {
1131   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
1132   PC_IS       *pcis = (PC_IS*)pc->data;
1133   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
1134   PetscInt    bs,ierr,i,vertex_size;
1135   PetscViewer viewer=pcbddc->dbg_viewer;
1136 
1137   PetscFunctionBegin;
1138   /* Init local Graph struct */
1139   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
1140 
1141   /* Check validity of the csr graph passed in by the user */
1142   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
1143     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
1144   }
1145   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
1146   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
1147     Mat mat_adj;
1148     const PetscInt *xadj,*adjncy;
1149     PetscBool flg_row=PETSC_TRUE;
1150 
1151     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
1152     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1153     if (!flg_row) {
1154       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
1155     }
1156     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
1157     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1158     if (!flg_row) {
1159       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
1160     }
1161     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
1162   }
1163 
1164   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
1165   vertex_size = 1;
1166   if (!pcbddc->n_ISForDofs) {
1167     IS *custom_ISForDofs;
1168 
1169     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
1170     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
1171     for (i=0;i<bs;i++) {
1172       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
1173     }
1174     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
1175     /* remove my references to IS objects */
1176     for (i=0;i<bs;i++) {
1177       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
1178     }
1179     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
1180   } else { /* mat block size as vertex size (used for elasticity) */
1181     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
1182   }
1183 
1184   /* Setup of Graph */
1185   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
1186 
1187   /* Graph's connected components analysis */
1188   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
1189 
1190   /* print some info to stdout */
1191   if (pcbddc->dbg_flag) {
1192     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
1193   }
1194   PetscFunctionReturn(0);
1195 }
1196 
1197 #undef __FUNCT__
1198 #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
1199 PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
1200 {
1201   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1202   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
1203   PetscErrorCode ierr;
1204 
1205   PetscFunctionBegin;
1206   n = 0;
1207   vertices = 0;
1208   if (pcbddc->ConstraintMatrix) {
1209     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
1210     for (i=0;i<local_primal_size;i++) {
1211       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1212       if (size_of_constraint == 1) n++;
1213       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1214     }
1215     ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
1216     n = 0;
1217     for (i=0;i<local_primal_size;i++) {
1218       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1219       if (size_of_constraint == 1) {
1220         vertices[n++]=row_cmat_indices[0];
1221       }
1222       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1223     }
1224   }
1225   *n_vertices = n;
1226   *vertices_idx = vertices;
1227   PetscFunctionReturn(0);
1228 }
1229 
1230 #undef __FUNCT__
1231 #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
1232 PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
1233 {
1234   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1235   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
1236   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
1237   PetscBool      *touched;
1238   PetscErrorCode ierr;
1239 
1240   PetscFunctionBegin;
1241   n = 0;
1242   constraints_index = 0;
1243   if (pcbddc->ConstraintMatrix) {
1244     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
1245     max_size_of_constraint = 0;
1246     for (i=0;i<local_primal_size;i++) {
1247       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1248       if (size_of_constraint > 1) {
1249         n++;
1250       }
1251       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1252       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1253     }
1254     ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1255     ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1256     ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1257     ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1258     n = 0;
1259     for (i=0;i<local_primal_size;i++) {
1260       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1261       if (size_of_constraint > 1) {
1262         ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
1263         min_index = row_cmat_global_indices[0];
1264         min_loc = 0;
1265         for (j=1;j<size_of_constraint;j++) {
1266           /* there can be more than one constraint on a single connected component */
1267           if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1268             min_index = row_cmat_global_indices[j];
1269             min_loc = j;
1270           }
1271         }
1272         touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1273         constraints_index[n++] = row_cmat_indices[min_loc];
1274       }
1275       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1276     }
1277   }
1278   ierr = PetscFree(touched);CHKERRQ(ierr);
1279   ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1280   *n_constraints = n;
1281   *constraints_idx = constraints_index;
1282   PetscFunctionReturn(0);
1283 }
1284 
1285 /* the next two functions has been adapted from pcis.c */
1286 #undef __FUNCT__
1287 #define __FUNCT__ "PCBDDCApplySchur"
1288 PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1289 {
1290   PetscErrorCode ierr;
1291   PC_IS          *pcis = (PC_IS*)(pc->data);
1292 
1293   PetscFunctionBegin;
1294   if (!vec2_B) { vec2_B = v; }
1295   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1296   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1297   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1298   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1299   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1300   PetscFunctionReturn(0);
1301 }
1302 
1303 #undef __FUNCT__
1304 #define __FUNCT__ "PCBDDCApplySchurTranspose"
1305 PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1306 {
1307   PetscErrorCode ierr;
1308   PC_IS          *pcis = (PC_IS*)(pc->data);
1309 
1310   PetscFunctionBegin;
1311   if (!vec2_B) { vec2_B = v; }
1312   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1313   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1314   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1315   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1316   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1317   PetscFunctionReturn(0);
1318 }
1319 
1320 #undef __FUNCT__
1321 #define __FUNCT__ "PCBDDCSubsetNumbering"
1322 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[])
1323 {
1324   Vec            local_vec,global_vec;
1325   IS             seqis,paris;
1326   VecScatter     scatter_ctx;
1327   PetscScalar    *array;
1328   PetscInt       *temp_global_dofs;
1329   PetscScalar    globalsum;
1330   PetscInt       i,j,s;
1331   PetscInt       nlocals,first_index,old_index,max_local;
1332   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1333   PetscMPIInt    *dof_sizes,*dof_displs;
1334   PetscBool      first_found;
1335   PetscErrorCode ierr;
1336 
1337   PetscFunctionBegin;
1338   /* mpi buffers */
1339   MPI_Comm_size(comm,&size_prec_comm);
1340   MPI_Comm_rank(comm,&rank_prec_comm);
1341   j = ( !rank_prec_comm ? size_prec_comm : 0);
1342   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1343   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1344   /* get maximum size of subset */
1345   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1346   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1347   max_local = 0;
1348   if (n_local_dofs) {
1349     max_local = temp_global_dofs[0];
1350     for (i=1;i<n_local_dofs;i++) {
1351       if (max_local < temp_global_dofs[i] ) {
1352         max_local = temp_global_dofs[i];
1353       }
1354     }
1355   }
1356   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1357   max_global++;
1358   max_local = 0;
1359   if (n_local_dofs) {
1360     max_local = local_dofs[0];
1361     for (i=1;i<n_local_dofs;i++) {
1362       if (max_local < local_dofs[i] ) {
1363         max_local = local_dofs[i];
1364       }
1365     }
1366   }
1367   max_local++;
1368   /* allocate workspace */
1369   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1370   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1371   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1372   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1373   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1374   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1375   /* create scatter */
1376   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1377   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1378   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1379   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1380   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1381   /* init array */
1382   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1383   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1384   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1385   if (local_dofs_mult) {
1386     for (i=0;i<n_local_dofs;i++) {
1387       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1388     }
1389   } else {
1390     for (i=0;i<n_local_dofs;i++) {
1391       array[local_dofs[i]]=1.0;
1392     }
1393   }
1394   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1395   /* scatter into global vec and get total number of global dofs */
1396   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1397   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1398   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
1399   *n_global_subset = (PetscInt)PetscRealPart(globalsum);
1400   /* Fill global_vec with cumulative function for global numbering */
1401   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1402   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1403   nlocals = 0;
1404   first_index = -1;
1405   first_found = PETSC_FALSE;
1406   for (i=0;i<s;i++) {
1407     if (!first_found && PetscRealPart(array[i]) > 0.0) {
1408       first_found = PETSC_TRUE;
1409       first_index = i;
1410     }
1411     nlocals += (PetscInt)PetscRealPart(array[i]);
1412   }
1413   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1414   if (!rank_prec_comm) {
1415     dof_displs[0]=0;
1416     for (i=1;i<size_prec_comm;i++) {
1417       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1418     }
1419   }
1420   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1421   if (first_found) {
1422     array[first_index] += (PetscScalar)nlocals;
1423     old_index = first_index;
1424     for (i=first_index+1;i<s;i++) {
1425       if (PetscRealPart(array[i]) > 0.0) {
1426         array[i] += array[old_index];
1427         old_index = i;
1428       }
1429     }
1430   }
1431   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1432   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1433   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1434   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1435   /* get global ordering of local dofs */
1436   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1437   if (local_dofs_mult) {
1438     for (i=0;i<n_local_dofs;i++) {
1439       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i];
1440     }
1441   } else {
1442     for (i=0;i<n_local_dofs;i++) {
1443       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1;
1444     }
1445   }
1446   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1447   /* free workspace */
1448   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1449   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1450   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1451   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1452   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1453   /* return pointer to global ordering of local dofs */
1454   *global_numbering_subset = temp_global_dofs;
1455   PetscFunctionReturn(0);
1456 }
1457