xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision f7c40c414c665fccbfb9d94fa7258db2e834b86f)
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 = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr);
51674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->coarse_vec);CHKERRQ(ierr);
52674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->coarse_rhs);CHKERRQ(ierr);
53674ae819SStefano Zampini   ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr);
54674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_mat);CHKERRQ(ierr);
55674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_phi_B);CHKERRQ(ierr);
56674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_phi_D);CHKERRQ(ierr);
5715aaf578SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_psi_B);CHKERRQ(ierr);
5815aaf578SStefano Zampini   ierr = MatDestroy(&pcbddc->coarse_psi_D);CHKERRQ(ierr);
59674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec1_P);CHKERRQ(ierr);
60674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec1_C);CHKERRQ(ierr);
61674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->local_auxmat1);CHKERRQ(ierr);
62674ae819SStefano Zampini   ierr = MatDestroy(&pcbddc->local_auxmat2);CHKERRQ(ierr);
63674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec1_R);CHKERRQ(ierr);
64674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec2_R);CHKERRQ(ierr);
65674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->vec4_D);CHKERRQ(ierr);
66674ae819SStefano Zampini   ierr = VecScatterDestroy(&pcbddc->R_to_B);CHKERRQ(ierr);
67674ae819SStefano Zampini   ierr = VecScatterDestroy(&pcbddc->R_to_D);CHKERRQ(ierr);
68674ae819SStefano Zampini   ierr = VecScatterDestroy(&pcbddc->coarse_loc_to_glob);CHKERRQ(ierr);
69674ae819SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr);
70674ae819SStefano Zampini   ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr);
71674ae819SStefano Zampini   ierr = PetscFree(pcbddc->local_primal_indices);CHKERRQ(ierr);
72674ae819SStefano Zampini   ierr = PetscFree(pcbddc->replicated_local_primal_indices);CHKERRQ(ierr);
73674ae819SStefano Zampini   ierr = PetscFree(pcbddc->replicated_local_primal_values);CHKERRQ(ierr);
74674ae819SStefano Zampini   ierr = PetscFree(pcbddc->local_primal_displacements);CHKERRQ(ierr);
75674ae819SStefano Zampini   ierr = PetscFree(pcbddc->local_primal_sizes);CHKERRQ(ierr);
76674ae819SStefano Zampini   PetscFunctionReturn(0);
77674ae819SStefano Zampini }
78674ae819SStefano Zampini 
79674ae819SStefano Zampini #undef __FUNCT__
80674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSolveSaddlePoint"
81674ae819SStefano Zampini static PetscErrorCode  PCBDDCSolveSaddlePoint(PC pc)
82674ae819SStefano Zampini {
83674ae819SStefano Zampini   PetscErrorCode ierr;
84674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
85674ae819SStefano Zampini 
86674ae819SStefano Zampini   PetscFunctionBegin;
87674ae819SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr);
88674ae819SStefano Zampini   if (pcbddc->local_auxmat1) {
89674ae819SStefano Zampini     ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr);
90674ae819SStefano Zampini     ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr);
91674ae819SStefano Zampini   }
92674ae819SStefano Zampini   PetscFunctionReturn(0);
93674ae819SStefano Zampini }
94674ae819SStefano Zampini 
95674ae819SStefano Zampini #undef __FUNCT__
96674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner"
97674ae819SStefano Zampini PetscErrorCode  PCBDDCApplyInterfacePreconditioner(PC pc)
98674ae819SStefano Zampini {
99674ae819SStefano Zampini   PetscErrorCode ierr;
100674ae819SStefano Zampini   PC_BDDC*        pcbddc = (PC_BDDC*)(pc->data);
101674ae819SStefano Zampini   PC_IS*            pcis = (PC_IS*)  (pc->data);
102674ae819SStefano Zampini   const PetscScalar zero = 0.0;
103674ae819SStefano Zampini 
104674ae819SStefano Zampini   PetscFunctionBegin;
10515aaf578SStefano Zampini   /* Application of PHI^T (or PSI^T)  */
10615aaf578SStefano Zampini   if (pcbddc->coarse_psi_B) {
10715aaf578SStefano Zampini     ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
10815aaf578SStefano Zampini     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
10915aaf578SStefano Zampini   } else {
110674ae819SStefano Zampini     ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
111674ae819SStefano Zampini     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
11215aaf578SStefano Zampini   }
113674ae819SStefano Zampini   /* Scatter data of coarse_rhs */
114674ae819SStefano Zampini   if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); }
115674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
116674ae819SStefano Zampini 
117674ae819SStefano Zampini   /* Local solution on R nodes */
118674ae819SStefano Zampini   ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr);
119674ae819SStefano Zampini   ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
120674ae819SStefano Zampini   ierr = VecScatterEnd  (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
121674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) {
122674ae819SStefano Zampini     ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
123674ae819SStefano Zampini     ierr = VecScatterEnd  (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
124674ae819SStefano Zampini   }
125674ae819SStefano Zampini   ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr);
126674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
127674ae819SStefano Zampini   ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
128674ae819SStefano Zampini   ierr = VecScatterEnd  (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
129674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) {
130674ae819SStefano Zampini     ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
131674ae819SStefano Zampini     ierr = VecScatterEnd  (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132674ae819SStefano Zampini   }
133674ae819SStefano Zampini 
134674ae819SStefano Zampini   /* Coarse solution */
135674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
136674ae819SStefano Zampini   if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */
137674ae819SStefano Zampini     ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr);
138674ae819SStefano Zampini   }
139674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
140674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataEnd  (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
141674ae819SStefano Zampini 
142674ae819SStefano Zampini   /* Sum contributions from two levels */
143674ae819SStefano Zampini   ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
144674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
145674ae819SStefano Zampini   PetscFunctionReturn(0);
146674ae819SStefano Zampini }
147674ae819SStefano Zampini 
148674ae819SStefano Zampini #undef __FUNCT__
149674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataBegin"
150674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
151674ae819SStefano Zampini {
152674ae819SStefano Zampini   PetscErrorCode ierr;
153674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
154674ae819SStefano Zampini 
155674ae819SStefano Zampini   PetscFunctionBegin;
156674ae819SStefano Zampini   switch (pcbddc->coarse_communications_type) {
157674ae819SStefano Zampini     case SCATTERS_BDDC:
158674ae819SStefano Zampini       ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
159674ae819SStefano Zampini       break;
160674ae819SStefano Zampini     case GATHERS_BDDC:
161674ae819SStefano Zampini       break;
162674ae819SStefano Zampini   }
163674ae819SStefano Zampini   PetscFunctionReturn(0);
164674ae819SStefano Zampini }
165674ae819SStefano Zampini 
166674ae819SStefano Zampini #undef __FUNCT__
167674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataEnd"
168674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
169674ae819SStefano Zampini {
170674ae819SStefano Zampini   PetscErrorCode ierr;
171674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
172674ae819SStefano Zampini   PetscScalar*   array_to;
173674ae819SStefano Zampini   PetscScalar*   array_from;
174674ae819SStefano Zampini   MPI_Comm       comm;
175674ae819SStefano Zampini   PetscInt       i;
176674ae819SStefano Zampini 
177674ae819SStefano Zampini   PetscFunctionBegin;
178674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
179674ae819SStefano Zampini   switch (pcbddc->coarse_communications_type) {
180674ae819SStefano Zampini     case SCATTERS_BDDC:
181674ae819SStefano Zampini       ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
182674ae819SStefano Zampini       break;
183674ae819SStefano Zampini     case GATHERS_BDDC:
184674ae819SStefano Zampini       if (vec_from) {
185674ae819SStefano Zampini         ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr);
186674ae819SStefano Zampini       }
187674ae819SStefano Zampini       if (vec_to) {
188674ae819SStefano Zampini         ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr);
189674ae819SStefano Zampini       }
190674ae819SStefano Zampini       switch(pcbddc->coarse_problem_type){
191674ae819SStefano Zampini         case SEQUENTIAL_BDDC:
192674ae819SStefano Zampini           if (smode == SCATTER_FORWARD) {
193674ae819SStefano 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);
194674ae819SStefano Zampini             if (vec_to) {
195674ae819SStefano Zampini               if (imode == ADD_VALUES) {
196674ae819SStefano Zampini                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
197674ae819SStefano Zampini                   array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
198674ae819SStefano Zampini                 }
199674ae819SStefano Zampini               } else {
200674ae819SStefano Zampini                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
201674ae819SStefano Zampini                   array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
202674ae819SStefano Zampini                 }
203674ae819SStefano Zampini               }
204674ae819SStefano Zampini             }
205674ae819SStefano Zampini           } else {
206674ae819SStefano Zampini             if (vec_from) {
207674ae819SStefano Zampini               if (imode == ADD_VALUES) {
208674ae819SStefano Zampini                 MPI_Comm vec_from_comm;
209674ae819SStefano Zampini                 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr);
210674ae819SStefano 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);
211674ae819SStefano Zampini               }
212674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
213674ae819SStefano Zampini                 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]];
214674ae819SStefano Zampini               }
215674ae819SStefano Zampini             }
216674ae819SStefano 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);
217674ae819SStefano Zampini           }
218674ae819SStefano Zampini           break;
219674ae819SStefano Zampini         case REPLICATED_BDDC:
220674ae819SStefano Zampini           if (smode == SCATTER_FORWARD) {
221674ae819SStefano 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);
222674ae819SStefano Zampini             if (imode == ADD_VALUES) {
223674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
224674ae819SStefano Zampini                 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
225674ae819SStefano Zampini               }
226674ae819SStefano Zampini             } else {
227674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
228674ae819SStefano Zampini                 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
229674ae819SStefano Zampini               }
230674ae819SStefano Zampini             }
231674ae819SStefano Zampini           } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */
232674ae819SStefano Zampini             if (imode == ADD_VALUES) {
233674ae819SStefano Zampini               for (i=0;i<pcbddc->local_primal_size;i++) {
234674ae819SStefano Zampini                 array_to[i]+=array_from[pcbddc->local_primal_indices[i]];
235674ae819SStefano Zampini               }
236674ae819SStefano Zampini             } else {
237674ae819SStefano Zampini               for (i=0;i<pcbddc->local_primal_size;i++) {
238674ae819SStefano Zampini                 array_to[i]=array_from[pcbddc->local_primal_indices[i]];
239674ae819SStefano Zampini               }
240674ae819SStefano Zampini             }
241674ae819SStefano Zampini           }
242674ae819SStefano Zampini           break;
243674ae819SStefano Zampini         case MULTILEVEL_BDDC:
244674ae819SStefano Zampini           break;
245674ae819SStefano Zampini         case PARALLEL_BDDC:
246674ae819SStefano Zampini           break;
247674ae819SStefano Zampini       }
248674ae819SStefano Zampini       if (vec_from) {
249674ae819SStefano Zampini         ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr);
250674ae819SStefano Zampini       }
251674ae819SStefano Zampini       if (vec_to) {
252674ae819SStefano Zampini         ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr);
253674ae819SStefano Zampini       }
254674ae819SStefano Zampini       break;
255674ae819SStefano Zampini   }
256674ae819SStefano Zampini   PetscFunctionReturn(0);
257674ae819SStefano Zampini }
258674ae819SStefano Zampini 
259674ae819SStefano Zampini #undef __FUNCT__
260674ae819SStefano Zampini #define __FUNCT__ "PCBDDCConstraintsSetUp"
261674ae819SStefano Zampini PetscErrorCode PCBDDCConstraintsSetUp(PC pc)
262674ae819SStefano Zampini {
263674ae819SStefano Zampini   PetscErrorCode ierr;
264674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)(pc->data);
265674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
266674ae819SStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
267674ae819SStefano Zampini   PetscInt       *nnz,*is_indices;
268674ae819SStefano Zampini   PetscScalar    *temp_quadrature_constraint;
269674ae819SStefano Zampini   PetscInt       *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B;
270674ae819SStefano Zampini   PetscInt       local_primal_size,i,j,k,total_counts,max_size_of_constraint;
271674ae819SStefano Zampini   PetscInt       n_vertices,size_of_constraint;
272674ae819SStefano Zampini   PetscScalar    quad_value;
273674ae819SStefano Zampini   PetscBool      nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true;
274674ae819SStefano Zampini   PetscInt       nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges;
275674ae819SStefano Zampini   IS             *used_IS,ISForVertices,*ISForFaces,*ISForEdges;
276674ae819SStefano Zampini   MatType        impMatType=MATSEQAIJ;
277674ae819SStefano Zampini   PetscBLASInt   Bs,Bt,lwork,lierr;
278674ae819SStefano Zampini   PetscReal      tol=1.0e-8;
279674ae819SStefano Zampini   MatNullSpace   nearnullsp;
280674ae819SStefano Zampini   const Vec      *nearnullvecs;
281674ae819SStefano Zampini   Vec            *localnearnullsp;
282674ae819SStefano Zampini   PetscScalar    *work,*temp_basis,*array_vector,*correlation_mat;
283674ae819SStefano Zampini   PetscReal      *rwork,*singular_vals;
284674ae819SStefano Zampini   PetscBLASInt   Bone=1,*ipiv;
285674ae819SStefano Zampini   Vec            temp_vec;
286674ae819SStefano Zampini   Mat            temp_mat;
287674ae819SStefano Zampini   KSP            temp_ksp;
288674ae819SStefano Zampini   PC             temp_pc;
289674ae819SStefano Zampini   PetscInt       s,start_constraint,dual_dofs;
290674ae819SStefano Zampini   PetscBool      compute_submatrix,useksp=PETSC_FALSE;
291674ae819SStefano Zampini   PetscInt       *aux_primal_permutation,*aux_primal_numbering;
29251b0f6cfSStefano Zampini   PetscBool      boolforchange,*change_basis;
293674ae819SStefano Zampini /* some ugly conditional declarations */
294674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
295674ae819SStefano Zampini   PetscScalar    one=1.0,zero=0.0;
296674ae819SStefano Zampini   PetscInt       ii;
297674ae819SStefano Zampini   PetscScalar    *singular_vectors;
298674ae819SStefano Zampini   PetscBLASInt   *iwork,*ifail;
299674ae819SStefano Zampini   PetscReal      dummy_real,abs_tol;
300674ae819SStefano Zampini   PetscBLASInt   eigs_found;
301674ae819SStefano Zampini #endif
302674ae819SStefano Zampini   PetscBLASInt   dummy_int;
303674ae819SStefano Zampini   PetscScalar    dummy_scalar;
304674ae819SStefano Zampini   PetscBool      used_vertex,get_faces,get_edges,get_vertices;
305674ae819SStefano Zampini 
306674ae819SStefano Zampini   PetscFunctionBegin;
307674ae819SStefano Zampini   /* Get index sets for faces, edges and vertices from graph */
308674ae819SStefano Zampini   get_faces = PETSC_TRUE;
309674ae819SStefano Zampini   get_edges = PETSC_TRUE;
310674ae819SStefano Zampini   get_vertices = PETSC_TRUE;
311674ae819SStefano Zampini   if (pcbddc->vertices_flag) {
312674ae819SStefano Zampini     get_faces = PETSC_FALSE;
313674ae819SStefano Zampini     get_edges = PETSC_FALSE;
314674ae819SStefano Zampini   }
315674ae819SStefano Zampini   if (pcbddc->constraints_flag) {
316674ae819SStefano Zampini     get_vertices = PETSC_FALSE;
317674ae819SStefano Zampini   }
318674ae819SStefano Zampini   if (pcbddc->faces_flag) {
319674ae819SStefano Zampini     get_edges = PETSC_FALSE;
320674ae819SStefano Zampini   }
321674ae819SStefano Zampini   if (pcbddc->edges_flag) {
322674ae819SStefano Zampini     get_faces = PETSC_FALSE;
323674ae819SStefano Zampini   }
324674ae819SStefano Zampini   /* default */
325674ae819SStefano Zampini   if (!get_faces && !get_edges && !get_vertices) {
326674ae819SStefano Zampini     get_vertices = PETSC_TRUE;
327674ae819SStefano Zampini   }
328674ae819SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
329674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
330674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
331674ae819SStefano Zampini     i = 0;
332674ae819SStefano Zampini     if (ISForVertices) {
333674ae819SStefano Zampini       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
334674ae819SStefano Zampini     }
335674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
336674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
33715aaf578SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
338674ae819SStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
339674ae819SStefano Zampini   }
340674ae819SStefano Zampini   /* check if near null space is attached to global mat */
341674ae819SStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
342674ae819SStefano Zampini   if (nearnullsp) {
343674ae819SStefano Zampini     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
344674ae819SStefano Zampini   } else { /* if near null space is not provided it uses constants */
345674ae819SStefano Zampini     nnsp_has_cnst = PETSC_TRUE;
346674ae819SStefano Zampini     use_nnsp_true = PETSC_TRUE;
347674ae819SStefano Zampini   }
348674ae819SStefano Zampini   if (nnsp_has_cnst) {
349674ae819SStefano Zampini     nnsp_addone = 1;
350674ae819SStefano Zampini   }
351674ae819SStefano Zampini   /*
352674ae819SStefano Zampini        Evaluate maximum storage size needed by the procedure
353674ae819SStefano Zampini        - temp_indices will contain start index of each constraint stored as follows
354674ae819SStefano 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
355674ae819SStefano 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
356674ae819SStefano Zampini        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
357674ae819SStefano Zampini                                                                                                                                                          */
358674ae819SStefano Zampini   total_counts = n_ISForFaces+n_ISForEdges;
359674ae819SStefano Zampini   total_counts *= (nnsp_addone+nnsp_size);
360674ae819SStefano Zampini   n_vertices = 0;
361674ae819SStefano Zampini   if (ISForVertices) {
362674ae819SStefano Zampini     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
363674ae819SStefano Zampini   }
364674ae819SStefano Zampini   total_counts += n_vertices;
365674ae819SStefano Zampini   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
366674ae819SStefano Zampini   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
367674ae819SStefano Zampini   total_counts = 0;
368674ae819SStefano Zampini   max_size_of_constraint = 0;
369674ae819SStefano Zampini   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
370674ae819SStefano Zampini     if (i<n_ISForEdges) {
371674ae819SStefano Zampini       used_IS = &ISForEdges[i];
372674ae819SStefano Zampini     } else {
373674ae819SStefano Zampini       used_IS = &ISForFaces[i-n_ISForEdges];
374674ae819SStefano Zampini     }
375674ae819SStefano Zampini     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
376674ae819SStefano Zampini     total_counts += j;
377674ae819SStefano Zampini     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
378674ae819SStefano Zampini   }
379674ae819SStefano Zampini   total_counts *= (nnsp_addone+nnsp_size);
380674ae819SStefano Zampini   total_counts += n_vertices;
381674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
382674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
383674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
384674ae819SStefano Zampini   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
385674ae819SStefano Zampini   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
386674ae819SStefano Zampini   for (i=0;i<pcis->n;i++) {
387674ae819SStefano Zampini     local_to_B[i]=-1;
388674ae819SStefano Zampini   }
389674ae819SStefano Zampini   for (i=0;i<pcis->n_B;i++) {
390674ae819SStefano Zampini     local_to_B[is_indices[i]]=i;
391674ae819SStefano Zampini   }
392674ae819SStefano Zampini   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
393674ae819SStefano Zampini 
394674ae819SStefano Zampini   /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */
395674ae819SStefano Zampini   rwork = 0;
396674ae819SStefano Zampini   work = 0;
397674ae819SStefano Zampini   singular_vals = 0;
398674ae819SStefano Zampini   temp_basis = 0;
399674ae819SStefano Zampini   correlation_mat = 0;
400674ae819SStefano Zampini   if (!pcbddc->use_nnsp_true) {
401674ae819SStefano Zampini     PetscScalar temp_work;
402674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
403674ae819SStefano Zampini     /* POD */
404674ae819SStefano Zampini     PetscInt max_n;
405674ae819SStefano Zampini     max_n = nnsp_addone+nnsp_size;
406674ae819SStefano Zampini     /* using some techniques borrowed from Proper Orthogonal Decomposition */
407674ae819SStefano Zampini     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
408674ae819SStefano Zampini     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr);
409674ae819SStefano Zampini     ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
410674ae819SStefano Zampini     ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
411674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX)
412674ae819SStefano Zampini     ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
413674ae819SStefano Zampini #endif
414674ae819SStefano Zampini     ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr);
415674ae819SStefano Zampini     ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr);
416674ae819SStefano Zampini     /* now we evaluate the optimal workspace using query with lwork=-1 */
417674ae819SStefano Zampini     ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr);
418674ae819SStefano Zampini     lwork=-1;
419674ae819SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
420674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
421674ae819SStefano Zampini     abs_tol=1.e-8;
422674ae819SStefano Zampini /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */
423*f7c40c41SStefano 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));
424674ae819SStefano Zampini #else
425674ae819SStefano Zampini /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */
426674ae819SStefano Zampini /*  LAPACK call is missing here! TODO */
427674ae819SStefano Zampini     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
428674ae819SStefano Zampini #endif
429674ae819SStefano Zampini     if ( lierr ) {
430674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr);
431674ae819SStefano Zampini     }
432674ae819SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
433674ae819SStefano Zampini #else /* on missing GESVD */
434674ae819SStefano Zampini     /* SVD */
435674ae819SStefano Zampini     PetscInt max_n,min_n;
436674ae819SStefano Zampini     max_n = max_size_of_constraint;
437674ae819SStefano Zampini     min_n = nnsp_addone+nnsp_size;
438674ae819SStefano Zampini     if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) {
439674ae819SStefano Zampini       min_n = max_size_of_constraint;
440674ae819SStefano Zampini       max_n = nnsp_addone+nnsp_size;
441674ae819SStefano Zampini     }
442674ae819SStefano Zampini     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
443674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX)
444674ae819SStefano Zampini     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
445674ae819SStefano Zampini #endif
446674ae819SStefano Zampini     /* now we evaluate the optimal workspace using query with lwork=-1 */
447674ae819SStefano Zampini     lwork=-1;
448674ae819SStefano Zampini     ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr);
449674ae819SStefano Zampini     ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr);
450674ae819SStefano Zampini     dummy_int = Bs;
451674ae819SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
452674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
453*f7c40c41SStefano 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));
454674ae819SStefano Zampini #else
455*f7c40c41SStefano 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));
456674ae819SStefano Zampini #endif
457674ae819SStefano Zampini     if ( lierr ) {
458674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr);
459674ae819SStefano Zampini     }
460674ae819SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
461674ae819SStefano Zampini #endif
462674ae819SStefano Zampini     /* Allocate optimal workspace */
463674ae819SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
464674ae819SStefano Zampini     total_counts = (PetscInt)lwork;
465674ae819SStefano Zampini     ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr);
466674ae819SStefano Zampini   }
467674ae819SStefano Zampini   /* get local part of global near null space vectors */
468674ae819SStefano Zampini   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
469674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
470674ae819SStefano Zampini     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
471674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
472674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
473674ae819SStefano Zampini   }
474674ae819SStefano Zampini   /* Now we can loop on constraining sets */
475674ae819SStefano Zampini   total_counts = 0;
476674ae819SStefano Zampini   temp_indices[0] = 0;
477674ae819SStefano Zampini   /* vertices */
478674ae819SStefano Zampini   if (ISForVertices) {
479674ae819SStefano Zampini     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
480674ae819SStefano Zampini     if (nnsp_has_cnst) { /* consider all vertices */
481674ae819SStefano Zampini       for (i=0;i<n_vertices;i++) {
482674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
483674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
484674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
485674ae819SStefano Zampini         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
486674ae819SStefano Zampini         change_basis[total_counts]=PETSC_FALSE;
487674ae819SStefano Zampini         total_counts++;
488674ae819SStefano Zampini       }
489674ae819SStefano Zampini     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
490674ae819SStefano Zampini       for (i=0;i<n_vertices;i++) {
491674ae819SStefano Zampini         used_vertex=PETSC_FALSE;
492674ae819SStefano Zampini         k=0;
493674ae819SStefano Zampini         while (!used_vertex && k<nnsp_size) {
494674ae819SStefano Zampini           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
495674ae819SStefano Zampini           if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) {
496674ae819SStefano Zampini             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
497674ae819SStefano Zampini             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
498674ae819SStefano Zampini             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
499674ae819SStefano Zampini             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
500674ae819SStefano Zampini             change_basis[total_counts]=PETSC_FALSE;
501674ae819SStefano Zampini             total_counts++;
502674ae819SStefano Zampini             used_vertex=PETSC_TRUE;
503674ae819SStefano Zampini           }
504674ae819SStefano Zampini           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
505674ae819SStefano Zampini           k++;
506674ae819SStefano Zampini         }
507674ae819SStefano Zampini       }
508674ae819SStefano Zampini     }
509674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
510674ae819SStefano Zampini     n_vertices = total_counts;
511674ae819SStefano Zampini   }
512674ae819SStefano Zampini   /* edges and faces */
513674ae819SStefano Zampini   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
514674ae819SStefano Zampini     if (i<n_ISForEdges) {
515674ae819SStefano Zampini       used_IS = &ISForEdges[i];
51651b0f6cfSStefano Zampini       boolforchange = pcbddc->use_change_of_basis;
517674ae819SStefano Zampini     } else {
518674ae819SStefano Zampini       used_IS = &ISForFaces[i-n_ISForEdges];
51951b0f6cfSStefano Zampini       boolforchange = pcbddc->use_change_on_faces;
520674ae819SStefano Zampini     }
521674ae819SStefano Zampini     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
522674ae819SStefano Zampini     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
523674ae819SStefano Zampini     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
524674ae819SStefano Zampini     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
52551b0f6cfSStefano Zampini     /* HACK: change of basis should not performed on local periodic nodes */
52651b0f6cfSStefano Zampini     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) {
52751b0f6cfSStefano Zampini       boolforchange = PETSC_FALSE;
52851b0f6cfSStefano Zampini     }
529674ae819SStefano Zampini     if (nnsp_has_cnst) {
530674ae819SStefano Zampini       temp_constraints++;
531674ae819SStefano Zampini       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
532674ae819SStefano Zampini       for (j=0;j<size_of_constraint;j++) {
533674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
534674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
535674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
536674ae819SStefano Zampini       }
537674ae819SStefano Zampini       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
53851b0f6cfSStefano Zampini       change_basis[total_counts]=boolforchange;
539674ae819SStefano Zampini       total_counts++;
540674ae819SStefano Zampini     }
541674ae819SStefano Zampini     for (k=0;k<nnsp_size;k++) {
542674ae819SStefano Zampini       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
543674ae819SStefano Zampini       for (j=0;j<size_of_constraint;j++) {
544674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
545674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
546674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]];
547674ae819SStefano Zampini       }
548674ae819SStefano Zampini       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
549674ae819SStefano Zampini       quad_value = 1.0;
550674ae819SStefano Zampini       if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */
551674ae819SStefano Zampini         ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
552*f7c40c41SStefano Zampini         PetscStackCallBLAS("BLASasum",quad_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone));
553674ae819SStefano Zampini       }
554674ae819SStefano Zampini       if (quad_value > 0.0) { /* keep indices and values */
555674ae819SStefano Zampini         temp_constraints++;
556674ae819SStefano Zampini         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
55751b0f6cfSStefano Zampini         change_basis[total_counts]=boolforchange;
558674ae819SStefano Zampini         total_counts++;
559674ae819SStefano Zampini       }
560674ae819SStefano Zampini     }
561674ae819SStefano Zampini     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
562674ae819SStefano Zampini     /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */
563674ae819SStefano Zampini     if (!use_nnsp_true) {
564674ae819SStefano Zampini       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
565674ae819SStefano Zampini       ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr);
566674ae819SStefano Zampini 
567674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
568674ae819SStefano Zampini       ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr);
569674ae819SStefano Zampini       /* Store upper triangular part of correlation matrix */
570674ae819SStefano Zampini       for (j=0;j<temp_constraints;j++) {
571674ae819SStefano Zampini         for (k=0;k<j+1;k++) {
572*f7c40c41SStefano 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));
573674ae819SStefano Zampini 
574674ae819SStefano Zampini         }
575674ae819SStefano Zampini       }
576674ae819SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
577674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
578674ae819SStefano Zampini /*      LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */
579*f7c40c41SStefano 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));
580674ae819SStefano Zampini #else
581674ae819SStefano Zampini /*  LAPACK call is missing here! TODO */
582674ae819SStefano Zampini       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
583674ae819SStefano Zampini #endif
584674ae819SStefano Zampini       if (lierr) {
585674ae819SStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr);
586674ae819SStefano Zampini       }
587674ae819SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
588674ae819SStefano Zampini       /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */
589674ae819SStefano Zampini       j=0;
590674ae819SStefano Zampini       while (j < Bt && singular_vals[j] < tol) j++;
591674ae819SStefano Zampini       total_counts=total_counts-j;
592674ae819SStefano Zampini       if (j<temp_constraints) {
593674ae819SStefano Zampini         for (k=j;k<Bt;k++) {
594674ae819SStefano Zampini           singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
595674ae819SStefano Zampini         }
596674ae819SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
597*f7c40c41SStefano 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));
598674ae819SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
599674ae819SStefano Zampini         /* copy POD basis into used quadrature memory */
600674ae819SStefano Zampini         for (k=0;k<Bt-j;k++) {
601674ae819SStefano Zampini           for (ii=0;ii<size_of_constraint;ii++) {
602674ae819SStefano 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];
603674ae819SStefano Zampini           }
604674ae819SStefano Zampini         }
605674ae819SStefano Zampini       }
606674ae819SStefano Zampini 
607674ae819SStefano Zampini #else  /* on missing GESVD */
608674ae819SStefano Zampini       PetscInt min_n = temp_constraints;
609674ae819SStefano Zampini       if (min_n > size_of_constraint) min_n = size_of_constraint;
610674ae819SStefano Zampini       dummy_int = Bs;
611674ae819SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
612674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
613*f7c40c41SStefano 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));
614674ae819SStefano Zampini #else
615*f7c40c41SStefano 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));
616674ae819SStefano Zampini #endif
617674ae819SStefano Zampini       if (lierr) {
618674ae819SStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr);
619674ae819SStefano Zampini       }
620674ae819SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
621674ae819SStefano Zampini       /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */
622674ae819SStefano Zampini       j = 0;
623674ae819SStefano Zampini       while (j < min_n && singular_vals[min_n-j-1] < tol) j++;
624674ae819SStefano Zampini       total_counts = total_counts-(PetscInt)Bt+(min_n-j);
625674ae819SStefano Zampini #endif
626674ae819SStefano Zampini     }
627674ae819SStefano Zampini   }
628674ae819SStefano Zampini   /* free index sets of faces, edges and vertices */
629674ae819SStefano Zampini   for (i=0;i<n_ISForFaces;i++) {
630674ae819SStefano Zampini     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
631674ae819SStefano Zampini   }
632674ae819SStefano Zampini   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
633674ae819SStefano Zampini   for (i=0;i<n_ISForEdges;i++) {
634674ae819SStefano Zampini     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
635674ae819SStefano Zampini   }
636674ae819SStefano Zampini   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
637674ae819SStefano Zampini   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
638674ae819SStefano Zampini 
639674ae819SStefano Zampini   /* set quantities in pcbddc data structure */
640674ae819SStefano Zampini   /* n_vertices defines the number of point primal dofs */
641674ae819SStefano Zampini   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
642674ae819SStefano Zampini   local_primal_size = total_counts;
643674ae819SStefano Zampini   pcbddc->n_vertices = n_vertices;
644674ae819SStefano Zampini   pcbddc->n_constraints = total_counts-n_vertices;
645674ae819SStefano Zampini   pcbddc->local_primal_size = local_primal_size;
646674ae819SStefano Zampini 
647674ae819SStefano Zampini   /* Create constraint matrix */
648674ae819SStefano Zampini   /* The constraint matrix is used to compute the l2g map of primal dofs */
649674ae819SStefano Zampini   /* so we need to set it up properly either with or without change of basis */
650674ae819SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
651674ae819SStefano Zampini   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
652674ae819SStefano Zampini   ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr);
653674ae819SStefano Zampini   /* compute a local numbering of constraints : vertices first then constraints */
654674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
655674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
656674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
657674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr);
658674ae819SStefano Zampini   total_counts=0;
659674ae819SStefano Zampini   /* find vertices: subdomain corners plus dofs with basis changed */
660674ae819SStefano Zampini   for (i=0;i<local_primal_size;i++) {
661674ae819SStefano Zampini     size_of_constraint=temp_indices[i+1]-temp_indices[i];
662674ae819SStefano Zampini     if (change_basis[i] || size_of_constraint == 1) {
663674ae819SStefano Zampini       k=0;
664674ae819SStefano Zampini       while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) {
665674ae819SStefano Zampini         k=k+1;
666674ae819SStefano Zampini       }
667674ae819SStefano Zampini       j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1];
668674ae819SStefano Zampini       array_vector[j] = 1.0;
669674ae819SStefano Zampini       aux_primal_numbering[total_counts]=j;
670674ae819SStefano Zampini       aux_primal_permutation[total_counts]=total_counts;
671674ae819SStefano Zampini       total_counts++;
672674ae819SStefano Zampini     }
673674ae819SStefano Zampini   }
674674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
675674ae819SStefano Zampini   /* permute indices in order to have a sorted set of vertices */
676674ae819SStefano Zampini   ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation);
677674ae819SStefano Zampini   /* nonzero structure */
678674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
679674ae819SStefano Zampini   for (i=0;i<total_counts;i++) {
680674ae819SStefano Zampini     nnz[i]=1;
681674ae819SStefano Zampini   }
682674ae819SStefano Zampini   j=total_counts;
683674ae819SStefano Zampini   for (i=n_vertices;i<local_primal_size;i++) {
684674ae819SStefano Zampini     if (!change_basis[i]) {
685674ae819SStefano Zampini       nnz[j]=temp_indices[i+1]-temp_indices[i];
686674ae819SStefano Zampini       j++;
687674ae819SStefano Zampini     }
688674ae819SStefano Zampini   }
689674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
690674ae819SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
691674ae819SStefano Zampini   /* set values in constraint matrix */
692674ae819SStefano Zampini   for (i=0;i<total_counts;i++) {
693674ae819SStefano Zampini     j = aux_primal_permutation[i];
694674ae819SStefano Zampini     k = aux_primal_numbering[j];
695674ae819SStefano Zampini     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr);
696674ae819SStefano Zampini   }
697674ae819SStefano Zampini   for (i=n_vertices;i<local_primal_size;i++) {
698674ae819SStefano Zampini     if (!change_basis[i]) {
699674ae819SStefano Zampini       size_of_constraint=temp_indices[i+1]-temp_indices[i];
700674ae819SStefano 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);
701674ae819SStefano Zampini       total_counts++;
702674ae819SStefano Zampini     }
703674ae819SStefano Zampini   }
704674ae819SStefano Zampini   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
705674ae819SStefano Zampini   ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr);
706674ae819SStefano Zampini   /* assembling */
707674ae819SStefano Zampini   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
708674ae819SStefano Zampini   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
709674ae819SStefano Zampini 
710674ae819SStefano Zampini   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
711674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
712674ae819SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
713674ae819SStefano Zampini     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
714674ae819SStefano Zampini     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
715674ae819SStefano Zampini     /* work arrays */
716674ae819SStefano Zampini     /* we need to reuse these arrays, so we free them */
717674ae819SStefano Zampini     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
718674ae819SStefano Zampini     ierr = PetscFree(work);CHKERRQ(ierr);
719674ae819SStefano Zampini     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
720674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
721674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr);
722674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr);
723674ae819SStefano Zampini     for (i=0;i<pcis->n_B;i++) {
724674ae819SStefano Zampini       nnz[i]=1;
725674ae819SStefano Zampini     }
726674ae819SStefano Zampini     /* Overestimated nonzeros per row */
727674ae819SStefano Zampini     k=1;
728674ae819SStefano Zampini     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
729674ae819SStefano Zampini       if (change_basis[i]) {
730674ae819SStefano Zampini         size_of_constraint = temp_indices[i+1]-temp_indices[i];
731674ae819SStefano Zampini         if (k < size_of_constraint) {
732674ae819SStefano Zampini           k = size_of_constraint;
733674ae819SStefano Zampini         }
734674ae819SStefano Zampini         for (j=0;j<size_of_constraint;j++) {
735674ae819SStefano Zampini           nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
736674ae819SStefano Zampini         }
737674ae819SStefano Zampini       }
738674ae819SStefano Zampini     }
739674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
740674ae819SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
741674ae819SStefano Zampini     /* Temporary array to store indices */
742674ae819SStefano Zampini     ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr);
743674ae819SStefano Zampini     /* Set initial identity in the matrix */
744674ae819SStefano Zampini     for (i=0;i<pcis->n_B;i++) {
745674ae819SStefano Zampini       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
746674ae819SStefano Zampini     }
747674ae819SStefano Zampini     /* Now we loop on the constraints which need a change of basis */
748674ae819SStefano Zampini     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
749674ae819SStefano Zampini     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */
750674ae819SStefano Zampini     temp_constraints = 0;
751674ae819SStefano Zampini     if (pcbddc->n_vertices < local_primal_size) {
752674ae819SStefano Zampini       temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]];
753674ae819SStefano Zampini     }
754674ae819SStefano Zampini     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
755674ae819SStefano Zampini       if (change_basis[i]) {
756674ae819SStefano Zampini         compute_submatrix = PETSC_FALSE;
757674ae819SStefano Zampini         useksp = PETSC_FALSE;
758674ae819SStefano Zampini         if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) {
759674ae819SStefano Zampini           temp_constraints++;
760674ae819SStefano Zampini           if (i == local_primal_size -1 ||  temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) {
761674ae819SStefano Zampini             compute_submatrix = PETSC_TRUE;
762674ae819SStefano Zampini           }
763674ae819SStefano Zampini         }
764674ae819SStefano Zampini         if (compute_submatrix) {
765674ae819SStefano Zampini           if (temp_constraints > 1 || pcbddc->use_nnsp_true) {
766674ae819SStefano Zampini             useksp = PETSC_TRUE;
767674ae819SStefano Zampini           }
768674ae819SStefano Zampini           size_of_constraint = temp_indices[i+1]-temp_indices[i];
769674ae819SStefano Zampini           if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */
770674ae819SStefano Zampini             ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr);
771674ae819SStefano Zampini             ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr);
772674ae819SStefano Zampini             ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr);
773674ae819SStefano Zampini             ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr);
774674ae819SStefano Zampini           }
775674ae819SStefano Zampini           /* First _size_of_constraint-temp_constraints_ columns */
776674ae819SStefano Zampini           dual_dofs = size_of_constraint-temp_constraints;
777674ae819SStefano Zampini           start_constraint = i+1-temp_constraints;
778674ae819SStefano Zampini           for (s=0;s<dual_dofs;s++) {
779674ae819SStefano Zampini             is_indices[0] = s;
780674ae819SStefano Zampini             for (j=0;j<temp_constraints;j++) {
781674ae819SStefano Zampini               for (k=0;k<temp_constraints;k++) {
782674ae819SStefano Zampini                 temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1];
783674ae819SStefano Zampini               }
784674ae819SStefano Zampini               work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s];
785674ae819SStefano Zampini               is_indices[j+1]=s+j+1;
786674ae819SStefano Zampini             }
787674ae819SStefano Zampini             Bt = temp_constraints;
788674ae819SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
789*f7c40c41SStefano Zampini             PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr));
790674ae819SStefano Zampini             if ( lierr ) {
791674ae819SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr);
792674ae819SStefano Zampini             }
793674ae819SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
794674ae819SStefano Zampini             j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s];
795674ae819SStefano 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);
796674ae819SStefano Zampini             if (useksp) {
797674ae819SStefano Zampini               /* temp mat with transposed rows and columns */
798674ae819SStefano Zampini               ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr);
799674ae819SStefano Zampini               ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr);
800674ae819SStefano Zampini             }
801674ae819SStefano Zampini           }
802674ae819SStefano Zampini           if (useksp) {
803674ae819SStefano Zampini             /* last rows of temp_mat */
804674ae819SStefano Zampini             for (j=0;j<size_of_constraint;j++) {
805674ae819SStefano Zampini               is_indices[j] = j;
806674ae819SStefano Zampini             }
807674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
808674ae819SStefano Zampini               k = s + dual_dofs;
809674ae819SStefano Zampini               ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr);
810674ae819SStefano Zampini             }
811674ae819SStefano Zampini             ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
812674ae819SStefano Zampini             ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
813674ae819SStefano Zampini             ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr);
814674ae819SStefano Zampini             ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr);
815674ae819SStefano Zampini             ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr);
816674ae819SStefano Zampini             ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr);
817674ae819SStefano Zampini             ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr);
818674ae819SStefano Zampini             ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr);
819674ae819SStefano Zampini             ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr);
820674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
821674ae819SStefano Zampini               ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr);
822674ae819SStefano Zampini               ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr);
823674ae819SStefano Zampini               ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr);
824674ae819SStefano Zampini               ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr);
825674ae819SStefano Zampini               ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr);
826674ae819SStefano Zampini               ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr);
827674ae819SStefano Zampini               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
828674ae819SStefano Zampini               /* last columns of change of basis matrix associated to new primal dofs */
829674ae819SStefano 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);
830674ae819SStefano Zampini               ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr);
831674ae819SStefano Zampini             }
832674ae819SStefano Zampini             ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
833674ae819SStefano Zampini             ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr);
834674ae819SStefano Zampini             ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
835674ae819SStefano Zampini           } else {
836674ae819SStefano Zampini             /* last columns of change of basis matrix associated to new primal dofs */
837674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
838674ae819SStefano Zampini               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
839674ae819SStefano 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);
840674ae819SStefano Zampini             }
841674ae819SStefano Zampini           }
842674ae819SStefano Zampini           /* prepare for the next cycle */
843674ae819SStefano Zampini           temp_constraints = 0;
844674ae819SStefano Zampini           if (i != local_primal_size -1 ) {
845674ae819SStefano Zampini             temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]];
846674ae819SStefano Zampini           }
847674ae819SStefano Zampini         }
848674ae819SStefano Zampini       }
849674ae819SStefano Zampini     }
850674ae819SStefano Zampini     /* assembling */
851674ae819SStefano Zampini     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
852674ae819SStefano Zampini     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
853674ae819SStefano Zampini     ierr = PetscFree(ipiv);CHKERRQ(ierr);
854674ae819SStefano Zampini     ierr = PetscFree(is_indices);CHKERRQ(ierr);
855674ae819SStefano Zampini   }
856674ae819SStefano Zampini   /* free workspace no longer needed */
857674ae819SStefano Zampini   ierr = PetscFree(rwork);CHKERRQ(ierr);
858674ae819SStefano Zampini   ierr = PetscFree(work);CHKERRQ(ierr);
859674ae819SStefano Zampini   ierr = PetscFree(temp_basis);CHKERRQ(ierr);
860674ae819SStefano Zampini   ierr = PetscFree(singular_vals);CHKERRQ(ierr);
861674ae819SStefano Zampini   ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
862674ae819SStefano Zampini   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
863674ae819SStefano Zampini   ierr = PetscFree(change_basis);CHKERRQ(ierr);
864674ae819SStefano Zampini   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
865674ae819SStefano Zampini   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
866674ae819SStefano Zampini   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
867674ae819SStefano Zampini   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
868674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
869674ae819SStefano Zampini   ierr = PetscFree(iwork);CHKERRQ(ierr);
870674ae819SStefano Zampini   ierr = PetscFree(ifail);CHKERRQ(ierr);
871674ae819SStefano Zampini   ierr = PetscFree(singular_vectors);CHKERRQ(ierr);
872674ae819SStefano Zampini #endif
873674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
874674ae819SStefano Zampini     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
875674ae819SStefano Zampini   }
876674ae819SStefano Zampini   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
877674ae819SStefano Zampini   PetscFunctionReturn(0);
878674ae819SStefano Zampini }
879674ae819SStefano Zampini 
880674ae819SStefano Zampini #undef __FUNCT__
881674ae819SStefano Zampini #define __FUNCT__ "PCBDDCAnalyzeInterface"
882674ae819SStefano Zampini PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
883674ae819SStefano Zampini {
884674ae819SStefano Zampini   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
885674ae819SStefano Zampini   PC_IS       *pcis = (PC_IS*)pc->data;
886674ae819SStefano Zampini   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
887674ae819SStefano Zampini   PetscInt    bs,ierr,i,vertex_size;
888674ae819SStefano Zampini   PetscViewer viewer=pcbddc->dbg_viewer;
889674ae819SStefano Zampini 
890674ae819SStefano Zampini   PetscFunctionBegin;
891674ae819SStefano Zampini   /* Init local Graph struct */
892674ae819SStefano Zampini   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
893674ae819SStefano Zampini 
894575ad6abSStefano Zampini   /* Check validity of the csr graph passed in by the user */
895575ad6abSStefano Zampini   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
896575ad6abSStefano Zampini     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
897575ad6abSStefano Zampini   }
898674ae819SStefano Zampini   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
899674ae819SStefano Zampini   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
900674ae819SStefano Zampini     Mat mat_adj;
901674ae819SStefano Zampini     const PetscInt *xadj,*adjncy;
902674ae819SStefano Zampini     PetscBool flg_row=PETSC_TRUE;
903674ae819SStefano Zampini 
904674ae819SStefano Zampini     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
905674ae819SStefano Zampini     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
906674ae819SStefano Zampini     if (!flg_row) {
907674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
908674ae819SStefano Zampini     }
909674ae819SStefano Zampini     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
910674ae819SStefano Zampini     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
911674ae819SStefano Zampini     if (!flg_row) {
912674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
913674ae819SStefano Zampini     }
914674ae819SStefano Zampini     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
915674ae819SStefano Zampini   }
916674ae819SStefano Zampini 
917674ae819SStefano Zampini   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
918674ae819SStefano Zampini   vertex_size = 1;
919674ae819SStefano Zampini   if (!pcbddc->n_ISForDofs) {
920674ae819SStefano Zampini     IS *custom_ISForDofs;
921674ae819SStefano Zampini 
922674ae819SStefano Zampini     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
923674ae819SStefano Zampini     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
924674ae819SStefano Zampini     for (i=0;i<bs;i++) {
925674ae819SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
926674ae819SStefano Zampini     }
927674ae819SStefano Zampini     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
928674ae819SStefano Zampini     /* remove my references to IS objects */
929674ae819SStefano Zampini     for (i=0;i<bs;i++) {
930674ae819SStefano Zampini       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
931674ae819SStefano Zampini     }
932674ae819SStefano Zampini     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
933674ae819SStefano Zampini   } else { /* mat block size as vertex size (used for elasticity) */
934674ae819SStefano Zampini     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
935674ae819SStefano Zampini   }
936674ae819SStefano Zampini 
937674ae819SStefano Zampini   /* Setup of Graph */
938674ae819SStefano Zampini   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
939674ae819SStefano Zampini 
940674ae819SStefano Zampini   /* Graph's connected components analysis */
941674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
942674ae819SStefano Zampini 
943674ae819SStefano Zampini   /* print some info to stdout */
944674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
945e49050b4SStefano Zampini     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
946674ae819SStefano Zampini   }
947674ae819SStefano Zampini   PetscFunctionReturn(0);
948674ae819SStefano Zampini }
949674ae819SStefano Zampini 
950674ae819SStefano Zampini #undef __FUNCT__
951674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
952674ae819SStefano Zampini PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
953674ae819SStefano Zampini {
954674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
955674ae819SStefano Zampini   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
956674ae819SStefano Zampini   PetscErrorCode ierr;
957674ae819SStefano Zampini 
958674ae819SStefano Zampini   PetscFunctionBegin;
959674ae819SStefano Zampini   n = 0;
960674ae819SStefano Zampini   vertices = 0;
961674ae819SStefano Zampini   if (pcbddc->ConstraintMatrix) {
962674ae819SStefano Zampini     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
963b120a5c6SStefano Zampini     for (i=0;i<local_primal_size;i++) {
964b120a5c6SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
965b120a5c6SStefano Zampini       if (size_of_constraint == 1) n++;
966b120a5c6SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
967b120a5c6SStefano Zampini     }
968b120a5c6SStefano Zampini     ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
969b120a5c6SStefano Zampini     n = 0;
970674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
971674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
972674ae819SStefano Zampini       if (size_of_constraint == 1) {
973674ae819SStefano Zampini         vertices[n++]=row_cmat_indices[0];
974674ae819SStefano Zampini       }
975674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
976674ae819SStefano Zampini     }
977674ae819SStefano Zampini   }
978674ae819SStefano Zampini   *n_vertices = n;
979674ae819SStefano Zampini   *vertices_idx = vertices;
980674ae819SStefano Zampini   PetscFunctionReturn(0);
981674ae819SStefano Zampini }
982674ae819SStefano Zampini 
983674ae819SStefano Zampini #undef __FUNCT__
984674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
985674ae819SStefano Zampini PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
986674ae819SStefano Zampini {
987674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
988674ae819SStefano Zampini   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
989674ae819SStefano Zampini   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
990674ae819SStefano Zampini   PetscBool      *touched;
991674ae819SStefano Zampini   PetscErrorCode ierr;
992674ae819SStefano Zampini 
993674ae819SStefano Zampini   PetscFunctionBegin;
994674ae819SStefano Zampini   n = 0;
995674ae819SStefano Zampini   constraints_index = 0;
996674ae819SStefano Zampini   if (pcbddc->ConstraintMatrix) {
997674ae819SStefano Zampini     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
998674ae819SStefano Zampini     max_size_of_constraint = 0;
999674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1000674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1001674ae819SStefano Zampini       if (size_of_constraint > 1) {
1002674ae819SStefano Zampini         n++;
1003674ae819SStefano Zampini       }
1004674ae819SStefano Zampini       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1005674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1006674ae819SStefano Zampini     }
1007674ae819SStefano Zampini     ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1008674ae819SStefano Zampini     ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1009674ae819SStefano Zampini     ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1010674ae819SStefano Zampini     ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1011674ae819SStefano Zampini     n = 0;
1012674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1013674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1014674ae819SStefano Zampini       if (size_of_constraint > 1) {
1015674ae819SStefano Zampini         ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
1016674ae819SStefano Zampini         min_index = row_cmat_global_indices[0];
1017674ae819SStefano Zampini         min_loc = 0;
1018674ae819SStefano Zampini         for (j=1;j<size_of_constraint;j++) {
1019674ae819SStefano Zampini           /* there can be more than one constraint on a single connected component */
1020674ae819SStefano Zampini           if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1021674ae819SStefano Zampini             min_index = row_cmat_global_indices[j];
1022674ae819SStefano Zampini             min_loc = j;
1023674ae819SStefano Zampini           }
1024674ae819SStefano Zampini         }
1025674ae819SStefano Zampini         touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1026674ae819SStefano Zampini         constraints_index[n++] = row_cmat_indices[min_loc];
1027674ae819SStefano Zampini       }
1028674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1029674ae819SStefano Zampini     }
1030674ae819SStefano Zampini   }
1031674ae819SStefano Zampini   ierr = PetscFree(touched);CHKERRQ(ierr);
1032674ae819SStefano Zampini   ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1033674ae819SStefano Zampini   *n_constraints = n;
1034674ae819SStefano Zampini   *constraints_idx = constraints_index;
1035674ae819SStefano Zampini   PetscFunctionReturn(0);
1036674ae819SStefano Zampini }
1037674ae819SStefano Zampini 
1038674ae819SStefano Zampini /* the next two functions has been adapted from pcis.c */
1039674ae819SStefano Zampini #undef __FUNCT__
1040674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchur"
1041674ae819SStefano Zampini PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1042674ae819SStefano Zampini {
1043674ae819SStefano Zampini   PetscErrorCode ierr;
1044674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
1045674ae819SStefano Zampini 
1046674ae819SStefano Zampini   PetscFunctionBegin;
1047674ae819SStefano Zampini   if (!vec2_B) { vec2_B = v; }
1048674ae819SStefano Zampini   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1049674ae819SStefano Zampini   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1050674ae819SStefano Zampini   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1051674ae819SStefano Zampini   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1052674ae819SStefano Zampini   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1053674ae819SStefano Zampini   PetscFunctionReturn(0);
1054674ae819SStefano Zampini }
1055674ae819SStefano Zampini 
1056674ae819SStefano Zampini #undef __FUNCT__
1057674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchurTranspose"
1058674ae819SStefano Zampini PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1059674ae819SStefano Zampini {
1060674ae819SStefano Zampini   PetscErrorCode ierr;
1061674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
1062674ae819SStefano Zampini 
1063674ae819SStefano Zampini   PetscFunctionBegin;
1064674ae819SStefano Zampini   if (!vec2_B) { vec2_B = v; }
1065674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1066674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1067674ae819SStefano Zampini   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1068674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1069674ae819SStefano Zampini   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1070674ae819SStefano Zampini   PetscFunctionReturn(0);
1071674ae819SStefano Zampini }
1072674ae819SStefano Zampini 
1073674ae819SStefano Zampini #undef __FUNCT__
1074674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSubsetNumbering"
1075674ae819SStefano 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[])
1076674ae819SStefano Zampini {
1077674ae819SStefano Zampini   Vec            local_vec,global_vec;
1078674ae819SStefano Zampini   IS             seqis,paris;
1079674ae819SStefano Zampini   VecScatter     scatter_ctx;
1080674ae819SStefano Zampini   PetscScalar    *array;
1081674ae819SStefano Zampini   PetscInt       *temp_global_dofs;
1082674ae819SStefano Zampini   PetscScalar    globalsum;
1083674ae819SStefano Zampini   PetscInt       i,j,s;
1084674ae819SStefano Zampini   PetscInt       nlocals,first_index,old_index,max_local;
1085674ae819SStefano Zampini   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1086674ae819SStefano Zampini   PetscMPIInt    *dof_sizes,*dof_displs;
1087674ae819SStefano Zampini   PetscBool      first_found;
1088674ae819SStefano Zampini   PetscErrorCode ierr;
1089674ae819SStefano Zampini 
1090674ae819SStefano Zampini   PetscFunctionBegin;
1091674ae819SStefano Zampini   /* mpi buffers */
1092674ae819SStefano Zampini   MPI_Comm_size(comm,&size_prec_comm);
1093674ae819SStefano Zampini   MPI_Comm_rank(comm,&rank_prec_comm);
1094674ae819SStefano Zampini   j = ( !rank_prec_comm ? size_prec_comm : 0);
1095674ae819SStefano Zampini   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1096674ae819SStefano Zampini   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1097674ae819SStefano Zampini   /* get maximum size of subset */
1098674ae819SStefano Zampini   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1099674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1100674ae819SStefano Zampini   max_local = 0;
1101674ae819SStefano Zampini   if (n_local_dofs) {
1102674ae819SStefano Zampini     max_local = temp_global_dofs[0];
1103674ae819SStefano Zampini     for (i=1;i<n_local_dofs;i++) {
1104674ae819SStefano Zampini       if (max_local < temp_global_dofs[i] ) {
1105674ae819SStefano Zampini         max_local = temp_global_dofs[i];
1106674ae819SStefano Zampini       }
1107674ae819SStefano Zampini     }
1108674ae819SStefano Zampini   }
1109674ae819SStefano Zampini   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1110674ae819SStefano Zampini   max_global++;
1111674ae819SStefano Zampini   max_local = 0;
1112674ae819SStefano Zampini   if (n_local_dofs) {
1113674ae819SStefano Zampini     max_local = local_dofs[0];
1114674ae819SStefano Zampini     for (i=1;i<n_local_dofs;i++) {
1115674ae819SStefano Zampini       if (max_local < local_dofs[i] ) {
1116674ae819SStefano Zampini         max_local = local_dofs[i];
1117674ae819SStefano Zampini       }
1118674ae819SStefano Zampini     }
1119674ae819SStefano Zampini   }
1120674ae819SStefano Zampini   max_local++;
1121674ae819SStefano Zampini   /* allocate workspace */
1122674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1123674ae819SStefano Zampini   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1124674ae819SStefano Zampini   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1125674ae819SStefano Zampini   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1126674ae819SStefano Zampini   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1127674ae819SStefano Zampini   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1128674ae819SStefano Zampini   /* create scatter */
1129674ae819SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1130674ae819SStefano Zampini   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1131674ae819SStefano Zampini   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1132674ae819SStefano Zampini   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1133674ae819SStefano Zampini   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1134674ae819SStefano Zampini   /* init array */
1135674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1136674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1137674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1138674ae819SStefano Zampini   if (local_dofs_mult) {
1139674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1140674ae819SStefano Zampini       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1141674ae819SStefano Zampini     }
1142674ae819SStefano Zampini   } else {
1143674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1144674ae819SStefano Zampini       array[local_dofs[i]]=1.0;
1145674ae819SStefano Zampini     }
1146674ae819SStefano Zampini   }
1147674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1148674ae819SStefano Zampini   /* scatter into global vec and get total number of global dofs */
1149674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1150674ae819SStefano Zampini   ierr = VecScatterEnd  (scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1151674ae819SStefano Zampini   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
1152674ae819SStefano Zampini   *n_global_subset = (PetscInt)globalsum;
1153674ae819SStefano Zampini   /* Fill global_vec with cumulative function for global numbering */
1154674ae819SStefano Zampini   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1155674ae819SStefano Zampini   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1156674ae819SStefano Zampini   nlocals = 0;
1157674ae819SStefano Zampini   first_index = -1;
1158674ae819SStefano Zampini   first_found = PETSC_FALSE;
1159674ae819SStefano Zampini   for (i=0;i<s;i++) {
1160674ae819SStefano Zampini     if (!first_found && array[i] > 0.0) {
1161674ae819SStefano Zampini       first_found = PETSC_TRUE;
1162674ae819SStefano Zampini       first_index = i;
1163674ae819SStefano Zampini     }
1164674ae819SStefano Zampini     nlocals += (PetscInt)array[i];
1165674ae819SStefano Zampini   }
1166674ae819SStefano Zampini   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1167674ae819SStefano Zampini   if (!rank_prec_comm) {
1168674ae819SStefano Zampini     dof_displs[0]=0;
1169674ae819SStefano Zampini     for (i=1;i<size_prec_comm;i++) {
1170674ae819SStefano Zampini       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1171674ae819SStefano Zampini     }
1172674ae819SStefano Zampini   }
1173674ae819SStefano Zampini   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1174674ae819SStefano Zampini   if (first_found) {
1175674ae819SStefano Zampini     array[first_index] += (PetscScalar)nlocals;
1176674ae819SStefano Zampini     old_index = first_index;
1177674ae819SStefano Zampini     for (i=first_index+1;i<s;i++) {
1178674ae819SStefano Zampini       if (array[i] > 0.0) {
1179674ae819SStefano Zampini         array[i] += array[old_index];
1180674ae819SStefano Zampini         old_index = i;
1181674ae819SStefano Zampini       }
1182674ae819SStefano Zampini     }
1183674ae819SStefano Zampini   }
1184674ae819SStefano Zampini   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1185674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1186674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1187674ae819SStefano Zampini   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1188674ae819SStefano Zampini   /* get global ordering of local dofs */
1189674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1190674ae819SStefano Zampini   if (local_dofs_mult) {
1191674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1192674ae819SStefano Zampini       temp_global_dofs[i] = (PetscInt)array[local_dofs[i]]-local_dofs_mult[i];
1193674ae819SStefano Zampini     }
1194674ae819SStefano Zampini   } else {
1195674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1196674ae819SStefano Zampini       temp_global_dofs[i] = (PetscInt)array[local_dofs[i]]-1;
1197674ae819SStefano Zampini     }
1198674ae819SStefano Zampini   }
1199674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1200674ae819SStefano Zampini   /* free workspace */
1201674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1202674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1203674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1204674ae819SStefano Zampini   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1205674ae819SStefano Zampini   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1206674ae819SStefano Zampini   /* return pointer to global ordering of local dofs */
1207674ae819SStefano Zampini   *global_numbering_subset = temp_global_dofs;
1208674ae819SStefano Zampini   PetscFunctionReturn(0);
1209674ae819SStefano Zampini }
1210