xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision 83b7ccabe1fbecb8bf298e2ac937c80b54d2b2de)
1674ae819SStefano Zampini #include "bddc.h"
2674ae819SStefano Zampini #include "bddcprivate.h"
3674ae819SStefano Zampini #include <petscblaslapack.h>
4674ae819SStefano Zampini 
5674ae819SStefano Zampini #undef __FUNCT__
6674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetCustomization"
7674ae819SStefano Zampini PetscErrorCode PCBDDCResetCustomization(PC pc)
8674ae819SStefano Zampini {
9674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
10674ae819SStefano Zampini   PetscInt       i;
11674ae819SStefano Zampini   PetscErrorCode ierr;
12674ae819SStefano Zampini 
13674ae819SStefano Zampini   PetscFunctionBegin;
14674ae819SStefano Zampini   ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
15674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr);
16674ae819SStefano Zampini   ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr);
17674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr);
18674ae819SStefano Zampini   ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr);
19674ae819SStefano Zampini   for (i=0;i<pcbddc->n_ISForDofs;i++) {
20674ae819SStefano Zampini     ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr);
21674ae819SStefano Zampini   }
22674ae819SStefano Zampini   ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr);
23674ae819SStefano Zampini   PetscFunctionReturn(0);
24674ae819SStefano Zampini }
25674ae819SStefano Zampini 
26674ae819SStefano Zampini #undef __FUNCT__
27674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetTopography"
28674ae819SStefano Zampini PetscErrorCode PCBDDCResetTopography(PC pc)
29674ae819SStefano Zampini {
30674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
31674ae819SStefano Zampini   PetscErrorCode ierr;
32674ae819SStefano Zampini 
33674ae819SStefano Zampini   PetscFunctionBegin;
34674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
35674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
36674ae819SStefano Zampini   ierr = PCBDDCGraphReset(pcbddc->mat_graph);CHKERRQ(ierr);
37674ae819SStefano Zampini   PetscFunctionReturn(0);
38674ae819SStefano Zampini }
39674ae819SStefano Zampini 
40674ae819SStefano Zampini #undef __FUNCT__
41674ae819SStefano Zampini #define __FUNCT__ "PCBDDCResetSolvers"
42674ae819SStefano Zampini PetscErrorCode PCBDDCResetSolvers(PC pc)
43674ae819SStefano Zampini {
44674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
45674ae819SStefano Zampini   PetscErrorCode ierr;
46674ae819SStefano Zampini 
47674ae819SStefano Zampini   PetscFunctionBegin;
48674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr);
49674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr);
50674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr);
51674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr);
52674ae819SStefano Zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
53674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr);
54674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr);
55674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr);
5615aaf578SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr);
5715aaf578SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr);
58674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr);
59674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr);
60674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr);
61674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr);
62674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr);
63674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr);
64674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr);
65674ae819SStefano Zampini   ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr);
66674ae819SStefano Zampini   ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr);
67674ae819SStefano Zampini   ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr);
68674ae819SStefano Zampini   ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr);
69674ae819SStefano Zampini   ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr);
70674ae819SStefano Zampini   ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr);
71674ae819SStefano Zampini   ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr);
72674ae819SStefano Zampini   ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr);
73674ae819SStefano Zampini   PetscFunctionReturn(0);
74674ae819SStefano Zampini }
75674ae819SStefano Zampini 
76674ae819SStefano Zampini #undef __FUNCT__
776bfb1811SStefano Zampini #define __FUNCT__ "PCBDDCCreateWorkVectors"
786bfb1811SStefano Zampini PetscErrorCode PCBDDCCreateWorkVectors(PC pc)
796bfb1811SStefano Zampini {
806bfb1811SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
816bfb1811SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
826bfb1811SStefano Zampini   VecType        impVecType;
836bfb1811SStefano Zampini   PetscInt       n_vertices,n_constraints,local_primal_size,n_R;
846bfb1811SStefano Zampini   PetscErrorCode ierr;
856bfb1811SStefano Zampini 
866bfb1811SStefano Zampini   PetscFunctionBegin;
876bfb1811SStefano Zampini   ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,NULL);CHKERRQ(ierr);
886bfb1811SStefano Zampini   ierr = PCBDDCGetPrimalConstraintsLocalIdx(pc,&n_constraints,NULL);CHKERRQ(ierr);
896bfb1811SStefano Zampini   local_primal_size = n_constraints+n_vertices;
906bfb1811SStefano Zampini   n_R = pcis->n-n_vertices;
91*83b7ccabSStefano Zampini   /* parallel work vectors used in presolve. TODO: move outside */
926bfb1811SStefano Zampini   ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr);
936bfb1811SStefano Zampini   ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr);
946bfb1811SStefano Zampini   /* local work vectors */
956bfb1811SStefano Zampini   ierr = VecGetType(pcis->vec1_N,&impVecType);CHKERRQ(ierr);
966bfb1811SStefano Zampini   ierr = VecDuplicate(pcis->vec1_D,&pcbddc->vec4_D);CHKERRQ(ierr);
976bfb1811SStefano Zampini   ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_R);CHKERRQ(ierr);
986bfb1811SStefano Zampini   ierr = VecSetSizes(pcbddc->vec1_R,PETSC_DECIDE,n_R);CHKERRQ(ierr);
996bfb1811SStefano Zampini   ierr = VecSetType(pcbddc->vec1_R,impVecType);CHKERRQ(ierr);
1006bfb1811SStefano Zampini   ierr = VecDuplicate(pcbddc->vec1_R,&pcbddc->vec2_R);CHKERRQ(ierr);
101*83b7ccabSStefano Zampini   ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_P);CHKERRQ(ierr);
1026bfb1811SStefano Zampini   ierr = VecSetSizes(pcbddc->vec1_P,PETSC_DECIDE,local_primal_size);CHKERRQ(ierr);
1036bfb1811SStefano Zampini   ierr = VecSetType(pcbddc->vec1_P,impVecType);CHKERRQ(ierr);
104*83b7ccabSStefano Zampini   if (n_constraints) {
105*83b7ccabSStefano Zampini     ierr = VecCreate(PetscObjectComm((PetscObject)pcis->vec1_N),&pcbddc->vec1_C);CHKERRQ(ierr);
106*83b7ccabSStefano Zampini     ierr = VecSetSizes(pcbddc->vec1_C,PETSC_DECIDE,n_constraints);CHKERRQ(ierr);
107*83b7ccabSStefano Zampini     ierr = VecSetType(pcbddc->vec1_C,impVecType);CHKERRQ(ierr);
108*83b7ccabSStefano Zampini   }
1096bfb1811SStefano Zampini   PetscFunctionReturn(0);
1106bfb1811SStefano Zampini }
1116bfb1811SStefano Zampini 
1126bfb1811SStefano Zampini #undef __FUNCT__
113aa0d41d4SStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalMatrices"
114aa0d41d4SStefano Zampini PetscErrorCode PCBDDCSetUpLocalMatrices(PC pc)
115aa0d41d4SStefano Zampini {
116aa0d41d4SStefano Zampini   PC_IS*            pcis = (PC_IS*)(pc->data);
117aa0d41d4SStefano Zampini   PC_BDDC*          pcbddc = (PC_BDDC*)pc->data;
118aa0d41d4SStefano Zampini   Mat_IS*           matis = (Mat_IS*)pc->pmat->data;
119aa0d41d4SStefano Zampini   /* manage repeated solves */
120aa0d41d4SStefano Zampini   MatReuse          reuse;
121aa0d41d4SStefano Zampini   MatStructure      matstruct;
122aa0d41d4SStefano Zampini   PetscErrorCode    ierr;
123aa0d41d4SStefano Zampini 
124aa0d41d4SStefano Zampini   PetscFunctionBegin;
125aa0d41d4SStefano Zampini   /* get mat flags */
126aa0d41d4SStefano Zampini   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
127aa0d41d4SStefano Zampini   reuse = MAT_INITIAL_MATRIX;
128aa0d41d4SStefano Zampini   if (pc->setupcalled) {
129aa0d41d4SStefano Zampini     /* when matstruct is SAME_PRECONDITIONER, we shouldn't be here */
130aa0d41d4SStefano Zampini     if (matstruct == SAME_PRECONDITIONER) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"This should not happen");
131aa0d41d4SStefano Zampini     if (matstruct == SAME_NONZERO_PATTERN) {
132aa0d41d4SStefano Zampini       reuse = MAT_REUSE_MATRIX;
133aa0d41d4SStefano Zampini     } else {
134aa0d41d4SStefano Zampini       reuse = MAT_INITIAL_MATRIX;
135aa0d41d4SStefano Zampini     }
136aa0d41d4SStefano Zampini   }
137aa0d41d4SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
138aa0d41d4SStefano Zampini     ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr);
139aa0d41d4SStefano Zampini     ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr);
140aa0d41d4SStefano Zampini     ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr);
141aa0d41d4SStefano Zampini     ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr);
142aa0d41d4SStefano Zampini     ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
143aa0d41d4SStefano Zampini   }
144aa0d41d4SStefano Zampini 
145aa0d41d4SStefano Zampini   /* transform local matrices if needed */
146aa0d41d4SStefano Zampini   if (pcbddc->use_change_of_basis) {
147aa0d41d4SStefano Zampini     Mat         change_mat_all;
148aa0d41d4SStefano Zampini     PetscScalar *row_cmat_values;
149aa0d41d4SStefano Zampini     PetscInt    *row_cmat_indices;
150aa0d41d4SStefano Zampini     PetscInt    *nnz,*is_indices,*temp_indices;
151aa0d41d4SStefano Zampini     PetscInt    i,j,k,n_D,n_B;
152aa0d41d4SStefano Zampini 
153aa0d41d4SStefano Zampini     /* Get Non-overlapping dimensions */
154aa0d41d4SStefano Zampini     n_B = pcis->n_B;
155aa0d41d4SStefano Zampini     n_D = pcis->n-n_B;
156aa0d41d4SStefano Zampini 
157aa0d41d4SStefano Zampini     /* compute nonzero structure of change of basis on all local nodes */
158aa0d41d4SStefano Zampini     ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
159aa0d41d4SStefano Zampini     ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
160aa0d41d4SStefano Zampini     for (i=0;i<n_D;i++) nnz[is_indices[i]] = 1;
161aa0d41d4SStefano Zampini     ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
162aa0d41d4SStefano Zampini     ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
163aa0d41d4SStefano Zampini     k=1;
164aa0d41d4SStefano Zampini     for (i=0;i<n_B;i++) {
165aa0d41d4SStefano Zampini       ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr);
166aa0d41d4SStefano Zampini       nnz[is_indices[i]]=j;
167aa0d41d4SStefano Zampini       if (k < j) k = j;
168aa0d41d4SStefano Zampini       ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,NULL,NULL);CHKERRQ(ierr);
169aa0d41d4SStefano Zampini     }
170aa0d41d4SStefano Zampini     ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
171aa0d41d4SStefano Zampini     /* assemble change of basis matrix on the whole set of local dofs */
172aa0d41d4SStefano Zampini     ierr = PetscMalloc(k*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
173aa0d41d4SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&change_mat_all);CHKERRQ(ierr);
174aa0d41d4SStefano Zampini     ierr = MatSetSizes(change_mat_all,pcis->n,pcis->n,pcis->n,pcis->n);CHKERRQ(ierr);
175aa0d41d4SStefano Zampini     ierr = MatSetType(change_mat_all,MATSEQAIJ);CHKERRQ(ierr);
176aa0d41d4SStefano Zampini     ierr = MatSeqAIJSetPreallocation(change_mat_all,0,nnz);CHKERRQ(ierr);
177aa0d41d4SStefano Zampini     ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
178aa0d41d4SStefano Zampini     for (i=0;i<n_D;i++) {
179aa0d41d4SStefano Zampini       ierr = MatSetValue(change_mat_all,is_indices[i],is_indices[i],1.0,INSERT_VALUES);CHKERRQ(ierr);
180aa0d41d4SStefano Zampini     }
181aa0d41d4SStefano Zampini     ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
182aa0d41d4SStefano Zampini     ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
183aa0d41d4SStefano Zampini     for (i=0;i<n_B;i++) {
184aa0d41d4SStefano Zampini       ierr = MatGetRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr);
185aa0d41d4SStefano Zampini       for (k=0; k<j; k++) temp_indices[k]=is_indices[row_cmat_indices[k]];
186aa0d41d4SStefano Zampini       ierr = MatSetValues(change_mat_all,1,&is_indices[i],j,temp_indices,row_cmat_values,INSERT_VALUES);CHKERRQ(ierr);
187aa0d41d4SStefano Zampini       ierr = MatRestoreRow(pcbddc->ChangeOfBasisMatrix,i,&j,(const PetscInt**)&row_cmat_indices,(const PetscScalar**)&row_cmat_values);CHKERRQ(ierr);
188aa0d41d4SStefano Zampini     }
189aa0d41d4SStefano Zampini     ierr = MatAssemblyBegin(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
190aa0d41d4SStefano Zampini     ierr = MatAssemblyEnd(change_mat_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191aa0d41d4SStefano Zampini     /* TODO: HOW TO WORK WITH BAIJ? PtAP not provided */
192aa0d41d4SStefano Zampini     ierr = MatGetBlockSize(matis->A,&i);CHKERRQ(ierr);
193aa0d41d4SStefano Zampini     if (i==1) {
194aa0d41d4SStefano Zampini       ierr = MatPtAP(matis->A,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
195aa0d41d4SStefano Zampini     } else {
196aa0d41d4SStefano Zampini       Mat work_mat;
197aa0d41d4SStefano Zampini       ierr = MatConvert(matis->A,MATSEQAIJ,MAT_INITIAL_MATRIX,&work_mat);CHKERRQ(ierr);
198aa0d41d4SStefano Zampini       ierr = MatPtAP(work_mat,change_mat_all,reuse,2.0,&pcbddc->local_mat);CHKERRQ(ierr);
199aa0d41d4SStefano Zampini       ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
200aa0d41d4SStefano Zampini     }
201aa0d41d4SStefano Zampini     ierr = MatDestroy(&change_mat_all);CHKERRQ(ierr);
202aa0d41d4SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
203aa0d41d4SStefano Zampini     ierr = PetscFree(temp_indices);CHKERRQ(ierr);
204aa0d41d4SStefano Zampini   } else {
205aa0d41d4SStefano Zampini     /* without change of basis, the local matrix is unchanged */
206aa0d41d4SStefano Zampini     if (!pcbddc->local_mat) {
207aa0d41d4SStefano Zampini       ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
208aa0d41d4SStefano Zampini       pcbddc->local_mat = matis->A;
209aa0d41d4SStefano Zampini     }
210aa0d41d4SStefano Zampini   }
211aa0d41d4SStefano Zampini 
212aa0d41d4SStefano Zampini   /* get submatrices */
213aa0d41d4SStefano Zampini   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr);
214aa0d41d4SStefano Zampini   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr);
215aa0d41d4SStefano Zampini   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr);
216aa0d41d4SStefano Zampini   ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr);
217aa0d41d4SStefano Zampini   PetscFunctionReturn(0);
218aa0d41d4SStefano Zampini }
219aa0d41d4SStefano Zampini 
220aa0d41d4SStefano Zampini #undef __FUNCT__
221a64d13efSStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalScatters"
222a64d13efSStefano Zampini PetscErrorCode PCBDDCSetUpLocalScatters(PC pc,IS* is_R_local_n)
223a64d13efSStefano Zampini {
224a64d13efSStefano Zampini   PC_IS*         pcis = (PC_IS*)(pc->data);
225a64d13efSStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
226a64d13efSStefano Zampini   IS             is_R_local,is_aux1,is_aux2;
227a64d13efSStefano Zampini   PetscInt       *vertices,*aux_array1,*aux_array2,*is_indices,*idx_R_local;
228a64d13efSStefano Zampini   PetscInt       n_vertices,n_constraints,i,j,n_R,n_D,n_B;
229a64d13efSStefano Zampini   PetscBool      *array_bool;
230a64d13efSStefano Zampini   PetscErrorCode ierr;
231a64d13efSStefano Zampini 
232a64d13efSStefano Zampini   PetscFunctionBegin;
233a64d13efSStefano Zampini   /* Set Non-overlapping dimensions */
234a64d13efSStefano Zampini   n_B = pcis->n_B; n_D = pcis->n - n_B;
235a64d13efSStefano Zampini   /* get vertex indices from constraint matrix */
236a64d13efSStefano Zampini   ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&n_vertices,&vertices);CHKERRQ(ierr);
237a64d13efSStefano Zampini   /* Set number of constraints */
238a64d13efSStefano Zampini   n_constraints = pcbddc->local_primal_size-n_vertices;
239a64d13efSStefano Zampini   /* Dohrmann's notation: dofs splitted in R (Remaining: all dofs but the vertices) and V (Vertices) */
240a64d13efSStefano Zampini   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&array_bool);CHKERRQ(ierr);
241a64d13efSStefano Zampini   for (i=0;i<pcis->n;i++) array_bool[i] = PETSC_TRUE;
242a64d13efSStefano Zampini   for (i=0;i<n_vertices;i++) array_bool[vertices[i]] = PETSC_FALSE;
243a64d13efSStefano Zampini   ierr = PetscMalloc((pcis->n-n_vertices)*sizeof(PetscInt),&idx_R_local);CHKERRQ(ierr);
244a64d13efSStefano Zampini   for (i=0, n_R=0; i<pcis->n; i++) {
245a64d13efSStefano Zampini     if (array_bool[i]) {
246a64d13efSStefano Zampini       idx_R_local[n_R] = i;
247a64d13efSStefano Zampini       n_R++;
248a64d13efSStefano Zampini     }
249a64d13efSStefano Zampini   }
250a64d13efSStefano Zampini   ierr = PetscFree(vertices);CHKERRQ(ierr);
251a64d13efSStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_R,idx_R_local,PETSC_OWN_POINTER,&is_R_local);CHKERRQ(ierr);
252a64d13efSStefano Zampini 
253a64d13efSStefano Zampini   /* print some info if requested */
254a64d13efSStefano Zampini   if (pcbddc->dbg_flag) {
255a64d13efSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
256a64d13efSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
257a64d13efSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d local dimensions\n",PetscGlobalRank);CHKERRQ(ierr);
258a64d13efSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"local_size = %d, dirichlet_size = %d, boundary_size = %d\n",pcis->n,n_D,n_B);CHKERRQ(ierr);
259a64d13efSStefano Zampini     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);
260a64d13efSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"pcbddc->n_vertices = %d, pcbddc->n_constraints = %d\n",pcbddc->n_vertices,pcbddc->n_constraints);CHKERRQ(ierr);
261a64d13efSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
262a64d13efSStefano Zampini   }
263a64d13efSStefano Zampini 
264a64d13efSStefano Zampini   /* VecScatters pcbddc->R_to_B and (optionally) pcbddc->R_to_D */
265a64d13efSStefano Zampini   ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr);
266a64d13efSStefano Zampini   ierr = PetscMalloc((pcis->n_B-n_vertices)*sizeof(PetscInt),&aux_array2);CHKERRQ(ierr);
267a64d13efSStefano Zampini   ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
268a64d13efSStefano Zampini   for (i=0; i<n_D; i++) array_bool[is_indices[i]] = PETSC_FALSE;
269a64d13efSStefano Zampini   ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
270a64d13efSStefano Zampini   for (i=0, j=0; i<n_R; i++) {
271a64d13efSStefano Zampini     if (array_bool[idx_R_local[i]]) {
272a64d13efSStefano Zampini       aux_array1[j] = i;
273a64d13efSStefano Zampini       j++;
274a64d13efSStefano Zampini     }
275a64d13efSStefano Zampini   }
276a64d13efSStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr);
277a64d13efSStefano Zampini   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
278a64d13efSStefano Zampini   for (i=0, j=0; i<n_B; i++) {
279a64d13efSStefano Zampini     if (array_bool[is_indices[i]]) {
280a64d13efSStefano Zampini       aux_array2[j] = i; j++;
281a64d13efSStefano Zampini     }
282a64d13efSStefano Zampini   }
283a64d13efSStefano Zampini   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
284a64d13efSStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array2,PETSC_OWN_POINTER,&is_aux2);CHKERRQ(ierr);
285a64d13efSStefano Zampini   ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_B,is_aux2,&pcbddc->R_to_B);CHKERRQ(ierr);
286a64d13efSStefano Zampini   ierr = ISDestroy(&is_aux1);CHKERRQ(ierr);
287a64d13efSStefano Zampini   ierr = ISDestroy(&is_aux2);CHKERRQ(ierr);
288a64d13efSStefano Zampini 
289a64d13efSStefano Zampini   if (pcbddc->inexact_prec_type || pcbddc->dbg_flag ) {
290a64d13efSStefano Zampini     ierr = PetscMalloc(n_D*sizeof(PetscInt),&aux_array1);CHKERRQ(ierr);
291a64d13efSStefano Zampini     for (i=0, j=0; i<n_R; i++) {
292a64d13efSStefano Zampini       if (!array_bool[idx_R_local[i]]) {
293a64d13efSStefano Zampini         aux_array1[j] = i;
294a64d13efSStefano Zampini         j++;
295a64d13efSStefano Zampini       }
296a64d13efSStefano Zampini     }
297a64d13efSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,j,aux_array1,PETSC_OWN_POINTER,&is_aux1);CHKERRQ(ierr);
298a64d13efSStefano Zampini     ierr = VecScatterCreate(pcbddc->vec1_R,is_aux1,pcis->vec1_D,(IS)0,&pcbddc->R_to_D);CHKERRQ(ierr);
299a64d13efSStefano Zampini     ierr = ISDestroy(&is_aux1);CHKERRQ(ierr);
300a64d13efSStefano Zampini   }
301a64d13efSStefano Zampini   ierr = PetscFree(array_bool);CHKERRQ(ierr);
302a64d13efSStefano Zampini   *is_R_local_n = is_R_local;
303a64d13efSStefano Zampini   PetscFunctionReturn(0);
304a64d13efSStefano Zampini }
305a64d13efSStefano Zampini 
306a64d13efSStefano Zampini #undef __FUNCT__
307304d26faSStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
308304d26faSStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool use)
309304d26faSStefano Zampini {
310304d26faSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
311304d26faSStefano Zampini 
312304d26faSStefano Zampini   PetscFunctionBegin;
313304d26faSStefano Zampini   pcbddc->use_exact_dirichlet=use;
314304d26faSStefano Zampini   PetscFunctionReturn(0);
315304d26faSStefano Zampini }
316304d26faSStefano Zampini 
317304d26faSStefano Zampini #undef __FUNCT__
318304d26faSStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalSolvers"
319304d26faSStefano Zampini PetscErrorCode PCBDDCSetUpLocalSolvers(PC pc, IS is_I_local, IS is_R_local)
320304d26faSStefano Zampini {
321304d26faSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
322304d26faSStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
323304d26faSStefano Zampini   PC             pc_temp;
324304d26faSStefano Zampini   Mat            A_RR;
325304d26faSStefano Zampini   Vec            vec1,vec2,vec3;
326304d26faSStefano Zampini   MatStructure   matstruct;
327304d26faSStefano Zampini   PetscScalar    m_one = -1.0;
328304d26faSStefano Zampini   PetscReal      value;
329304d26faSStefano Zampini   PetscInt       n_D,n_R,use_exact,use_exact_reduced;
330304d26faSStefano Zampini   PetscErrorCode ierr;
331304d26faSStefano Zampini 
332304d26faSStefano Zampini   PetscFunctionBegin;
333304d26faSStefano Zampini   /* Creating PC contexts for local Dirichlet and Neumann problems */
334304d26faSStefano Zampini   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
335304d26faSStefano Zampini 
336304d26faSStefano Zampini   /* DIRICHLET PROBLEM */
337ac78edfcSStefano Zampini   /* Matrix for Dirichlet problem is pcis->A_II */
338304d26faSStefano Zampini   ierr = ISGetSize(is_I_local,&n_D);CHKERRQ(ierr);
339304d26faSStefano Zampini   if (!pcbddc->ksp_D) { /* create object if not yet build */
340304d26faSStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_D);CHKERRQ(ierr);
341304d26faSStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
342304d26faSStefano Zampini     /* default */
343304d26faSStefano Zampini     ierr = KSPSetType(pcbddc->ksp_D,KSPPREONLY);CHKERRQ(ierr);
344304d26faSStefano Zampini     ierr = KSPSetOptionsPrefix(pcbddc->ksp_D,"dirichlet_");CHKERRQ(ierr);
345304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
346304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
347304d26faSStefano Zampini     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
348304d26faSStefano Zampini   }
349304d26faSStefano Zampini   ierr = KSPSetOperators(pcbddc->ksp_D,pcis->A_II,pcis->A_II,matstruct);CHKERRQ(ierr);
350304d26faSStefano Zampini   /* Allow user's customization */
351304d26faSStefano Zampini   ierr = KSPSetFromOptions(pcbddc->ksp_D);CHKERRQ(ierr);
352304d26faSStefano Zampini   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
353304d26faSStefano Zampini   if (!n_D) {
354304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
355304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
356304d26faSStefano Zampini   }
357304d26faSStefano Zampini   /* Set Up KSP for Dirichlet problem of BDDC */
358304d26faSStefano Zampini   ierr = KSPSetUp(pcbddc->ksp_D);CHKERRQ(ierr);
359304d26faSStefano Zampini   /* set ksp_D into pcis data */
360304d26faSStefano Zampini   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
361304d26faSStefano Zampini   ierr = PetscObjectReference((PetscObject)pcbddc->ksp_D);CHKERRQ(ierr);
362304d26faSStefano Zampini   pcis->ksp_D = pcbddc->ksp_D;
363304d26faSStefano Zampini 
364304d26faSStefano Zampini   /* NEUMANN PROBLEM */
365304d26faSStefano Zampini   /* Matrix for Neumann problem is A_RR -> we need to create it */
366304d26faSStefano Zampini   ierr = ISGetSize(is_R_local,&n_R);CHKERRQ(ierr);
367304d26faSStefano Zampini   ierr = MatGetSubMatrix(pcbddc->local_mat,is_R_local,is_R_local,MAT_INITIAL_MATRIX,&A_RR);CHKERRQ(ierr);
368304d26faSStefano Zampini   if (!pcbddc->ksp_R) { /* create object if not yet build */
369304d26faSStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_R);CHKERRQ(ierr);
370304d26faSStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_R,(PetscObject)pc,1);CHKERRQ(ierr);
371304d26faSStefano Zampini     /* default */
372304d26faSStefano Zampini     ierr = KSPSetType(pcbddc->ksp_R,KSPPREONLY);CHKERRQ(ierr);
373304d26faSStefano Zampini     ierr = KSPSetOptionsPrefix(pcbddc->ksp_R,"neumann_");CHKERRQ(ierr);
374304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
375304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
376304d26faSStefano Zampini     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
377304d26faSStefano Zampini   }
378304d26faSStefano Zampini   ierr = KSPSetOperators(pcbddc->ksp_R,A_RR,A_RR,matstruct);CHKERRQ(ierr);
379304d26faSStefano Zampini   /* Allow user's customization */
380304d26faSStefano Zampini   ierr = KSPSetFromOptions(pcbddc->ksp_R);CHKERRQ(ierr);
381304d26faSStefano Zampini   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
382304d26faSStefano Zampini   if (!n_R) {
383304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
384304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
385304d26faSStefano Zampini   }
386304d26faSStefano Zampini   /* Set Up KSP for Neumann problem of BDDC */
387304d26faSStefano Zampini   ierr = KSPSetUp(pcbddc->ksp_R);CHKERRQ(ierr);
388304d26faSStefano Zampini 
389304d26faSStefano Zampini   /* check Dirichlet and Neumann solvers and adapt them if a nullspace correction is needed */
390304d26faSStefano Zampini 
391304d26faSStefano Zampini   /* Dirichlet */
392304d26faSStefano Zampini   ierr = MatGetVecs(pcis->A_II,&vec1,&vec2);CHKERRQ(ierr);
393304d26faSStefano Zampini   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
394304d26faSStefano Zampini   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
395304d26faSStefano Zampini   ierr = MatMult(pcis->A_II,vec1,vec2);CHKERRQ(ierr);
396304d26faSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,vec2,vec3);CHKERRQ(ierr);
397304d26faSStefano Zampini   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
398304d26faSStefano Zampini   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
399304d26faSStefano Zampini   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
400304d26faSStefano Zampini   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
401304d26faSStefano Zampini   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
402304d26faSStefano Zampini   /* need to be adapted? */
403304d26faSStefano Zampini   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
404304d26faSStefano Zampini   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
405304d26faSStefano Zampini   ierr = PCBDDCSetUseExactDirichlet(pc,(PetscBool)use_exact_reduced);CHKERRQ(ierr);
406304d26faSStefano Zampini   /* print info */
407304d26faSStefano Zampini   if (pcbddc->dbg_flag) {
408304d26faSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
409304d26faSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
410304d26faSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Checking solution of Dirichlet and Neumann problems\n");CHKERRQ(ierr);
411304d26faSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
412304d26faSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
413304d26faSStefano Zampini   }
414304d26faSStefano Zampini   if (n_D && pcbddc->NullSpace && !use_exact_reduced && !pcbddc->inexact_prec_type) {
415304d26faSStefano Zampini     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_I_local);CHKERRQ(ierr);
416304d26faSStefano Zampini   }
417304d26faSStefano Zampini 
418304d26faSStefano Zampini   /* Neumann */
419304d26faSStefano Zampini   ierr = MatGetVecs(A_RR,&vec1,&vec2);CHKERRQ(ierr);
420304d26faSStefano Zampini   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
421304d26faSStefano Zampini   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
422304d26faSStefano Zampini   ierr = MatMult(A_RR,vec1,vec2);CHKERRQ(ierr);
423304d26faSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_R,vec2,vec3);CHKERRQ(ierr);
424304d26faSStefano Zampini   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
425304d26faSStefano Zampini   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
426304d26faSStefano Zampini   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
427304d26faSStefano Zampini   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
428304d26faSStefano Zampini   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
429304d26faSStefano Zampini   /* need to be adapted? */
430304d26faSStefano Zampini   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
431304d26faSStefano Zampini   if (PetscAbsReal(value) > 1.e-4) use_exact = 0;
432304d26faSStefano Zampini   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
433304d26faSStefano Zampini   /* print info */
434304d26faSStefano Zampini   if (pcbddc->dbg_flag) {
435304d26faSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for  Neumann  solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
436304d26faSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
437304d26faSStefano Zampini   }
438304d26faSStefano Zampini   if (n_R && pcbddc->NullSpace && !use_exact_reduced) { /* is it the right logic? */
439304d26faSStefano Zampini     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_R_local);CHKERRQ(ierr);
440304d26faSStefano Zampini   }
441304d26faSStefano Zampini 
442304d26faSStefano Zampini   /* free Neumann problem's matrix */
443304d26faSStefano Zampini   ierr = MatDestroy(&A_RR);CHKERRQ(ierr);
444304d26faSStefano Zampini   PetscFunctionReturn(0);
445304d26faSStefano Zampini }
446304d26faSStefano Zampini 
447304d26faSStefano Zampini #undef __FUNCT__
448674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSolveSaddlePoint"
449674ae819SStefano Zampini static PetscErrorCode  PCBDDCSolveSaddlePoint(PC pc)
450674ae819SStefano Zampini {
451674ae819SStefano Zampini   PetscErrorCode ierr;
452674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
453674ae819SStefano Zampini 
454674ae819SStefano Zampini   PetscFunctionBegin;
455674ae819SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr);
456674ae819SStefano Zampini   if (pcbddc->local_auxmat1) {
457674ae819SStefano Zampini     ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr);
458674ae819SStefano Zampini     ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr);
459674ae819SStefano Zampini   }
460674ae819SStefano Zampini   PetscFunctionReturn(0);
461674ae819SStefano Zampini }
462674ae819SStefano Zampini 
463674ae819SStefano Zampini #undef __FUNCT__
464674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner"
465674ae819SStefano Zampini PetscErrorCode  PCBDDCApplyInterfacePreconditioner(PC pc)
466674ae819SStefano Zampini {
467674ae819SStefano Zampini   PetscErrorCode ierr;
468674ae819SStefano Zampini   PC_BDDC*        pcbddc = (PC_BDDC*)(pc->data);
469674ae819SStefano Zampini   PC_IS*            pcis = (PC_IS*)  (pc->data);
470674ae819SStefano Zampini   const PetscScalar zero = 0.0;
471674ae819SStefano Zampini 
472674ae819SStefano Zampini   PetscFunctionBegin;
47315aaf578SStefano Zampini   /* Application of PHI^T (or PSI^T)  */
47415aaf578SStefano Zampini   if (pcbddc->coarse_psi_B) {
47515aaf578SStefano Zampini     ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
47615aaf578SStefano Zampini     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
47715aaf578SStefano Zampini   } else {
478674ae819SStefano Zampini     ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
479674ae819SStefano Zampini     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
48015aaf578SStefano Zampini   }
481674ae819SStefano Zampini   /* Scatter data of coarse_rhs */
482674ae819SStefano Zampini   if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); }
483674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
484674ae819SStefano Zampini 
485674ae819SStefano Zampini   /* Local solution on R nodes */
486674ae819SStefano Zampini   ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr);
487674ae819SStefano Zampini   ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
488674ae819SStefano Zampini   ierr = VecScatterEnd  (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
489674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) {
490674ae819SStefano Zampini     ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
491674ae819SStefano Zampini     ierr = VecScatterEnd  (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
492674ae819SStefano Zampini   }
493674ae819SStefano Zampini   ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr);
494674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
495674ae819SStefano Zampini   ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
496674ae819SStefano Zampini   ierr = VecScatterEnd  (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
497674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) {
498674ae819SStefano Zampini     ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
499674ae819SStefano Zampini     ierr = VecScatterEnd  (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
500674ae819SStefano Zampini   }
501674ae819SStefano Zampini 
502674ae819SStefano Zampini   /* Coarse solution */
503674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
504674ae819SStefano Zampini   if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */
505674ae819SStefano Zampini     ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr);
506674ae819SStefano Zampini   }
507674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
508674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataEnd  (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
509674ae819SStefano Zampini 
510674ae819SStefano Zampini   /* Sum contributions from two levels */
511674ae819SStefano Zampini   ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
512674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
513674ae819SStefano Zampini   PetscFunctionReturn(0);
514674ae819SStefano Zampini }
515674ae819SStefano Zampini 
516674ae819SStefano Zampini #undef __FUNCT__
517674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataBegin"
518674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
519674ae819SStefano Zampini {
520674ae819SStefano Zampini   PetscErrorCode ierr;
521674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
522674ae819SStefano Zampini 
523674ae819SStefano Zampini   PetscFunctionBegin;
524674ae819SStefano Zampini   switch (pcbddc->coarse_communications_type) {
525674ae819SStefano Zampini     case SCATTERS_BDDC:
526674ae819SStefano Zampini       ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
527674ae819SStefano Zampini       break;
528674ae819SStefano Zampini     case GATHERS_BDDC:
529674ae819SStefano Zampini       break;
530674ae819SStefano Zampini   }
531674ae819SStefano Zampini   PetscFunctionReturn(0);
532674ae819SStefano Zampini }
533674ae819SStefano Zampini 
534674ae819SStefano Zampini #undef __FUNCT__
535674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataEnd"
536674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
537674ae819SStefano Zampini {
538674ae819SStefano Zampini   PetscErrorCode ierr;
539674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
540674ae819SStefano Zampini   PetscScalar*   array_to;
541674ae819SStefano Zampini   PetscScalar*   array_from;
542674ae819SStefano Zampini   MPI_Comm       comm;
543674ae819SStefano Zampini   PetscInt       i;
544674ae819SStefano Zampini 
545674ae819SStefano Zampini   PetscFunctionBegin;
546674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
547674ae819SStefano Zampini   switch (pcbddc->coarse_communications_type) {
548674ae819SStefano Zampini     case SCATTERS_BDDC:
549674ae819SStefano Zampini       ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
550674ae819SStefano Zampini       break;
551674ae819SStefano Zampini     case GATHERS_BDDC:
552674ae819SStefano Zampini       if (vec_from) {
553674ae819SStefano Zampini         ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr);
554674ae819SStefano Zampini       }
555674ae819SStefano Zampini       if (vec_to) {
556674ae819SStefano Zampini         ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr);
557674ae819SStefano Zampini       }
558674ae819SStefano Zampini       switch(pcbddc->coarse_problem_type){
559674ae819SStefano Zampini         case SEQUENTIAL_BDDC:
560674ae819SStefano Zampini           if (smode == SCATTER_FORWARD) {
561674ae819SStefano Zampini             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);
562674ae819SStefano Zampini             if (vec_to) {
563674ae819SStefano Zampini               if (imode == ADD_VALUES) {
564674ae819SStefano Zampini                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
565674ae819SStefano Zampini                   array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
566674ae819SStefano Zampini                 }
567674ae819SStefano Zampini               } else {
568674ae819SStefano Zampini                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
569674ae819SStefano Zampini                   array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
570674ae819SStefano Zampini                 }
571674ae819SStefano Zampini               }
572674ae819SStefano Zampini             }
573674ae819SStefano Zampini           } else {
574674ae819SStefano Zampini             if (vec_from) {
575674ae819SStefano Zampini               if (imode == ADD_VALUES) {
576674ae819SStefano Zampini                 MPI_Comm vec_from_comm;
577674ae819SStefano Zampini                 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr);
578674ae819SStefano Zampini                 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);
579674ae819SStefano Zampini               }
580674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
581674ae819SStefano Zampini                 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]];
582674ae819SStefano Zampini               }
583674ae819SStefano Zampini             }
584674ae819SStefano Zampini             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);
585674ae819SStefano Zampini           }
586674ae819SStefano Zampini           break;
587674ae819SStefano Zampini         case REPLICATED_BDDC:
588674ae819SStefano Zampini           if (smode == SCATTER_FORWARD) {
589674ae819SStefano Zampini             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);
590674ae819SStefano Zampini             if (imode == ADD_VALUES) {
591674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
592674ae819SStefano Zampini                 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
593674ae819SStefano Zampini               }
594674ae819SStefano Zampini             } else {
595674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
596674ae819SStefano Zampini                 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
597674ae819SStefano Zampini               }
598674ae819SStefano Zampini             }
599674ae819SStefano Zampini           } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */
600674ae819SStefano Zampini             if (imode == ADD_VALUES) {
601674ae819SStefano Zampini               for (i=0;i<pcbddc->local_primal_size;i++) {
602674ae819SStefano Zampini                 array_to[i]+=array_from[pcbddc->local_primal_indices[i]];
603674ae819SStefano Zampini               }
604674ae819SStefano Zampini             } else {
605674ae819SStefano Zampini               for (i=0;i<pcbddc->local_primal_size;i++) {
606674ae819SStefano Zampini                 array_to[i]=array_from[pcbddc->local_primal_indices[i]];
607674ae819SStefano Zampini               }
608674ae819SStefano Zampini             }
609674ae819SStefano Zampini           }
610674ae819SStefano Zampini           break;
611674ae819SStefano Zampini         case MULTILEVEL_BDDC:
612674ae819SStefano Zampini           break;
613674ae819SStefano Zampini         case PARALLEL_BDDC:
614674ae819SStefano Zampini           break;
615674ae819SStefano Zampini       }
616674ae819SStefano Zampini       if (vec_from) {
617674ae819SStefano Zampini         ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr);
618674ae819SStefano Zampini       }
619674ae819SStefano Zampini       if (vec_to) {
620674ae819SStefano Zampini         ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr);
621674ae819SStefano Zampini       }
622674ae819SStefano Zampini       break;
623674ae819SStefano Zampini   }
624674ae819SStefano Zampini   PetscFunctionReturn(0);
625674ae819SStefano Zampini }
626674ae819SStefano Zampini 
627674ae819SStefano Zampini #undef __FUNCT__
628674ae819SStefano Zampini #define __FUNCT__ "PCBDDCConstraintsSetUp"
629674ae819SStefano Zampini PetscErrorCode PCBDDCConstraintsSetUp(PC pc)
630674ae819SStefano Zampini {
631674ae819SStefano Zampini   PetscErrorCode ierr;
632674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)(pc->data);
633674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
634674ae819SStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
635674ae819SStefano Zampini   PetscInt       *nnz,*is_indices;
636674ae819SStefano Zampini   PetscScalar    *temp_quadrature_constraint;
637674ae819SStefano Zampini   PetscInt       *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B;
638674ae819SStefano Zampini   PetscInt       local_primal_size,i,j,k,total_counts,max_size_of_constraint;
639674ae819SStefano Zampini   PetscInt       n_vertices,size_of_constraint;
6405b08dc53SStefano Zampini   PetscReal      real_value;
641674ae819SStefano Zampini   PetscBool      nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true;
642674ae819SStefano Zampini   PetscInt       nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges;
643674ae819SStefano Zampini   IS             *used_IS,ISForVertices,*ISForFaces,*ISForEdges;
644674ae819SStefano Zampini   MatType        impMatType=MATSEQAIJ;
645674ae819SStefano Zampini   PetscBLASInt   Bs,Bt,lwork,lierr;
646674ae819SStefano Zampini   PetscReal      tol=1.0e-8;
647674ae819SStefano Zampini   MatNullSpace   nearnullsp;
648674ae819SStefano Zampini   const Vec      *nearnullvecs;
649674ae819SStefano Zampini   Vec            *localnearnullsp;
650674ae819SStefano Zampini   PetscScalar    *work,*temp_basis,*array_vector,*correlation_mat;
651674ae819SStefano Zampini   PetscReal      *rwork,*singular_vals;
652674ae819SStefano Zampini   PetscBLASInt   Bone=1,*ipiv;
653674ae819SStefano Zampini   Vec            temp_vec;
654674ae819SStefano Zampini   Mat            temp_mat;
655674ae819SStefano Zampini   KSP            temp_ksp;
656674ae819SStefano Zampini   PC             temp_pc;
657674ae819SStefano Zampini   PetscInt       s,start_constraint,dual_dofs;
658674ae819SStefano Zampini   PetscBool      compute_submatrix,useksp=PETSC_FALSE;
659674ae819SStefano Zampini   PetscInt       *aux_primal_permutation,*aux_primal_numbering;
66051b0f6cfSStefano Zampini   PetscBool      boolforchange,*change_basis;
661674ae819SStefano Zampini /* some ugly conditional declarations */
662674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
663674ae819SStefano Zampini   PetscScalar    one=1.0,zero=0.0;
664674ae819SStefano Zampini   PetscInt       ii;
665674ae819SStefano Zampini   PetscScalar    *singular_vectors;
666674ae819SStefano Zampini   PetscBLASInt   *iwork,*ifail;
667674ae819SStefano Zampini   PetscReal      dummy_real,abs_tol;
668674ae819SStefano Zampini   PetscBLASInt   eigs_found;
669674ae819SStefano Zampini #endif
670674ae819SStefano Zampini   PetscBLASInt   dummy_int;
671674ae819SStefano Zampini   PetscScalar    dummy_scalar;
672674ae819SStefano Zampini   PetscBool      used_vertex,get_faces,get_edges,get_vertices;
673674ae819SStefano Zampini 
674674ae819SStefano Zampini   PetscFunctionBegin;
675674ae819SStefano Zampini   /* Get index sets for faces, edges and vertices from graph */
676674ae819SStefano Zampini   get_faces = PETSC_TRUE;
677674ae819SStefano Zampini   get_edges = PETSC_TRUE;
678674ae819SStefano Zampini   get_vertices = PETSC_TRUE;
679674ae819SStefano Zampini   if (pcbddc->vertices_flag) {
680674ae819SStefano Zampini     get_faces = PETSC_FALSE;
681674ae819SStefano Zampini     get_edges = PETSC_FALSE;
682674ae819SStefano Zampini   }
683674ae819SStefano Zampini   if (pcbddc->constraints_flag) {
684674ae819SStefano Zampini     get_vertices = PETSC_FALSE;
685674ae819SStefano Zampini   }
686674ae819SStefano Zampini   if (pcbddc->faces_flag) {
687674ae819SStefano Zampini     get_edges = PETSC_FALSE;
688674ae819SStefano Zampini   }
689674ae819SStefano Zampini   if (pcbddc->edges_flag) {
690674ae819SStefano Zampini     get_faces = PETSC_FALSE;
691674ae819SStefano Zampini   }
692674ae819SStefano Zampini   /* default */
693674ae819SStefano Zampini   if (!get_faces && !get_edges && !get_vertices) {
694674ae819SStefano Zampini     get_vertices = PETSC_TRUE;
695674ae819SStefano Zampini   }
696674ae819SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
697674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
698674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
699674ae819SStefano Zampini     i = 0;
700674ae819SStefano Zampini     if (ISForVertices) {
701674ae819SStefano Zampini       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
702674ae819SStefano Zampini     }
703674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
704674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
70515aaf578SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
706674ae819SStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
707674ae819SStefano Zampini   }
708674ae819SStefano Zampini   /* check if near null space is attached to global mat */
709674ae819SStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
710674ae819SStefano Zampini   if (nearnullsp) {
711674ae819SStefano Zampini     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
712674ae819SStefano Zampini   } else { /* if near null space is not provided it uses constants */
713674ae819SStefano Zampini     nnsp_has_cnst = PETSC_TRUE;
714674ae819SStefano Zampini     use_nnsp_true = PETSC_TRUE;
715674ae819SStefano Zampini   }
716674ae819SStefano Zampini   if (nnsp_has_cnst) {
717674ae819SStefano Zampini     nnsp_addone = 1;
718674ae819SStefano Zampini   }
719674ae819SStefano Zampini   /*
720674ae819SStefano Zampini        Evaluate maximum storage size needed by the procedure
721674ae819SStefano Zampini        - temp_indices will contain start index of each constraint stored as follows
722674ae819SStefano Zampini        - temp_indices_to_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the indices (in local numbering) on which the constraint acts
723674ae819SStefano Zampini        - 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
724674ae819SStefano Zampini        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
725674ae819SStefano Zampini                                                                                                                                                          */
726674ae819SStefano Zampini   total_counts = n_ISForFaces+n_ISForEdges;
727674ae819SStefano Zampini   total_counts *= (nnsp_addone+nnsp_size);
728674ae819SStefano Zampini   n_vertices = 0;
729674ae819SStefano Zampini   if (ISForVertices) {
730674ae819SStefano Zampini     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
731674ae819SStefano Zampini   }
732674ae819SStefano Zampini   total_counts += n_vertices;
733674ae819SStefano Zampini   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
734674ae819SStefano Zampini   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
735674ae819SStefano Zampini   total_counts = 0;
736674ae819SStefano Zampini   max_size_of_constraint = 0;
737674ae819SStefano Zampini   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
738674ae819SStefano Zampini     if (i<n_ISForEdges) {
739674ae819SStefano Zampini       used_IS = &ISForEdges[i];
740674ae819SStefano Zampini     } else {
741674ae819SStefano Zampini       used_IS = &ISForFaces[i-n_ISForEdges];
742674ae819SStefano Zampini     }
743674ae819SStefano Zampini     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
744674ae819SStefano Zampini     total_counts += j;
745674ae819SStefano Zampini     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
746674ae819SStefano Zampini   }
747674ae819SStefano Zampini   total_counts *= (nnsp_addone+nnsp_size);
748674ae819SStefano Zampini   total_counts += n_vertices;
749674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
750674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
751674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
752674ae819SStefano Zampini   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
753674ae819SStefano Zampini   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
754674ae819SStefano Zampini   for (i=0;i<pcis->n;i++) {
755674ae819SStefano Zampini     local_to_B[i]=-1;
756674ae819SStefano Zampini   }
757674ae819SStefano Zampini   for (i=0;i<pcis->n_B;i++) {
758674ae819SStefano Zampini     local_to_B[is_indices[i]]=i;
759674ae819SStefano Zampini   }
760674ae819SStefano Zampini   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
761674ae819SStefano Zampini 
762674ae819SStefano Zampini   /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */
763674ae819SStefano Zampini   rwork = 0;
764674ae819SStefano Zampini   work = 0;
765674ae819SStefano Zampini   singular_vals = 0;
766674ae819SStefano Zampini   temp_basis = 0;
767674ae819SStefano Zampini   correlation_mat = 0;
768674ae819SStefano Zampini   if (!pcbddc->use_nnsp_true) {
769674ae819SStefano Zampini     PetscScalar temp_work;
770674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
771674ae819SStefano Zampini     /* POD */
772674ae819SStefano Zampini     PetscInt max_n;
773674ae819SStefano Zampini     max_n = nnsp_addone+nnsp_size;
774674ae819SStefano Zampini     /* using some techniques borrowed from Proper Orthogonal Decomposition */
775674ae819SStefano Zampini     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
776674ae819SStefano Zampini     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr);
777674ae819SStefano Zampini     ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
778674ae819SStefano Zampini     ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
779674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX)
780674ae819SStefano Zampini     ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
781674ae819SStefano Zampini #endif
782674ae819SStefano Zampini     ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr);
783674ae819SStefano Zampini     ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr);
784674ae819SStefano Zampini     /* now we evaluate the optimal workspace using query with lwork=-1 */
785674ae819SStefano Zampini     ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr);
786674ae819SStefano Zampini     lwork=-1;
787674ae819SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
788674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
789674ae819SStefano Zampini     abs_tol=1.e-8;
790674ae819SStefano Zampini /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */
791f7c40c41SStefano Zampini     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));
792674ae819SStefano Zampini #else
793674ae819SStefano Zampini /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */
794674ae819SStefano Zampini /*  LAPACK call is missing here! TODO */
795674ae819SStefano Zampini     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
796674ae819SStefano Zampini #endif
797674ae819SStefano Zampini     if ( lierr ) {
798674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr);
799674ae819SStefano Zampini     }
800674ae819SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
801674ae819SStefano Zampini #else /* on missing GESVD */
802674ae819SStefano Zampini     /* SVD */
803674ae819SStefano Zampini     PetscInt max_n,min_n;
804674ae819SStefano Zampini     max_n = max_size_of_constraint;
805674ae819SStefano Zampini     min_n = nnsp_addone+nnsp_size;
806674ae819SStefano Zampini     if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) {
807674ae819SStefano Zampini       min_n = max_size_of_constraint;
808674ae819SStefano Zampini       max_n = nnsp_addone+nnsp_size;
809674ae819SStefano Zampini     }
810674ae819SStefano Zampini     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
811674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX)
812674ae819SStefano Zampini     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
813674ae819SStefano Zampini #endif
814674ae819SStefano Zampini     /* now we evaluate the optimal workspace using query with lwork=-1 */
815674ae819SStefano Zampini     lwork=-1;
816674ae819SStefano Zampini     ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr);
817674ae819SStefano Zampini     ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr);
818674ae819SStefano Zampini     dummy_int = Bs;
819674ae819SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
820674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
821f7c40c41SStefano Zampini     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));
822674ae819SStefano Zampini #else
823f7c40c41SStefano Zampini     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));
824674ae819SStefano Zampini #endif
825674ae819SStefano Zampini     if ( lierr ) {
826674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr);
827674ae819SStefano Zampini     }
828674ae819SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
829674ae819SStefano Zampini #endif
830674ae819SStefano Zampini     /* Allocate optimal workspace */
831674ae819SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
832674ae819SStefano Zampini     total_counts = (PetscInt)lwork;
833674ae819SStefano Zampini     ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr);
834674ae819SStefano Zampini   }
835674ae819SStefano Zampini   /* get local part of global near null space vectors */
836674ae819SStefano Zampini   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
837674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
838674ae819SStefano Zampini     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
839674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
840674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
841674ae819SStefano Zampini   }
842674ae819SStefano Zampini   /* Now we can loop on constraining sets */
843674ae819SStefano Zampini   total_counts = 0;
844674ae819SStefano Zampini   temp_indices[0] = 0;
845674ae819SStefano Zampini   /* vertices */
846674ae819SStefano Zampini   if (ISForVertices) {
847674ae819SStefano Zampini     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
848674ae819SStefano Zampini     if (nnsp_has_cnst) { /* consider all vertices */
849674ae819SStefano Zampini       for (i=0;i<n_vertices;i++) {
850674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
851674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
852674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
853674ae819SStefano Zampini         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
854674ae819SStefano Zampini         change_basis[total_counts]=PETSC_FALSE;
855674ae819SStefano Zampini         total_counts++;
856674ae819SStefano Zampini       }
857674ae819SStefano Zampini     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
858674ae819SStefano Zampini       for (i=0;i<n_vertices;i++) {
859674ae819SStefano Zampini         used_vertex=PETSC_FALSE;
860674ae819SStefano Zampini         k=0;
861674ae819SStefano Zampini         while (!used_vertex && k<nnsp_size) {
862674ae819SStefano Zampini           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
863674ae819SStefano Zampini           if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) {
864674ae819SStefano Zampini             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
865674ae819SStefano Zampini             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
866674ae819SStefano Zampini             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
867674ae819SStefano Zampini             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
868674ae819SStefano Zampini             change_basis[total_counts]=PETSC_FALSE;
869674ae819SStefano Zampini             total_counts++;
870674ae819SStefano Zampini             used_vertex=PETSC_TRUE;
871674ae819SStefano Zampini           }
872674ae819SStefano Zampini           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
873674ae819SStefano Zampini           k++;
874674ae819SStefano Zampini         }
875674ae819SStefano Zampini       }
876674ae819SStefano Zampini     }
877674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
878674ae819SStefano Zampini     n_vertices = total_counts;
879674ae819SStefano Zampini   }
880674ae819SStefano Zampini   /* edges and faces */
881674ae819SStefano Zampini   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
882674ae819SStefano Zampini     if (i<n_ISForEdges) {
883674ae819SStefano Zampini       used_IS = &ISForEdges[i];
88451b0f6cfSStefano Zampini       boolforchange = pcbddc->use_change_of_basis;
885674ae819SStefano Zampini     } else {
886674ae819SStefano Zampini       used_IS = &ISForFaces[i-n_ISForEdges];
88751b0f6cfSStefano Zampini       boolforchange = pcbddc->use_change_on_faces;
888674ae819SStefano Zampini     }
889674ae819SStefano Zampini     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
890674ae819SStefano Zampini     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
891674ae819SStefano Zampini     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
892674ae819SStefano Zampini     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
89351b0f6cfSStefano Zampini     /* HACK: change of basis should not performed on local periodic nodes */
89451b0f6cfSStefano Zampini     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) {
89551b0f6cfSStefano Zampini       boolforchange = PETSC_FALSE;
89651b0f6cfSStefano Zampini     }
897674ae819SStefano Zampini     if (nnsp_has_cnst) {
8985b08dc53SStefano Zampini       PetscScalar quad_value;
899674ae819SStefano Zampini       temp_constraints++;
900674ae819SStefano Zampini       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
901674ae819SStefano Zampini       for (j=0;j<size_of_constraint;j++) {
902674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
903674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
904674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
905674ae819SStefano Zampini       }
906674ae819SStefano Zampini       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
90751b0f6cfSStefano Zampini       change_basis[total_counts]=boolforchange;
908674ae819SStefano Zampini       total_counts++;
909674ae819SStefano Zampini     }
910674ae819SStefano Zampini     for (k=0;k<nnsp_size;k++) {
911674ae819SStefano Zampini       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
912674ae819SStefano Zampini       for (j=0;j<size_of_constraint;j++) {
913674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
914674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
915674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]];
916674ae819SStefano Zampini       }
917674ae819SStefano Zampini       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
9185b08dc53SStefano Zampini       real_value = 1.0;
919674ae819SStefano Zampini       if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */
920674ae819SStefano Zampini         ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
9215b08dc53SStefano Zampini         PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone));
922674ae819SStefano Zampini       }
9235b08dc53SStefano Zampini       if (real_value > 0.0) { /* keep indices and values */
924674ae819SStefano Zampini         temp_constraints++;
925674ae819SStefano Zampini         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
92651b0f6cfSStefano Zampini         change_basis[total_counts]=boolforchange;
927674ae819SStefano Zampini         total_counts++;
928674ae819SStefano Zampini       }
929674ae819SStefano Zampini     }
930674ae819SStefano Zampini     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
931674ae819SStefano Zampini     /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */
932674ae819SStefano Zampini     if (!use_nnsp_true) {
933674ae819SStefano Zampini       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
934674ae819SStefano Zampini       ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr);
935674ae819SStefano Zampini 
936674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
937674ae819SStefano Zampini       ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr);
938674ae819SStefano Zampini       /* Store upper triangular part of correlation matrix */
939674ae819SStefano Zampini       for (j=0;j<temp_constraints;j++) {
940674ae819SStefano Zampini         for (k=0;k<j+1;k++) {
9415b08dc53SStefano Zampini           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));
942674ae819SStefano Zampini 
943674ae819SStefano Zampini         }
944674ae819SStefano Zampini       }
945674ae819SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
946674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
947674ae819SStefano Zampini /*      LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */
948f7c40c41SStefano Zampini       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));
949674ae819SStefano Zampini #else
950674ae819SStefano Zampini /*  LAPACK call is missing here! TODO */
951674ae819SStefano Zampini       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
952674ae819SStefano Zampini #endif
953674ae819SStefano Zampini       if (lierr) {
954674ae819SStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr);
955674ae819SStefano Zampini       }
956674ae819SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
957674ae819SStefano Zampini       /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */
958674ae819SStefano Zampini       j=0;
959674ae819SStefano Zampini       while (j < Bt && singular_vals[j] < tol) j++;
960674ae819SStefano Zampini       total_counts=total_counts-j;
961674ae819SStefano Zampini       if (j<temp_constraints) {
962674ae819SStefano Zampini         for (k=j;k<Bt;k++) {
963674ae819SStefano Zampini           singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
964674ae819SStefano Zampini         }
965674ae819SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
966f7c40c41SStefano Zampini         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));
967674ae819SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
968674ae819SStefano Zampini         /* copy POD basis into used quadrature memory */
969674ae819SStefano Zampini         for (k=0;k<Bt-j;k++) {
970674ae819SStefano Zampini           for (ii=0;ii<size_of_constraint;ii++) {
971674ae819SStefano Zampini             temp_quadrature_constraint[temp_indices[temp_start_ptr+k]+ii]=singular_vals[Bt-1-k]*temp_basis[(Bt-1-k)*size_of_constraint+ii];
972674ae819SStefano Zampini           }
973674ae819SStefano Zampini         }
974674ae819SStefano Zampini       }
975674ae819SStefano Zampini 
976674ae819SStefano Zampini #else  /* on missing GESVD */
977674ae819SStefano Zampini       PetscInt min_n = temp_constraints;
978674ae819SStefano Zampini       if (min_n > size_of_constraint) min_n = size_of_constraint;
979674ae819SStefano Zampini       dummy_int = Bs;
980674ae819SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
981674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
982f7c40c41SStefano Zampini       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));
983674ae819SStefano Zampini #else
984f7c40c41SStefano Zampini       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));
985674ae819SStefano Zampini #endif
986674ae819SStefano Zampini       if (lierr) {
987674ae819SStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr);
988674ae819SStefano Zampini       }
989674ae819SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
990674ae819SStefano Zampini       /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */
991674ae819SStefano Zampini       j = 0;
992674ae819SStefano Zampini       while (j < min_n && singular_vals[min_n-j-1] < tol) j++;
993674ae819SStefano Zampini       total_counts = total_counts-(PetscInt)Bt+(min_n-j);
994674ae819SStefano Zampini #endif
995674ae819SStefano Zampini     }
996674ae819SStefano Zampini   }
997674ae819SStefano Zampini   /* free index sets of faces, edges and vertices */
998674ae819SStefano Zampini   for (i=0;i<n_ISForFaces;i++) {
999674ae819SStefano Zampini     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
1000674ae819SStefano Zampini   }
1001674ae819SStefano Zampini   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
1002674ae819SStefano Zampini   for (i=0;i<n_ISForEdges;i++) {
1003674ae819SStefano Zampini     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
1004674ae819SStefano Zampini   }
1005674ae819SStefano Zampini   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
1006674ae819SStefano Zampini   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
1007674ae819SStefano Zampini 
1008674ae819SStefano Zampini   /* set quantities in pcbddc data structure */
1009674ae819SStefano Zampini   /* n_vertices defines the number of point primal dofs */
1010674ae819SStefano Zampini   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
1011674ae819SStefano Zampini   local_primal_size = total_counts;
1012674ae819SStefano Zampini   pcbddc->n_vertices = n_vertices;
1013674ae819SStefano Zampini   pcbddc->n_constraints = total_counts-n_vertices;
1014674ae819SStefano Zampini   pcbddc->local_primal_size = local_primal_size;
1015674ae819SStefano Zampini 
1016674ae819SStefano Zampini   /* Create constraint matrix */
1017674ae819SStefano Zampini   /* The constraint matrix is used to compute the l2g map of primal dofs */
1018674ae819SStefano Zampini   /* so we need to set it up properly either with or without change of basis */
1019674ae819SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
1020674ae819SStefano Zampini   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
1021674ae819SStefano Zampini   ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr);
1022674ae819SStefano Zampini   /* compute a local numbering of constraints : vertices first then constraints */
1023674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
1024674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
1025674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
1026674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr);
1027674ae819SStefano Zampini   total_counts=0;
1028674ae819SStefano Zampini   /* find vertices: subdomain corners plus dofs with basis changed */
1029674ae819SStefano Zampini   for (i=0;i<local_primal_size;i++) {
1030674ae819SStefano Zampini     size_of_constraint=temp_indices[i+1]-temp_indices[i];
1031674ae819SStefano Zampini     if (change_basis[i] || size_of_constraint == 1) {
1032674ae819SStefano Zampini       k=0;
1033674ae819SStefano Zampini       while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) {
1034674ae819SStefano Zampini         k=k+1;
1035674ae819SStefano Zampini       }
1036674ae819SStefano Zampini       j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1];
1037674ae819SStefano Zampini       array_vector[j] = 1.0;
1038674ae819SStefano Zampini       aux_primal_numbering[total_counts]=j;
1039674ae819SStefano Zampini       aux_primal_permutation[total_counts]=total_counts;
1040674ae819SStefano Zampini       total_counts++;
1041674ae819SStefano Zampini     }
1042674ae819SStefano Zampini   }
1043674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
1044674ae819SStefano Zampini   /* permute indices in order to have a sorted set of vertices */
1045674ae819SStefano Zampini   ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation);
1046674ae819SStefano Zampini   /* nonzero structure */
1047674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
1048674ae819SStefano Zampini   for (i=0;i<total_counts;i++) {
1049674ae819SStefano Zampini     nnz[i]=1;
1050674ae819SStefano Zampini   }
1051674ae819SStefano Zampini   j=total_counts;
1052674ae819SStefano Zampini   for (i=n_vertices;i<local_primal_size;i++) {
1053674ae819SStefano Zampini     if (!change_basis[i]) {
1054674ae819SStefano Zampini       nnz[j]=temp_indices[i+1]-temp_indices[i];
1055674ae819SStefano Zampini       j++;
1056674ae819SStefano Zampini     }
1057674ae819SStefano Zampini   }
1058674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
1059674ae819SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1060674ae819SStefano Zampini   /* set values in constraint matrix */
1061674ae819SStefano Zampini   for (i=0;i<total_counts;i++) {
1062674ae819SStefano Zampini     j = aux_primal_permutation[i];
1063674ae819SStefano Zampini     k = aux_primal_numbering[j];
1064674ae819SStefano Zampini     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr);
1065674ae819SStefano Zampini   }
1066674ae819SStefano Zampini   for (i=n_vertices;i<local_primal_size;i++) {
1067674ae819SStefano Zampini     if (!change_basis[i]) {
1068674ae819SStefano Zampini       size_of_constraint=temp_indices[i+1]-temp_indices[i];
1069674ae819SStefano Zampini       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);
1070674ae819SStefano Zampini       total_counts++;
1071674ae819SStefano Zampini     }
1072674ae819SStefano Zampini   }
1073674ae819SStefano Zampini   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
1074674ae819SStefano Zampini   ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr);
1075674ae819SStefano Zampini   /* assembling */
1076674ae819SStefano Zampini   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1077674ae819SStefano Zampini   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1078674ae819SStefano Zampini 
1079674ae819SStefano Zampini   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
1080674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
1081674ae819SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
1082674ae819SStefano Zampini     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
1083674ae819SStefano Zampini     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
1084674ae819SStefano Zampini     /* work arrays */
1085674ae819SStefano Zampini     /* we need to reuse these arrays, so we free them */
1086674ae819SStefano Zampini     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
1087674ae819SStefano Zampini     ierr = PetscFree(work);CHKERRQ(ierr);
1088674ae819SStefano Zampini     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
1089674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
1090674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr);
1091674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr);
1092674ae819SStefano Zampini     for (i=0;i<pcis->n_B;i++) {
1093674ae819SStefano Zampini       nnz[i]=1;
1094674ae819SStefano Zampini     }
1095674ae819SStefano Zampini     /* Overestimated nonzeros per row */
1096674ae819SStefano Zampini     k=1;
1097674ae819SStefano Zampini     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
1098674ae819SStefano Zampini       if (change_basis[i]) {
1099674ae819SStefano Zampini         size_of_constraint = temp_indices[i+1]-temp_indices[i];
1100674ae819SStefano Zampini         if (k < size_of_constraint) {
1101674ae819SStefano Zampini           k = size_of_constraint;
1102674ae819SStefano Zampini         }
1103674ae819SStefano Zampini         for (j=0;j<size_of_constraint;j++) {
1104674ae819SStefano Zampini           nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
1105674ae819SStefano Zampini         }
1106674ae819SStefano Zampini       }
1107674ae819SStefano Zampini     }
1108674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
1109674ae819SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
1110674ae819SStefano Zampini     /* Temporary array to store indices */
1111674ae819SStefano Zampini     ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr);
1112674ae819SStefano Zampini     /* Set initial identity in the matrix */
1113674ae819SStefano Zampini     for (i=0;i<pcis->n_B;i++) {
1114674ae819SStefano Zampini       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
1115674ae819SStefano Zampini     }
1116674ae819SStefano Zampini     /* Now we loop on the constraints which need a change of basis */
1117674ae819SStefano Zampini     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
1118674ae819SStefano Zampini     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */
1119674ae819SStefano Zampini     temp_constraints = 0;
1120674ae819SStefano Zampini     if (pcbddc->n_vertices < local_primal_size) {
1121674ae819SStefano Zampini       temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]];
1122674ae819SStefano Zampini     }
1123674ae819SStefano Zampini     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
1124674ae819SStefano Zampini       if (change_basis[i]) {
1125674ae819SStefano Zampini         compute_submatrix = PETSC_FALSE;
1126674ae819SStefano Zampini         useksp = PETSC_FALSE;
1127674ae819SStefano Zampini         if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) {
1128674ae819SStefano Zampini           temp_constraints++;
1129674ae819SStefano Zampini           if (i == local_primal_size -1 ||  temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) {
1130674ae819SStefano Zampini             compute_submatrix = PETSC_TRUE;
1131674ae819SStefano Zampini           }
1132674ae819SStefano Zampini         }
1133674ae819SStefano Zampini         if (compute_submatrix) {
1134674ae819SStefano Zampini           if (temp_constraints > 1 || pcbddc->use_nnsp_true) {
1135674ae819SStefano Zampini             useksp = PETSC_TRUE;
1136674ae819SStefano Zampini           }
1137674ae819SStefano Zampini           size_of_constraint = temp_indices[i+1]-temp_indices[i];
1138674ae819SStefano Zampini           if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */
1139674ae819SStefano Zampini             ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr);
1140674ae819SStefano Zampini             ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr);
1141674ae819SStefano Zampini             ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr);
1142674ae819SStefano Zampini             ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr);
1143674ae819SStefano Zampini           }
1144674ae819SStefano Zampini           /* First _size_of_constraint-temp_constraints_ columns */
1145674ae819SStefano Zampini           dual_dofs = size_of_constraint-temp_constraints;
1146674ae819SStefano Zampini           start_constraint = i+1-temp_constraints;
1147674ae819SStefano Zampini           for (s=0;s<dual_dofs;s++) {
1148674ae819SStefano Zampini             is_indices[0] = s;
1149674ae819SStefano Zampini             for (j=0;j<temp_constraints;j++) {
1150674ae819SStefano Zampini               for (k=0;k<temp_constraints;k++) {
1151674ae819SStefano Zampini                 temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1];
1152674ae819SStefano Zampini               }
1153674ae819SStefano Zampini               work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s];
1154674ae819SStefano Zampini               is_indices[j+1]=s+j+1;
1155674ae819SStefano Zampini             }
1156674ae819SStefano Zampini             Bt = temp_constraints;
1157674ae819SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1158f7c40c41SStefano Zampini             PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr));
1159674ae819SStefano Zampini             if ( lierr ) {
1160674ae819SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr);
1161674ae819SStefano Zampini             }
1162674ae819SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1163674ae819SStefano Zampini             j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s];
1164674ae819SStefano Zampini             ierr = MatSetValues(pcbddc->ChangeOfBasisMatrix,temp_constraints,&temp_indices_to_constraint_B[temp_indices[start_constraint]+s+1],1,&j,work,INSERT_VALUES);CHKERRQ(ierr);
1165674ae819SStefano Zampini             if (useksp) {
1166674ae819SStefano Zampini               /* temp mat with transposed rows and columns */
1167674ae819SStefano Zampini               ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr);
1168674ae819SStefano Zampini               ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr);
1169674ae819SStefano Zampini             }
1170674ae819SStefano Zampini           }
1171674ae819SStefano Zampini           if (useksp) {
1172674ae819SStefano Zampini             /* last rows of temp_mat */
1173674ae819SStefano Zampini             for (j=0;j<size_of_constraint;j++) {
1174674ae819SStefano Zampini               is_indices[j] = j;
1175674ae819SStefano Zampini             }
1176674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
1177674ae819SStefano Zampini               k = s + dual_dofs;
1178674ae819SStefano Zampini               ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr);
1179674ae819SStefano Zampini             }
1180674ae819SStefano Zampini             ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1181674ae819SStefano Zampini             ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1182674ae819SStefano Zampini             ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr);
1183674ae819SStefano Zampini             ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr);
1184674ae819SStefano Zampini             ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr);
1185674ae819SStefano Zampini             ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr);
1186674ae819SStefano Zampini             ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr);
1187674ae819SStefano Zampini             ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr);
1188674ae819SStefano Zampini             ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr);
1189674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
1190674ae819SStefano Zampini               ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr);
1191674ae819SStefano Zampini               ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr);
1192674ae819SStefano Zampini               ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr);
1193674ae819SStefano Zampini               ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr);
1194674ae819SStefano Zampini               ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr);
1195674ae819SStefano Zampini               ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr);
1196674ae819SStefano Zampini               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
1197674ae819SStefano Zampini               /* last columns of change of basis matrix associated to new primal dofs */
1198674ae819SStefano Zampini               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);
1199674ae819SStefano Zampini               ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr);
1200674ae819SStefano Zampini             }
1201674ae819SStefano Zampini             ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
1202674ae819SStefano Zampini             ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr);
1203674ae819SStefano Zampini             ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1204674ae819SStefano Zampini           } else {
1205674ae819SStefano Zampini             /* last columns of change of basis matrix associated to new primal dofs */
1206674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
1207674ae819SStefano Zampini               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
1208674ae819SStefano Zampini               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);
1209674ae819SStefano Zampini             }
1210674ae819SStefano Zampini           }
1211674ae819SStefano Zampini           /* prepare for the next cycle */
1212674ae819SStefano Zampini           temp_constraints = 0;
1213674ae819SStefano Zampini           if (i != local_primal_size -1 ) {
1214674ae819SStefano Zampini             temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]];
1215674ae819SStefano Zampini           }
1216674ae819SStefano Zampini         }
1217674ae819SStefano Zampini       }
1218674ae819SStefano Zampini     }
1219674ae819SStefano Zampini     /* assembling */
1220674ae819SStefano Zampini     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1221674ae819SStefano Zampini     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1222674ae819SStefano Zampini     ierr = PetscFree(ipiv);CHKERRQ(ierr);
1223674ae819SStefano Zampini     ierr = PetscFree(is_indices);CHKERRQ(ierr);
1224674ae819SStefano Zampini   }
1225674ae819SStefano Zampini   /* free workspace no longer needed */
1226674ae819SStefano Zampini   ierr = PetscFree(rwork);CHKERRQ(ierr);
1227674ae819SStefano Zampini   ierr = PetscFree(work);CHKERRQ(ierr);
1228674ae819SStefano Zampini   ierr = PetscFree(temp_basis);CHKERRQ(ierr);
1229674ae819SStefano Zampini   ierr = PetscFree(singular_vals);CHKERRQ(ierr);
1230674ae819SStefano Zampini   ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
1231674ae819SStefano Zampini   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
1232674ae819SStefano Zampini   ierr = PetscFree(change_basis);CHKERRQ(ierr);
1233674ae819SStefano Zampini   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
1234674ae819SStefano Zampini   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
1235674ae819SStefano Zampini   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
1236674ae819SStefano Zampini   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
1237674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
1238674ae819SStefano Zampini   ierr = PetscFree(iwork);CHKERRQ(ierr);
1239674ae819SStefano Zampini   ierr = PetscFree(ifail);CHKERRQ(ierr);
1240674ae819SStefano Zampini   ierr = PetscFree(singular_vectors);CHKERRQ(ierr);
1241674ae819SStefano Zampini #endif
1242674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
1243674ae819SStefano Zampini     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
1244674ae819SStefano Zampini   }
1245674ae819SStefano Zampini   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
1246674ae819SStefano Zampini   PetscFunctionReturn(0);
1247674ae819SStefano Zampini }
1248674ae819SStefano Zampini 
1249674ae819SStefano Zampini #undef __FUNCT__
1250674ae819SStefano Zampini #define __FUNCT__ "PCBDDCAnalyzeInterface"
1251674ae819SStefano Zampini PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
1252674ae819SStefano Zampini {
1253674ae819SStefano Zampini   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
1254674ae819SStefano Zampini   PC_IS       *pcis = (PC_IS*)pc->data;
1255674ae819SStefano Zampini   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
1256674ae819SStefano Zampini   PetscInt    bs,ierr,i,vertex_size;
1257674ae819SStefano Zampini   PetscViewer viewer=pcbddc->dbg_viewer;
1258674ae819SStefano Zampini 
1259674ae819SStefano Zampini   PetscFunctionBegin;
1260674ae819SStefano Zampini   /* Init local Graph struct */
1261674ae819SStefano Zampini   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
1262674ae819SStefano Zampini 
1263575ad6abSStefano Zampini   /* Check validity of the csr graph passed in by the user */
1264575ad6abSStefano Zampini   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
1265575ad6abSStefano Zampini     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
1266575ad6abSStefano Zampini   }
1267674ae819SStefano Zampini   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
1268674ae819SStefano Zampini   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
1269674ae819SStefano Zampini     Mat mat_adj;
1270674ae819SStefano Zampini     const PetscInt *xadj,*adjncy;
1271674ae819SStefano Zampini     PetscBool flg_row=PETSC_TRUE;
1272674ae819SStefano Zampini 
1273674ae819SStefano Zampini     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
1274674ae819SStefano Zampini     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1275674ae819SStefano Zampini     if (!flg_row) {
1276674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
1277674ae819SStefano Zampini     }
1278674ae819SStefano Zampini     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
1279674ae819SStefano Zampini     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1280674ae819SStefano Zampini     if (!flg_row) {
1281674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
1282674ae819SStefano Zampini     }
1283674ae819SStefano Zampini     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
1284674ae819SStefano Zampini   }
1285674ae819SStefano Zampini 
1286674ae819SStefano Zampini   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
1287674ae819SStefano Zampini   vertex_size = 1;
1288674ae819SStefano Zampini   if (!pcbddc->n_ISForDofs) {
1289674ae819SStefano Zampini     IS *custom_ISForDofs;
1290674ae819SStefano Zampini 
1291674ae819SStefano Zampini     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
1292674ae819SStefano Zampini     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
1293674ae819SStefano Zampini     for (i=0;i<bs;i++) {
1294674ae819SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
1295674ae819SStefano Zampini     }
1296674ae819SStefano Zampini     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
1297674ae819SStefano Zampini     /* remove my references to IS objects */
1298674ae819SStefano Zampini     for (i=0;i<bs;i++) {
1299674ae819SStefano Zampini       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
1300674ae819SStefano Zampini     }
1301674ae819SStefano Zampini     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
1302674ae819SStefano Zampini   } else { /* mat block size as vertex size (used for elasticity) */
1303674ae819SStefano Zampini     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
1304674ae819SStefano Zampini   }
1305674ae819SStefano Zampini 
1306674ae819SStefano Zampini   /* Setup of Graph */
1307674ae819SStefano Zampini   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
1308674ae819SStefano Zampini 
1309674ae819SStefano Zampini   /* Graph's connected components analysis */
1310674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
1311674ae819SStefano Zampini 
1312674ae819SStefano Zampini   /* print some info to stdout */
1313674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
1314e49050b4SStefano Zampini     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
1315674ae819SStefano Zampini   }
1316674ae819SStefano Zampini   PetscFunctionReturn(0);
1317674ae819SStefano Zampini }
1318674ae819SStefano Zampini 
1319674ae819SStefano Zampini #undef __FUNCT__
1320674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
1321674ae819SStefano Zampini PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
1322674ae819SStefano Zampini {
1323674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1324674ae819SStefano Zampini   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
1325674ae819SStefano Zampini   PetscErrorCode ierr;
1326674ae819SStefano Zampini 
1327674ae819SStefano Zampini   PetscFunctionBegin;
1328674ae819SStefano Zampini   n = 0;
1329674ae819SStefano Zampini   vertices = 0;
1330674ae819SStefano Zampini   if (pcbddc->ConstraintMatrix) {
1331674ae819SStefano Zampini     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
1332b120a5c6SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1333b120a5c6SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1334b120a5c6SStefano Zampini       if (size_of_constraint == 1) n++;
1335b120a5c6SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1336b120a5c6SStefano Zampini     }
1337811e8ca2SStefano Zampini     if (vertices_idx) {
1338b120a5c6SStefano Zampini       ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
1339b120a5c6SStefano Zampini       n = 0;
1340674ae819SStefano Zampini       for (i=0;i<local_primal_size;i++) {
1341674ae819SStefano Zampini         ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1342674ae819SStefano Zampini         if (size_of_constraint == 1) {
1343674ae819SStefano Zampini           vertices[n++]=row_cmat_indices[0];
1344674ae819SStefano Zampini         }
1345674ae819SStefano Zampini         ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1346674ae819SStefano Zampini       }
1347674ae819SStefano Zampini     }
1348811e8ca2SStefano Zampini   }
1349674ae819SStefano Zampini   *n_vertices = n;
1350811e8ca2SStefano Zampini   if (vertices_idx) *vertices_idx = vertices;
1351674ae819SStefano Zampini   PetscFunctionReturn(0);
1352674ae819SStefano Zampini }
1353674ae819SStefano Zampini 
1354674ae819SStefano Zampini #undef __FUNCT__
1355674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
1356674ae819SStefano Zampini PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
1357674ae819SStefano Zampini {
1358674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1359674ae819SStefano Zampini   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
1360674ae819SStefano Zampini   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
1361674ae819SStefano Zampini   PetscBool      *touched;
1362674ae819SStefano Zampini   PetscErrorCode ierr;
1363674ae819SStefano Zampini 
1364674ae819SStefano Zampini   PetscFunctionBegin;
1365674ae819SStefano Zampini   n = 0;
1366674ae819SStefano Zampini   constraints_index = 0;
1367674ae819SStefano Zampini   if (pcbddc->ConstraintMatrix) {
1368674ae819SStefano Zampini     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
1369674ae819SStefano Zampini     max_size_of_constraint = 0;
1370674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1371674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1372674ae819SStefano Zampini       if (size_of_constraint > 1) {
1373674ae819SStefano Zampini         n++;
1374674ae819SStefano Zampini       }
1375674ae819SStefano Zampini       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1376674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1377674ae819SStefano Zampini     }
1378811e8ca2SStefano Zampini     if (constraints_idx) {
1379674ae819SStefano Zampini       ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1380674ae819SStefano Zampini       ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1381674ae819SStefano Zampini       ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1382674ae819SStefano Zampini       ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1383674ae819SStefano Zampini       n = 0;
1384674ae819SStefano Zampini       for (i=0;i<local_primal_size;i++) {
1385674ae819SStefano Zampini         ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1386674ae819SStefano Zampini         if (size_of_constraint > 1) {
1387674ae819SStefano Zampini           ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
138882d3d8afSStefano Zampini           /* find first untouched local node */
138982d3d8afSStefano Zampini           j = 0;
139082d3d8afSStefano Zampini           while(touched[row_cmat_indices[j]]) j++;
139182d3d8afSStefano Zampini           min_index = row_cmat_global_indices[j];
139282d3d8afSStefano Zampini           min_loc = j;
139382d3d8afSStefano Zampini           /* search the minimum among nodes not yet touched on the connected component
139482d3d8afSStefano Zampini              since there can be more than one constraint on a single cc */
1395674ae819SStefano Zampini           for (j=1;j<size_of_constraint;j++) {
1396674ae819SStefano Zampini             if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1397674ae819SStefano Zampini               min_index = row_cmat_global_indices[j];
1398674ae819SStefano Zampini               min_loc = j;
1399674ae819SStefano Zampini             }
1400674ae819SStefano Zampini           }
1401674ae819SStefano Zampini           touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1402674ae819SStefano Zampini           constraints_index[n++] = row_cmat_indices[min_loc];
1403674ae819SStefano Zampini         }
1404674ae819SStefano Zampini         ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1405674ae819SStefano Zampini       }
1406674ae819SStefano Zampini       ierr = PetscFree(touched);CHKERRQ(ierr);
1407674ae819SStefano Zampini       ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1408811e8ca2SStefano Zampini     }
1409811e8ca2SStefano Zampini   }
1410674ae819SStefano Zampini   *n_constraints = n;
1411811e8ca2SStefano Zampini   if (constraints_idx) *constraints_idx = constraints_index;
1412674ae819SStefano Zampini   PetscFunctionReturn(0);
1413674ae819SStefano Zampini }
1414674ae819SStefano Zampini 
1415674ae819SStefano Zampini /* the next two functions has been adapted from pcis.c */
1416674ae819SStefano Zampini #undef __FUNCT__
1417674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchur"
1418674ae819SStefano Zampini PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1419674ae819SStefano Zampini {
1420674ae819SStefano Zampini   PetscErrorCode ierr;
1421674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
1422674ae819SStefano Zampini 
1423674ae819SStefano Zampini   PetscFunctionBegin;
1424674ae819SStefano Zampini   if (!vec2_B) { vec2_B = v; }
1425674ae819SStefano Zampini   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1426674ae819SStefano Zampini   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1427674ae819SStefano Zampini   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1428674ae819SStefano Zampini   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1429674ae819SStefano Zampini   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1430674ae819SStefano Zampini   PetscFunctionReturn(0);
1431674ae819SStefano Zampini }
1432674ae819SStefano Zampini 
1433674ae819SStefano Zampini #undef __FUNCT__
1434674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchurTranspose"
1435674ae819SStefano Zampini PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1436674ae819SStefano Zampini {
1437674ae819SStefano Zampini   PetscErrorCode ierr;
1438674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
1439674ae819SStefano Zampini 
1440674ae819SStefano Zampini   PetscFunctionBegin;
1441674ae819SStefano Zampini   if (!vec2_B) { vec2_B = v; }
1442674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1443674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1444674ae819SStefano Zampini   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1445674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1446674ae819SStefano Zampini   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1447674ae819SStefano Zampini   PetscFunctionReturn(0);
1448674ae819SStefano Zampini }
1449674ae819SStefano Zampini 
1450674ae819SStefano Zampini #undef __FUNCT__
1451674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSubsetNumbering"
1452674ae819SStefano Zampini 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[])
1453674ae819SStefano Zampini {
1454674ae819SStefano Zampini   Vec            local_vec,global_vec;
1455674ae819SStefano Zampini   IS             seqis,paris;
1456674ae819SStefano Zampini   VecScatter     scatter_ctx;
1457674ae819SStefano Zampini   PetscScalar    *array;
1458674ae819SStefano Zampini   PetscInt       *temp_global_dofs;
1459674ae819SStefano Zampini   PetscScalar    globalsum;
1460674ae819SStefano Zampini   PetscInt       i,j,s;
1461674ae819SStefano Zampini   PetscInt       nlocals,first_index,old_index,max_local;
1462674ae819SStefano Zampini   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1463674ae819SStefano Zampini   PetscMPIInt    *dof_sizes,*dof_displs;
1464674ae819SStefano Zampini   PetscBool      first_found;
1465674ae819SStefano Zampini   PetscErrorCode ierr;
1466674ae819SStefano Zampini 
1467674ae819SStefano Zampini   PetscFunctionBegin;
1468674ae819SStefano Zampini   /* mpi buffers */
1469674ae819SStefano Zampini   MPI_Comm_size(comm,&size_prec_comm);
1470674ae819SStefano Zampini   MPI_Comm_rank(comm,&rank_prec_comm);
1471674ae819SStefano Zampini   j = ( !rank_prec_comm ? size_prec_comm : 0);
1472674ae819SStefano Zampini   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1473674ae819SStefano Zampini   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1474674ae819SStefano Zampini   /* get maximum size of subset */
1475674ae819SStefano Zampini   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1476674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1477674ae819SStefano Zampini   max_local = 0;
1478674ae819SStefano Zampini   if (n_local_dofs) {
1479674ae819SStefano Zampini     max_local = temp_global_dofs[0];
1480674ae819SStefano Zampini     for (i=1;i<n_local_dofs;i++) {
1481674ae819SStefano Zampini       if (max_local < temp_global_dofs[i] ) {
1482674ae819SStefano Zampini         max_local = temp_global_dofs[i];
1483674ae819SStefano Zampini       }
1484674ae819SStefano Zampini     }
1485674ae819SStefano Zampini   }
1486674ae819SStefano Zampini   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1487674ae819SStefano Zampini   max_global++;
1488674ae819SStefano Zampini   max_local = 0;
1489674ae819SStefano Zampini   if (n_local_dofs) {
1490674ae819SStefano Zampini     max_local = local_dofs[0];
1491674ae819SStefano Zampini     for (i=1;i<n_local_dofs;i++) {
1492674ae819SStefano Zampini       if (max_local < local_dofs[i] ) {
1493674ae819SStefano Zampini         max_local = local_dofs[i];
1494674ae819SStefano Zampini       }
1495674ae819SStefano Zampini     }
1496674ae819SStefano Zampini   }
1497674ae819SStefano Zampini   max_local++;
1498674ae819SStefano Zampini   /* allocate workspace */
1499674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1500674ae819SStefano Zampini   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1501674ae819SStefano Zampini   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1502674ae819SStefano Zampini   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1503674ae819SStefano Zampini   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1504674ae819SStefano Zampini   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1505674ae819SStefano Zampini   /* create scatter */
1506674ae819SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1507674ae819SStefano Zampini   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1508674ae819SStefano Zampini   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1509674ae819SStefano Zampini   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1510674ae819SStefano Zampini   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1511674ae819SStefano Zampini   /* init array */
1512674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1513674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1514674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1515674ae819SStefano Zampini   if (local_dofs_mult) {
1516674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1517674ae819SStefano Zampini       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1518674ae819SStefano Zampini     }
1519674ae819SStefano Zampini   } else {
1520674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1521674ae819SStefano Zampini       array[local_dofs[i]]=1.0;
1522674ae819SStefano Zampini     }
1523674ae819SStefano Zampini   }
1524674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1525674ae819SStefano Zampini   /* scatter into global vec and get total number of global dofs */
1526674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1527674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1528674ae819SStefano Zampini   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
15295b08dc53SStefano Zampini   *n_global_subset = (PetscInt)PetscRealPart(globalsum);
1530674ae819SStefano Zampini   /* Fill global_vec with cumulative function for global numbering */
1531674ae819SStefano Zampini   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1532674ae819SStefano Zampini   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1533674ae819SStefano Zampini   nlocals = 0;
1534674ae819SStefano Zampini   first_index = -1;
1535674ae819SStefano Zampini   first_found = PETSC_FALSE;
1536674ae819SStefano Zampini   for (i=0;i<s;i++) {
15375b08dc53SStefano Zampini     if (!first_found && PetscRealPart(array[i]) > 0.0) {
1538674ae819SStefano Zampini       first_found = PETSC_TRUE;
1539674ae819SStefano Zampini       first_index = i;
1540674ae819SStefano Zampini     }
15415b08dc53SStefano Zampini     nlocals += (PetscInt)PetscRealPart(array[i]);
1542674ae819SStefano Zampini   }
1543674ae819SStefano Zampini   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1544674ae819SStefano Zampini   if (!rank_prec_comm) {
1545674ae819SStefano Zampini     dof_displs[0]=0;
1546674ae819SStefano Zampini     for (i=1;i<size_prec_comm;i++) {
1547674ae819SStefano Zampini       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1548674ae819SStefano Zampini     }
1549674ae819SStefano Zampini   }
1550674ae819SStefano Zampini   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1551674ae819SStefano Zampini   if (first_found) {
1552674ae819SStefano Zampini     array[first_index] += (PetscScalar)nlocals;
1553674ae819SStefano Zampini     old_index = first_index;
1554674ae819SStefano Zampini     for (i=first_index+1;i<s;i++) {
15555b08dc53SStefano Zampini       if (PetscRealPart(array[i]) > 0.0) {
1556674ae819SStefano Zampini         array[i] += array[old_index];
1557674ae819SStefano Zampini         old_index = i;
1558674ae819SStefano Zampini       }
1559674ae819SStefano Zampini     }
1560674ae819SStefano Zampini   }
1561674ae819SStefano Zampini   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1562674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1563674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1564674ae819SStefano Zampini   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1565674ae819SStefano Zampini   /* get global ordering of local dofs */
1566674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1567674ae819SStefano Zampini   if (local_dofs_mult) {
1568674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
15695b08dc53SStefano Zampini       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i];
1570674ae819SStefano Zampini     }
1571674ae819SStefano Zampini   } else {
1572674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
15735b08dc53SStefano Zampini       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1;
1574674ae819SStefano Zampini     }
1575674ae819SStefano Zampini   }
1576674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1577674ae819SStefano Zampini   /* free workspace */
1578674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1579674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1580674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1581674ae819SStefano Zampini   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1582674ae819SStefano Zampini   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1583674ae819SStefano Zampini   /* return pointer to global ordering of local dofs */
1584674ae819SStefano Zampini   *global_numbering_subset = temp_global_dofs;
1585674ae819SStefano Zampini   PetscFunctionReturn(0);
1586674ae819SStefano Zampini }
1587