xref: /petsc/src/ksp/pc/impls/bddc/bddcprivate.c (revision ac78edfc57d8f19d48ae2e55e4bd407c89988fef)
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__
77304d26faSStefano Zampini #define __FUNCT__ "PCBDDCSetUseExactDirichlet"
78304d26faSStefano Zampini PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool use)
79304d26faSStefano Zampini {
80304d26faSStefano Zampini   PC_BDDC  *pcbddc = (PC_BDDC*)pc->data;
81304d26faSStefano Zampini 
82304d26faSStefano Zampini   PetscFunctionBegin;
83304d26faSStefano Zampini   pcbddc->use_exact_dirichlet=use;
84304d26faSStefano Zampini   PetscFunctionReturn(0);
85304d26faSStefano Zampini }
86304d26faSStefano Zampini 
87304d26faSStefano Zampini #undef __FUNCT__
88304d26faSStefano Zampini #define __FUNCT__ "PCBDDCSetUpLocalSolvers"
89304d26faSStefano Zampini PetscErrorCode PCBDDCSetUpLocalSolvers(PC pc, IS is_I_local, IS is_R_local)
90304d26faSStefano Zampini {
91304d26faSStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)pc->data;
92304d26faSStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
93304d26faSStefano Zampini   PC             pc_temp;
94304d26faSStefano Zampini   Mat            A_RR;
95304d26faSStefano Zampini   Vec            vec1,vec2,vec3;
96304d26faSStefano Zampini   MatStructure   matstruct;
97304d26faSStefano Zampini   PetscScalar    m_one = -1.0;
98304d26faSStefano Zampini   PetscReal      value;
99304d26faSStefano Zampini   PetscInt       n_D,n_R,use_exact,use_exact_reduced;
100304d26faSStefano Zampini   PetscErrorCode ierr;
101304d26faSStefano Zampini 
102304d26faSStefano Zampini   PetscFunctionBegin;
103304d26faSStefano Zampini   /* Creating PC contexts for local Dirichlet and Neumann problems */
104304d26faSStefano Zampini   ierr = PCGetOperators(pc,NULL,NULL,&matstruct);CHKERRQ(ierr);
105304d26faSStefano Zampini 
106304d26faSStefano Zampini   /* DIRICHLET PROBLEM */
107*ac78edfcSStefano Zampini   /* Matrix for Dirichlet problem is pcis->A_II */
108304d26faSStefano Zampini   ierr = ISGetSize(is_I_local,&n_D);CHKERRQ(ierr);
109304d26faSStefano Zampini   if (!pcbddc->ksp_D) { /* create object if not yet build */
110304d26faSStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_D);CHKERRQ(ierr);
111304d26faSStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr);
112304d26faSStefano Zampini     /* default */
113304d26faSStefano Zampini     ierr = KSPSetType(pcbddc->ksp_D,KSPPREONLY);CHKERRQ(ierr);
114304d26faSStefano Zampini     ierr = KSPSetOptionsPrefix(pcbddc->ksp_D,"dirichlet_");CHKERRQ(ierr);
115304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
116304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
117304d26faSStefano Zampini     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
118304d26faSStefano Zampini   }
119304d26faSStefano Zampini   ierr = KSPSetOperators(pcbddc->ksp_D,pcis->A_II,pcis->A_II,matstruct);CHKERRQ(ierr);
120304d26faSStefano Zampini   /* Allow user's customization */
121304d26faSStefano Zampini   ierr = KSPSetFromOptions(pcbddc->ksp_D);CHKERRQ(ierr);
122304d26faSStefano Zampini   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
123304d26faSStefano Zampini   if (!n_D) {
124304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
125304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
126304d26faSStefano Zampini   }
127304d26faSStefano Zampini   /* Set Up KSP for Dirichlet problem of BDDC */
128304d26faSStefano Zampini   ierr = KSPSetUp(pcbddc->ksp_D);CHKERRQ(ierr);
129304d26faSStefano Zampini   /* set ksp_D into pcis data */
130304d26faSStefano Zampini   ierr = KSPDestroy(&pcis->ksp_D);CHKERRQ(ierr);
131304d26faSStefano Zampini   ierr = PetscObjectReference((PetscObject)pcbddc->ksp_D);CHKERRQ(ierr);
132304d26faSStefano Zampini   pcis->ksp_D = pcbddc->ksp_D;
133304d26faSStefano Zampini 
134304d26faSStefano Zampini   /* NEUMANN PROBLEM */
135304d26faSStefano Zampini   /* Matrix for Neumann problem is A_RR -> we need to create it */
136304d26faSStefano Zampini   ierr = ISGetSize(is_R_local,&n_R);CHKERRQ(ierr);
137304d26faSStefano Zampini   ierr = MatGetSubMatrix(pcbddc->local_mat,is_R_local,is_R_local,MAT_INITIAL_MATRIX,&A_RR);CHKERRQ(ierr);
138304d26faSStefano Zampini   if (!pcbddc->ksp_R) { /* create object if not yet build */
139304d26faSStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&pcbddc->ksp_R);CHKERRQ(ierr);
140304d26faSStefano Zampini     ierr = PetscObjectIncrementTabLevel((PetscObject)pcbddc->ksp_R,(PetscObject)pc,1);CHKERRQ(ierr);
141304d26faSStefano Zampini     /* default */
142304d26faSStefano Zampini     ierr = KSPSetType(pcbddc->ksp_R,KSPPREONLY);CHKERRQ(ierr);
143304d26faSStefano Zampini     ierr = KSPSetOptionsPrefix(pcbddc->ksp_R,"neumann_");CHKERRQ(ierr);
144304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
145304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
146304d26faSStefano Zampini     ierr = PCFactorSetReuseFill(pc_temp,PETSC_TRUE);CHKERRQ(ierr);
147304d26faSStefano Zampini   }
148304d26faSStefano Zampini   ierr = KSPSetOperators(pcbddc->ksp_R,A_RR,A_RR,matstruct);CHKERRQ(ierr);
149304d26faSStefano Zampini   /* Allow user's customization */
150304d26faSStefano Zampini   ierr = KSPSetFromOptions(pcbddc->ksp_R);CHKERRQ(ierr);
151304d26faSStefano Zampini   /* umfpack interface has a bug when matrix dimension is zero. TODO solve from umfpack interface */
152304d26faSStefano Zampini   if (!n_R) {
153304d26faSStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_R,&pc_temp);CHKERRQ(ierr);
154304d26faSStefano Zampini     ierr = PCSetType(pc_temp,PCNONE);CHKERRQ(ierr);
155304d26faSStefano Zampini   }
156304d26faSStefano Zampini   /* Set Up KSP for Neumann problem of BDDC */
157304d26faSStefano Zampini   ierr = KSPSetUp(pcbddc->ksp_R);CHKERRQ(ierr);
158304d26faSStefano Zampini 
159304d26faSStefano Zampini   /* check Dirichlet and Neumann solvers and adapt them if a nullspace correction is needed */
160304d26faSStefano Zampini 
161304d26faSStefano Zampini   /* Dirichlet */
162304d26faSStefano Zampini   ierr = MatGetVecs(pcis->A_II,&vec1,&vec2);CHKERRQ(ierr);
163304d26faSStefano Zampini   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
164304d26faSStefano Zampini   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
165304d26faSStefano Zampini   ierr = MatMult(pcis->A_II,vec1,vec2);CHKERRQ(ierr);
166304d26faSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_D,vec2,vec3);CHKERRQ(ierr);
167304d26faSStefano Zampini   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
168304d26faSStefano Zampini   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
169304d26faSStefano Zampini   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
170304d26faSStefano Zampini   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
171304d26faSStefano Zampini   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
172304d26faSStefano Zampini   /* need to be adapted? */
173304d26faSStefano Zampini   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
174304d26faSStefano Zampini   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
175304d26faSStefano Zampini   ierr = PCBDDCSetUseExactDirichlet(pc,(PetscBool)use_exact_reduced);CHKERRQ(ierr);
176304d26faSStefano Zampini   /* print info */
177304d26faSStefano Zampini   if (pcbddc->dbg_flag) {
178304d26faSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
179304d26faSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
180304d26faSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Checking solution of Dirichlet and Neumann problems\n");CHKERRQ(ierr);
181304d26faSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for Dirichlet solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
182304d26faSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
183304d26faSStefano Zampini   }
184304d26faSStefano Zampini   if (n_D && pcbddc->NullSpace && !use_exact_reduced && !pcbddc->inexact_prec_type) {
185304d26faSStefano Zampini     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_I_local);CHKERRQ(ierr);
186304d26faSStefano Zampini   }
187304d26faSStefano Zampini 
188304d26faSStefano Zampini   /* Neumann */
189304d26faSStefano Zampini   ierr = MatGetVecs(A_RR,&vec1,&vec2);CHKERRQ(ierr);
190304d26faSStefano Zampini   ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
191304d26faSStefano Zampini   ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
192304d26faSStefano Zampini   ierr = MatMult(A_RR,vec1,vec2);CHKERRQ(ierr);
193304d26faSStefano Zampini   ierr = KSPSolve(pcbddc->ksp_R,vec2,vec3);CHKERRQ(ierr);
194304d26faSStefano Zampini   ierr = VecAXPY(vec3,m_one,vec1);CHKERRQ(ierr);
195304d26faSStefano Zampini   ierr = VecNorm(vec3,NORM_INFINITY,&value);CHKERRQ(ierr);
196304d26faSStefano Zampini   ierr = VecDestroy(&vec1);CHKERRQ(ierr);
197304d26faSStefano Zampini   ierr = VecDestroy(&vec2);CHKERRQ(ierr);
198304d26faSStefano Zampini   ierr = VecDestroy(&vec3);CHKERRQ(ierr);
199304d26faSStefano Zampini   /* need to be adapted? */
200304d26faSStefano Zampini   use_exact = (PetscAbsReal(value) > 1.e-4 ? 0 : 1);
201304d26faSStefano Zampini   if (PetscAbsReal(value) > 1.e-4) use_exact = 0;
202304d26faSStefano Zampini   ierr = MPI_Allreduce(&use_exact,&use_exact_reduced,1,MPIU_INT,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
203304d26faSStefano Zampini   /* print info */
204304d26faSStefano Zampini   if (pcbddc->dbg_flag) {
205304d26faSStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d infinity error for  Neumann  solve = % 1.14e \n",PetscGlobalRank,value);CHKERRQ(ierr);
206304d26faSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
207304d26faSStefano Zampini   }
208304d26faSStefano Zampini   if (n_R && pcbddc->NullSpace && !use_exact_reduced) { /* is it the right logic? */
209304d26faSStefano Zampini     ierr = PCBDDCNullSpaceAssembleCorrection(pc,is_R_local);CHKERRQ(ierr);
210304d26faSStefano Zampini   }
211304d26faSStefano Zampini 
212304d26faSStefano Zampini   /* free Neumann problem's matrix */
213304d26faSStefano Zampini   ierr = MatDestroy(&A_RR);CHKERRQ(ierr);
214304d26faSStefano Zampini   PetscFunctionReturn(0);
215304d26faSStefano Zampini }
216304d26faSStefano Zampini 
217304d26faSStefano Zampini #undef __FUNCT__
218674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSolveSaddlePoint"
219674ae819SStefano Zampini static PetscErrorCode  PCBDDCSolveSaddlePoint(PC pc)
220674ae819SStefano Zampini {
221674ae819SStefano Zampini   PetscErrorCode ierr;
222674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
223674ae819SStefano Zampini 
224674ae819SStefano Zampini   PetscFunctionBegin;
225674ae819SStefano Zampini   ierr = KSPSolve(pcbddc->ksp_R,pcbddc->vec1_R,pcbddc->vec2_R);CHKERRQ(ierr);
226674ae819SStefano Zampini   if (pcbddc->local_auxmat1) {
227674ae819SStefano Zampini     ierr = MatMult(pcbddc->local_auxmat1,pcbddc->vec2_R,pcbddc->vec1_C);CHKERRQ(ierr);
228674ae819SStefano Zampini     ierr = MatMultAdd(pcbddc->local_auxmat2,pcbddc->vec1_C,pcbddc->vec2_R,pcbddc->vec2_R);CHKERRQ(ierr);
229674ae819SStefano Zampini   }
230674ae819SStefano Zampini   PetscFunctionReturn(0);
231674ae819SStefano Zampini }
232674ae819SStefano Zampini 
233674ae819SStefano Zampini #undef __FUNCT__
234674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplyInterfacePreconditioner"
235674ae819SStefano Zampini PetscErrorCode  PCBDDCApplyInterfacePreconditioner(PC pc)
236674ae819SStefano Zampini {
237674ae819SStefano Zampini   PetscErrorCode ierr;
238674ae819SStefano Zampini   PC_BDDC*        pcbddc = (PC_BDDC*)(pc->data);
239674ae819SStefano Zampini   PC_IS*            pcis = (PC_IS*)  (pc->data);
240674ae819SStefano Zampini   const PetscScalar zero = 0.0;
241674ae819SStefano Zampini 
242674ae819SStefano Zampini   PetscFunctionBegin;
24315aaf578SStefano Zampini   /* Application of PHI^T (or PSI^T)  */
24415aaf578SStefano Zampini   if (pcbddc->coarse_psi_B) {
24515aaf578SStefano Zampini     ierr = MatMultTranspose(pcbddc->coarse_psi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
24615aaf578SStefano Zampini     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_psi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
24715aaf578SStefano Zampini   } else {
248674ae819SStefano Zampini     ierr = MatMultTranspose(pcbddc->coarse_phi_B,pcis->vec1_B,pcbddc->vec1_P);CHKERRQ(ierr);
249674ae819SStefano Zampini     if (pcbddc->inexact_prec_type) { ierr = MatMultTransposeAdd(pcbddc->coarse_phi_D,pcis->vec1_D,pcbddc->vec1_P,pcbddc->vec1_P);CHKERRQ(ierr); }
25015aaf578SStefano Zampini   }
251674ae819SStefano Zampini   /* Scatter data of coarse_rhs */
252674ae819SStefano Zampini   if (pcbddc->coarse_rhs) { ierr = VecSet(pcbddc->coarse_rhs,zero);CHKERRQ(ierr); }
253674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
254674ae819SStefano Zampini 
255674ae819SStefano Zampini   /* Local solution on R nodes */
256674ae819SStefano Zampini   ierr = VecSet(pcbddc->vec1_R,zero);CHKERRQ(ierr);
257674ae819SStefano Zampini   ierr = VecScatterBegin(pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
258674ae819SStefano Zampini   ierr = VecScatterEnd  (pcbddc->R_to_B,pcis->vec1_B,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
259674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) {
260674ae819SStefano Zampini     ierr = VecScatterBegin(pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
261674ae819SStefano Zampini     ierr = VecScatterEnd  (pcbddc->R_to_D,pcis->vec1_D,pcbddc->vec1_R,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
262674ae819SStefano Zampini   }
263674ae819SStefano Zampini   ierr = PCBDDCSolveSaddlePoint(pc);CHKERRQ(ierr);
264674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr);
265674ae819SStefano Zampini   ierr = VecScatterBegin(pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
266674ae819SStefano Zampini   ierr = VecScatterEnd  (pcbddc->R_to_B,pcbddc->vec2_R,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
267674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) {
268674ae819SStefano Zampini     ierr = VecScatterBegin(pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
269674ae819SStefano Zampini     ierr = VecScatterEnd  (pcbddc->R_to_D,pcbddc->vec2_R,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
270674ae819SStefano Zampini   }
271674ae819SStefano Zampini 
272674ae819SStefano Zampini   /* Coarse solution */
273674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataEnd(pc,pcbddc->vec1_P,pcbddc->coarse_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
274674ae819SStefano Zampini   if (pcbddc->coarse_rhs) { /* TODO remove null space when doing multilevel */
275674ae819SStefano Zampini     ierr = KSPSolve(pcbddc->coarse_ksp,pcbddc->coarse_rhs,pcbddc->coarse_vec);CHKERRQ(ierr);
276674ae819SStefano Zampini   }
277674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataBegin(pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
278674ae819SStefano Zampini   ierr = PCBDDCScatterCoarseDataEnd  (pc,pcbddc->coarse_vec,pcbddc->vec1_P,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
279674ae819SStefano Zampini 
280674ae819SStefano Zampini   /* Sum contributions from two levels */
281674ae819SStefano Zampini   ierr = MatMultAdd(pcbddc->coarse_phi_B,pcbddc->vec1_P,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
282674ae819SStefano Zampini   if (pcbddc->inexact_prec_type) { ierr = MatMultAdd(pcbddc->coarse_phi_D,pcbddc->vec1_P,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); }
283674ae819SStefano Zampini   PetscFunctionReturn(0);
284674ae819SStefano Zampini }
285674ae819SStefano Zampini 
286674ae819SStefano Zampini #undef __FUNCT__
287674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataBegin"
288674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataBegin(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
289674ae819SStefano Zampini {
290674ae819SStefano Zampini   PetscErrorCode ierr;
291674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
292674ae819SStefano Zampini 
293674ae819SStefano Zampini   PetscFunctionBegin;
294674ae819SStefano Zampini   switch (pcbddc->coarse_communications_type) {
295674ae819SStefano Zampini     case SCATTERS_BDDC:
296674ae819SStefano Zampini       ierr = VecScatterBegin(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
297674ae819SStefano Zampini       break;
298674ae819SStefano Zampini     case GATHERS_BDDC:
299674ae819SStefano Zampini       break;
300674ae819SStefano Zampini   }
301674ae819SStefano Zampini   PetscFunctionReturn(0);
302674ae819SStefano Zampini }
303674ae819SStefano Zampini 
304674ae819SStefano Zampini #undef __FUNCT__
305674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScatterCoarseDataEnd"
306674ae819SStefano Zampini PetscErrorCode PCBDDCScatterCoarseDataEnd(PC pc,Vec vec_from, Vec vec_to, InsertMode imode, ScatterMode smode)
307674ae819SStefano Zampini {
308674ae819SStefano Zampini   PetscErrorCode ierr;
309674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)(pc->data);
310674ae819SStefano Zampini   PetscScalar*   array_to;
311674ae819SStefano Zampini   PetscScalar*   array_from;
312674ae819SStefano Zampini   MPI_Comm       comm;
313674ae819SStefano Zampini   PetscInt       i;
314674ae819SStefano Zampini 
315674ae819SStefano Zampini   PetscFunctionBegin;
316674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
317674ae819SStefano Zampini   switch (pcbddc->coarse_communications_type) {
318674ae819SStefano Zampini     case SCATTERS_BDDC:
319674ae819SStefano Zampini       ierr = VecScatterEnd(pcbddc->coarse_loc_to_glob,vec_from,vec_to,imode,smode);CHKERRQ(ierr);
320674ae819SStefano Zampini       break;
321674ae819SStefano Zampini     case GATHERS_BDDC:
322674ae819SStefano Zampini       if (vec_from) {
323674ae819SStefano Zampini         ierr = VecGetArray(vec_from,&array_from);CHKERRQ(ierr);
324674ae819SStefano Zampini       }
325674ae819SStefano Zampini       if (vec_to) {
326674ae819SStefano Zampini         ierr = VecGetArray(vec_to,&array_to);CHKERRQ(ierr);
327674ae819SStefano Zampini       }
328674ae819SStefano Zampini       switch(pcbddc->coarse_problem_type){
329674ae819SStefano Zampini         case SEQUENTIAL_BDDC:
330674ae819SStefano Zampini           if (smode == SCATTER_FORWARD) {
331674ae819SStefano 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);
332674ae819SStefano Zampini             if (vec_to) {
333674ae819SStefano Zampini               if (imode == ADD_VALUES) {
334674ae819SStefano Zampini                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
335674ae819SStefano Zampini                   array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
336674ae819SStefano Zampini                 }
337674ae819SStefano Zampini               } else {
338674ae819SStefano Zampini                 for (i=0;i<pcbddc->replicated_primal_size;i++) {
339674ae819SStefano Zampini                   array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
340674ae819SStefano Zampini                 }
341674ae819SStefano Zampini               }
342674ae819SStefano Zampini             }
343674ae819SStefano Zampini           } else {
344674ae819SStefano Zampini             if (vec_from) {
345674ae819SStefano Zampini               if (imode == ADD_VALUES) {
346674ae819SStefano Zampini                 MPI_Comm vec_from_comm;
347674ae819SStefano Zampini                 ierr = PetscObjectGetComm((PetscObject)(vec_from),&vec_from_comm);CHKERRQ(ierr);
348674ae819SStefano 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);
349674ae819SStefano Zampini               }
350674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
351674ae819SStefano Zampini                 pcbddc->replicated_local_primal_values[i]=array_from[pcbddc->replicated_local_primal_indices[i]];
352674ae819SStefano Zampini               }
353674ae819SStefano Zampini             }
354674ae819SStefano 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);
355674ae819SStefano Zampini           }
356674ae819SStefano Zampini           break;
357674ae819SStefano Zampini         case REPLICATED_BDDC:
358674ae819SStefano Zampini           if (smode == SCATTER_FORWARD) {
359674ae819SStefano 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);
360674ae819SStefano Zampini             if (imode == ADD_VALUES) {
361674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
362674ae819SStefano Zampini                 array_to[pcbddc->replicated_local_primal_indices[i]]+=pcbddc->replicated_local_primal_values[i];
363674ae819SStefano Zampini               }
364674ae819SStefano Zampini             } else {
365674ae819SStefano Zampini               for (i=0;i<pcbddc->replicated_primal_size;i++) {
366674ae819SStefano Zampini                 array_to[pcbddc->replicated_local_primal_indices[i]]=pcbddc->replicated_local_primal_values[i];
367674ae819SStefano Zampini               }
368674ae819SStefano Zampini             }
369674ae819SStefano Zampini           } else { /* no communications needed for SCATTER_REVERSE since needed data is already present */
370674ae819SStefano Zampini             if (imode == ADD_VALUES) {
371674ae819SStefano Zampini               for (i=0;i<pcbddc->local_primal_size;i++) {
372674ae819SStefano Zampini                 array_to[i]+=array_from[pcbddc->local_primal_indices[i]];
373674ae819SStefano Zampini               }
374674ae819SStefano Zampini             } else {
375674ae819SStefano Zampini               for (i=0;i<pcbddc->local_primal_size;i++) {
376674ae819SStefano Zampini                 array_to[i]=array_from[pcbddc->local_primal_indices[i]];
377674ae819SStefano Zampini               }
378674ae819SStefano Zampini             }
379674ae819SStefano Zampini           }
380674ae819SStefano Zampini           break;
381674ae819SStefano Zampini         case MULTILEVEL_BDDC:
382674ae819SStefano Zampini           break;
383674ae819SStefano Zampini         case PARALLEL_BDDC:
384674ae819SStefano Zampini           break;
385674ae819SStefano Zampini       }
386674ae819SStefano Zampini       if (vec_from) {
387674ae819SStefano Zampini         ierr = VecRestoreArray(vec_from,&array_from);CHKERRQ(ierr);
388674ae819SStefano Zampini       }
389674ae819SStefano Zampini       if (vec_to) {
390674ae819SStefano Zampini         ierr = VecRestoreArray(vec_to,&array_to);CHKERRQ(ierr);
391674ae819SStefano Zampini       }
392674ae819SStefano Zampini       break;
393674ae819SStefano Zampini   }
394674ae819SStefano Zampini   PetscFunctionReturn(0);
395674ae819SStefano Zampini }
396674ae819SStefano Zampini 
397674ae819SStefano Zampini #undef __FUNCT__
398674ae819SStefano Zampini #define __FUNCT__ "PCBDDCConstraintsSetUp"
399674ae819SStefano Zampini PetscErrorCode PCBDDCConstraintsSetUp(PC pc)
400674ae819SStefano Zampini {
401674ae819SStefano Zampini   PetscErrorCode ierr;
402674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)(pc->data);
403674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
404674ae819SStefano Zampini   Mat_IS         *matis = (Mat_IS*)pc->pmat->data;
405674ae819SStefano Zampini   PetscInt       *nnz,*is_indices;
406674ae819SStefano Zampini   PetscScalar    *temp_quadrature_constraint;
407674ae819SStefano Zampini   PetscInt       *temp_indices,*temp_indices_to_constraint,*temp_indices_to_constraint_B,*local_to_B;
408674ae819SStefano Zampini   PetscInt       local_primal_size,i,j,k,total_counts,max_size_of_constraint;
409674ae819SStefano Zampini   PetscInt       n_vertices,size_of_constraint;
4105b08dc53SStefano Zampini   PetscReal      real_value;
411674ae819SStefano Zampini   PetscBool      nnsp_has_cnst=PETSC_FALSE,use_nnsp_true=pcbddc->use_nnsp_true;
412674ae819SStefano Zampini   PetscInt       nnsp_size=0,nnsp_addone=0,temp_constraints,temp_start_ptr,n_ISForFaces,n_ISForEdges;
413674ae819SStefano Zampini   IS             *used_IS,ISForVertices,*ISForFaces,*ISForEdges;
414674ae819SStefano Zampini   MatType        impMatType=MATSEQAIJ;
415674ae819SStefano Zampini   PetscBLASInt   Bs,Bt,lwork,lierr;
416674ae819SStefano Zampini   PetscReal      tol=1.0e-8;
417674ae819SStefano Zampini   MatNullSpace   nearnullsp;
418674ae819SStefano Zampini   const Vec      *nearnullvecs;
419674ae819SStefano Zampini   Vec            *localnearnullsp;
420674ae819SStefano Zampini   PetscScalar    *work,*temp_basis,*array_vector,*correlation_mat;
421674ae819SStefano Zampini   PetscReal      *rwork,*singular_vals;
422674ae819SStefano Zampini   PetscBLASInt   Bone=1,*ipiv;
423674ae819SStefano Zampini   Vec            temp_vec;
424674ae819SStefano Zampini   Mat            temp_mat;
425674ae819SStefano Zampini   KSP            temp_ksp;
426674ae819SStefano Zampini   PC             temp_pc;
427674ae819SStefano Zampini   PetscInt       s,start_constraint,dual_dofs;
428674ae819SStefano Zampini   PetscBool      compute_submatrix,useksp=PETSC_FALSE;
429674ae819SStefano Zampini   PetscInt       *aux_primal_permutation,*aux_primal_numbering;
43051b0f6cfSStefano Zampini   PetscBool      boolforchange,*change_basis;
431674ae819SStefano Zampini /* some ugly conditional declarations */
432674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
433674ae819SStefano Zampini   PetscScalar    one=1.0,zero=0.0;
434674ae819SStefano Zampini   PetscInt       ii;
435674ae819SStefano Zampini   PetscScalar    *singular_vectors;
436674ae819SStefano Zampini   PetscBLASInt   *iwork,*ifail;
437674ae819SStefano Zampini   PetscReal      dummy_real,abs_tol;
438674ae819SStefano Zampini   PetscBLASInt   eigs_found;
439674ae819SStefano Zampini #endif
440674ae819SStefano Zampini   PetscBLASInt   dummy_int;
441674ae819SStefano Zampini   PetscScalar    dummy_scalar;
442674ae819SStefano Zampini   PetscBool      used_vertex,get_faces,get_edges,get_vertices;
443674ae819SStefano Zampini 
444674ae819SStefano Zampini   PetscFunctionBegin;
445674ae819SStefano Zampini   /* Get index sets for faces, edges and vertices from graph */
446674ae819SStefano Zampini   get_faces = PETSC_TRUE;
447674ae819SStefano Zampini   get_edges = PETSC_TRUE;
448674ae819SStefano Zampini   get_vertices = PETSC_TRUE;
449674ae819SStefano Zampini   if (pcbddc->vertices_flag) {
450674ae819SStefano Zampini     get_faces = PETSC_FALSE;
451674ae819SStefano Zampini     get_edges = PETSC_FALSE;
452674ae819SStefano Zampini   }
453674ae819SStefano Zampini   if (pcbddc->constraints_flag) {
454674ae819SStefano Zampini     get_vertices = PETSC_FALSE;
455674ae819SStefano Zampini   }
456674ae819SStefano Zampini   if (pcbddc->faces_flag) {
457674ae819SStefano Zampini     get_edges = PETSC_FALSE;
458674ae819SStefano Zampini   }
459674ae819SStefano Zampini   if (pcbddc->edges_flag) {
460674ae819SStefano Zampini     get_faces = PETSC_FALSE;
461674ae819SStefano Zampini   }
462674ae819SStefano Zampini   /* default */
463674ae819SStefano Zampini   if (!get_faces && !get_edges && !get_vertices) {
464674ae819SStefano Zampini     get_vertices = PETSC_TRUE;
465674ae819SStefano Zampini   }
466674ae819SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,get_faces,get_edges,get_vertices,&n_ISForFaces,&ISForFaces,&n_ISForEdges,&ISForEdges,&ISForVertices);
467674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
468674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"--------------------------------------------------------------\n");CHKERRQ(ierr);
469674ae819SStefano Zampini     i = 0;
470674ae819SStefano Zampini     if (ISForVertices) {
471674ae819SStefano Zampini       ierr = ISGetSize(ISForVertices,&i);CHKERRQ(ierr);
472674ae819SStefano Zampini     }
473674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate vertices\n",PetscGlobalRank,i);CHKERRQ(ierr);
474674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate edges\n",PetscGlobalRank,n_ISForEdges);CHKERRQ(ierr);
47515aaf578SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Subdomain %04d got %02d local candidate faces\n",PetscGlobalRank,n_ISForFaces);CHKERRQ(ierr);
476674ae819SStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
477674ae819SStefano Zampini   }
478674ae819SStefano Zampini   /* check if near null space is attached to global mat */
479674ae819SStefano Zampini   ierr = MatGetNearNullSpace(pc->pmat,&nearnullsp);CHKERRQ(ierr);
480674ae819SStefano Zampini   if (nearnullsp) {
481674ae819SStefano Zampini     ierr = MatNullSpaceGetVecs(nearnullsp,&nnsp_has_cnst,&nnsp_size,&nearnullvecs);CHKERRQ(ierr);
482674ae819SStefano Zampini   } else { /* if near null space is not provided it uses constants */
483674ae819SStefano Zampini     nnsp_has_cnst = PETSC_TRUE;
484674ae819SStefano Zampini     use_nnsp_true = PETSC_TRUE;
485674ae819SStefano Zampini   }
486674ae819SStefano Zampini   if (nnsp_has_cnst) {
487674ae819SStefano Zampini     nnsp_addone = 1;
488674ae819SStefano Zampini   }
489674ae819SStefano Zampini   /*
490674ae819SStefano Zampini        Evaluate maximum storage size needed by the procedure
491674ae819SStefano Zampini        - temp_indices will contain start index of each constraint stored as follows
492674ae819SStefano 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
493674ae819SStefano 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
494674ae819SStefano Zampini        - temp_quadrature_constraint  [temp_indices[i],...,temp[indices[i+1]-1] will contain the scalars representing the constraint itself
495674ae819SStefano Zampini                                                                                                                                                          */
496674ae819SStefano Zampini   total_counts = n_ISForFaces+n_ISForEdges;
497674ae819SStefano Zampini   total_counts *= (nnsp_addone+nnsp_size);
498674ae819SStefano Zampini   n_vertices = 0;
499674ae819SStefano Zampini   if (ISForVertices) {
500674ae819SStefano Zampini     ierr = ISGetSize(ISForVertices,&n_vertices);CHKERRQ(ierr);
501674ae819SStefano Zampini   }
502674ae819SStefano Zampini   total_counts += n_vertices;
503674ae819SStefano Zampini   ierr = PetscMalloc((total_counts+1)*sizeof(PetscInt),&temp_indices);CHKERRQ(ierr);
504674ae819SStefano Zampini   ierr = PetscMalloc((total_counts+1)*sizeof(PetscBool),&change_basis);CHKERRQ(ierr);
505674ae819SStefano Zampini   total_counts = 0;
506674ae819SStefano Zampini   max_size_of_constraint = 0;
507674ae819SStefano Zampini   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
508674ae819SStefano Zampini     if (i<n_ISForEdges) {
509674ae819SStefano Zampini       used_IS = &ISForEdges[i];
510674ae819SStefano Zampini     } else {
511674ae819SStefano Zampini       used_IS = &ISForFaces[i-n_ISForEdges];
512674ae819SStefano Zampini     }
513674ae819SStefano Zampini     ierr = ISGetSize(*used_IS,&j);CHKERRQ(ierr);
514674ae819SStefano Zampini     total_counts += j;
515674ae819SStefano Zampini     max_size_of_constraint = PetscMax(j,max_size_of_constraint);
516674ae819SStefano Zampini   }
517674ae819SStefano Zampini   total_counts *= (nnsp_addone+nnsp_size);
518674ae819SStefano Zampini   total_counts += n_vertices;
519674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&temp_quadrature_constraint);CHKERRQ(ierr);
520674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint);CHKERRQ(ierr);
521674ae819SStefano Zampini   ierr = PetscMalloc(total_counts*sizeof(PetscInt),&temp_indices_to_constraint_B);CHKERRQ(ierr);
522674ae819SStefano Zampini   ierr = PetscMalloc(pcis->n*sizeof(PetscInt),&local_to_B);CHKERRQ(ierr);
523674ae819SStefano Zampini   ierr = ISGetIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
524674ae819SStefano Zampini   for (i=0;i<pcis->n;i++) {
525674ae819SStefano Zampini     local_to_B[i]=-1;
526674ae819SStefano Zampini   }
527674ae819SStefano Zampini   for (i=0;i<pcis->n_B;i++) {
528674ae819SStefano Zampini     local_to_B[is_indices[i]]=i;
529674ae819SStefano Zampini   }
530674ae819SStefano Zampini   ierr = ISRestoreIndices(pcis->is_B_local,(const PetscInt**)&is_indices);CHKERRQ(ierr);
531674ae819SStefano Zampini 
532674ae819SStefano Zampini   /* First we issue queries to allocate optimal workspace for LAPACKgesvd or LAPACKsyev/LAPACKheev */
533674ae819SStefano Zampini   rwork = 0;
534674ae819SStefano Zampini   work = 0;
535674ae819SStefano Zampini   singular_vals = 0;
536674ae819SStefano Zampini   temp_basis = 0;
537674ae819SStefano Zampini   correlation_mat = 0;
538674ae819SStefano Zampini   if (!pcbddc->use_nnsp_true) {
539674ae819SStefano Zampini     PetscScalar temp_work;
540674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
541674ae819SStefano Zampini     /* POD */
542674ae819SStefano Zampini     PetscInt max_n;
543674ae819SStefano Zampini     max_n = nnsp_addone+nnsp_size;
544674ae819SStefano Zampini     /* using some techniques borrowed from Proper Orthogonal Decomposition */
545674ae819SStefano Zampini     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&correlation_mat);CHKERRQ(ierr);
546674ae819SStefano Zampini     ierr = PetscMalloc(max_n*max_n*sizeof(PetscScalar),&singular_vectors);CHKERRQ(ierr);
547674ae819SStefano Zampini     ierr = PetscMalloc(max_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
548674ae819SStefano Zampini     ierr = PetscMalloc(max_size_of_constraint*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
549674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX)
550674ae819SStefano Zampini     ierr = PetscMalloc(3*max_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
551674ae819SStefano Zampini #endif
552674ae819SStefano Zampini     ierr = PetscMalloc(5*max_n*sizeof(PetscBLASInt),&iwork);CHKERRQ(ierr);
553674ae819SStefano Zampini     ierr = PetscMalloc(max_n*sizeof(PetscBLASInt),&ifail);CHKERRQ(ierr);
554674ae819SStefano Zampini     /* now we evaluate the optimal workspace using query with lwork=-1 */
555674ae819SStefano Zampini     ierr = PetscBLASIntCast(max_n,&Bt);CHKERRQ(ierr);
556674ae819SStefano Zampini     lwork=-1;
557674ae819SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
558674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
559674ae819SStefano Zampini     abs_tol=1.e-8;
560674ae819SStefano Zampini /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,&lierr); */
561f7c40c41SStefano 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));
562674ae819SStefano Zampini #else
563674ae819SStefano Zampini /*    LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,&temp_work,&lwork,rwork,&lierr); */
564674ae819SStefano Zampini /*  LAPACK call is missing here! TODO */
565674ae819SStefano Zampini     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
566674ae819SStefano Zampini #endif
567674ae819SStefano Zampini     if ( lierr ) {
568674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYEVX Lapack routine %d",(int)lierr);
569674ae819SStefano Zampini     }
570674ae819SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
571674ae819SStefano Zampini #else /* on missing GESVD */
572674ae819SStefano Zampini     /* SVD */
573674ae819SStefano Zampini     PetscInt max_n,min_n;
574674ae819SStefano Zampini     max_n = max_size_of_constraint;
575674ae819SStefano Zampini     min_n = nnsp_addone+nnsp_size;
576674ae819SStefano Zampini     if (max_size_of_constraint < ( nnsp_addone+nnsp_size ) ) {
577674ae819SStefano Zampini       min_n = max_size_of_constraint;
578674ae819SStefano Zampini       max_n = nnsp_addone+nnsp_size;
579674ae819SStefano Zampini     }
580674ae819SStefano Zampini     ierr = PetscMalloc(min_n*sizeof(PetscReal),&singular_vals);CHKERRQ(ierr);
581674ae819SStefano Zampini #if defined(PETSC_USE_COMPLEX)
582674ae819SStefano Zampini     ierr = PetscMalloc(5*min_n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
583674ae819SStefano Zampini #endif
584674ae819SStefano Zampini     /* now we evaluate the optimal workspace using query with lwork=-1 */
585674ae819SStefano Zampini     lwork=-1;
586674ae819SStefano Zampini     ierr = PetscBLASIntCast(max_n,&Bs);CHKERRQ(ierr);
587674ae819SStefano Zampini     ierr = PetscBLASIntCast(min_n,&Bt);CHKERRQ(ierr);
588674ae819SStefano Zampini     dummy_int = Bs;
589674ae819SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
590674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
591f7c40c41SStefano 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));
592674ae819SStefano Zampini #else
593f7c40c41SStefano 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));
594674ae819SStefano Zampini #endif
595674ae819SStefano Zampini     if ( lierr ) {
596674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SVD Lapack routine %d",(int)lierr);
597674ae819SStefano Zampini     }
598674ae819SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
599674ae819SStefano Zampini #endif
600674ae819SStefano Zampini     /* Allocate optimal workspace */
601674ae819SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(temp_work),&lwork);CHKERRQ(ierr);
602674ae819SStefano Zampini     total_counts = (PetscInt)lwork;
603674ae819SStefano Zampini     ierr = PetscMalloc(total_counts*sizeof(PetscScalar),&work);CHKERRQ(ierr);
604674ae819SStefano Zampini   }
605674ae819SStefano Zampini   /* get local part of global near null space vectors */
606674ae819SStefano Zampini   ierr = PetscMalloc(nnsp_size*sizeof(Vec),&localnearnullsp);CHKERRQ(ierr);
607674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
608674ae819SStefano Zampini     ierr = VecDuplicate(pcis->vec1_N,&localnearnullsp[k]);CHKERRQ(ierr);
609674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
610674ae819SStefano Zampini     ierr = VecScatterEnd(matis->ctx,nearnullvecs[k],localnearnullsp[k],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
611674ae819SStefano Zampini   }
612674ae819SStefano Zampini   /* Now we can loop on constraining sets */
613674ae819SStefano Zampini   total_counts = 0;
614674ae819SStefano Zampini   temp_indices[0] = 0;
615674ae819SStefano Zampini   /* vertices */
616674ae819SStefano Zampini   if (ISForVertices) {
617674ae819SStefano Zampini     ierr = ISGetIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
618674ae819SStefano Zampini     if (nnsp_has_cnst) { /* consider all vertices */
619674ae819SStefano Zampini       for (i=0;i<n_vertices;i++) {
620674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
621674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
622674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
623674ae819SStefano Zampini         temp_indices[total_counts+1]=temp_indices[total_counts]+1;
624674ae819SStefano Zampini         change_basis[total_counts]=PETSC_FALSE;
625674ae819SStefano Zampini         total_counts++;
626674ae819SStefano Zampini       }
627674ae819SStefano Zampini     } else { /* consider vertices for which exist at least a localnearnullsp which is not null there */
628674ae819SStefano Zampini       for (i=0;i<n_vertices;i++) {
629674ae819SStefano Zampini         used_vertex=PETSC_FALSE;
630674ae819SStefano Zampini         k=0;
631674ae819SStefano Zampini         while (!used_vertex && k<nnsp_size) {
632674ae819SStefano Zampini           ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
633674ae819SStefano Zampini           if (PetscAbsScalar(array_vector[is_indices[i]])>0.0) {
634674ae819SStefano Zampini             temp_indices_to_constraint[temp_indices[total_counts]]=is_indices[i];
635674ae819SStefano Zampini             temp_indices_to_constraint_B[temp_indices[total_counts]]=local_to_B[is_indices[i]];
636674ae819SStefano Zampini             temp_quadrature_constraint[temp_indices[total_counts]]=1.0;
637674ae819SStefano Zampini             temp_indices[total_counts+1]=temp_indices[total_counts]+1;
638674ae819SStefano Zampini             change_basis[total_counts]=PETSC_FALSE;
639674ae819SStefano Zampini             total_counts++;
640674ae819SStefano Zampini             used_vertex=PETSC_TRUE;
641674ae819SStefano Zampini           }
642674ae819SStefano Zampini           ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
643674ae819SStefano Zampini           k++;
644674ae819SStefano Zampini         }
645674ae819SStefano Zampini       }
646674ae819SStefano Zampini     }
647674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForVertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
648674ae819SStefano Zampini     n_vertices = total_counts;
649674ae819SStefano Zampini   }
650674ae819SStefano Zampini   /* edges and faces */
651674ae819SStefano Zampini   for (i=0;i<n_ISForEdges+n_ISForFaces;i++) {
652674ae819SStefano Zampini     if (i<n_ISForEdges) {
653674ae819SStefano Zampini       used_IS = &ISForEdges[i];
65451b0f6cfSStefano Zampini       boolforchange = pcbddc->use_change_of_basis;
655674ae819SStefano Zampini     } else {
656674ae819SStefano Zampini       used_IS = &ISForFaces[i-n_ISForEdges];
65751b0f6cfSStefano Zampini       boolforchange = pcbddc->use_change_on_faces;
658674ae819SStefano Zampini     }
659674ae819SStefano Zampini     temp_constraints = 0;          /* zero the number of constraints I have on this conn comp */
660674ae819SStefano Zampini     temp_start_ptr = total_counts; /* need to know the starting index of constraints stored */
661674ae819SStefano Zampini     ierr = ISGetSize(*used_IS,&size_of_constraint);CHKERRQ(ierr);
662674ae819SStefano Zampini     ierr = ISGetIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
66351b0f6cfSStefano Zampini     /* HACK: change of basis should not performed on local periodic nodes */
66451b0f6cfSStefano Zampini     if (pcbddc->mat_graph->mirrors && pcbddc->mat_graph->mirrors[is_indices[0]]) {
66551b0f6cfSStefano Zampini       boolforchange = PETSC_FALSE;
66651b0f6cfSStefano Zampini     }
667674ae819SStefano Zampini     if (nnsp_has_cnst) {
6685b08dc53SStefano Zampini       PetscScalar quad_value;
669674ae819SStefano Zampini       temp_constraints++;
670674ae819SStefano Zampini       quad_value = (PetscScalar)(1.0/PetscSqrtReal((PetscReal)size_of_constraint));
671674ae819SStefano Zampini       for (j=0;j<size_of_constraint;j++) {
672674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
673674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
674674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]+j]=quad_value;
675674ae819SStefano Zampini       }
676674ae819SStefano Zampini       temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
67751b0f6cfSStefano Zampini       change_basis[total_counts]=boolforchange;
678674ae819SStefano Zampini       total_counts++;
679674ae819SStefano Zampini     }
680674ae819SStefano Zampini     for (k=0;k<nnsp_size;k++) {
681674ae819SStefano Zampini       ierr = VecGetArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
682674ae819SStefano Zampini       for (j=0;j<size_of_constraint;j++) {
683674ae819SStefano Zampini         temp_indices_to_constraint[temp_indices[total_counts]+j]=is_indices[j];
684674ae819SStefano Zampini         temp_indices_to_constraint_B[temp_indices[total_counts]+j]=local_to_B[is_indices[j]];
685674ae819SStefano Zampini         temp_quadrature_constraint[temp_indices[total_counts]+j]=array_vector[is_indices[j]];
686674ae819SStefano Zampini       }
687674ae819SStefano Zampini       ierr = VecRestoreArrayRead(localnearnullsp[k],(const PetscScalar**)&array_vector);CHKERRQ(ierr);
6885b08dc53SStefano Zampini       real_value = 1.0;
689674ae819SStefano Zampini       if (use_nnsp_true) { /* check if array is null on the connected component in case use_nnsp_true has been requested */
690674ae819SStefano Zampini         ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
6915b08dc53SStefano Zampini         PetscStackCallBLAS("BLASasum",real_value = BLASasum_(&Bs,&temp_quadrature_constraint[temp_indices[total_counts]],&Bone));
692674ae819SStefano Zampini       }
6935b08dc53SStefano Zampini       if (real_value > 0.0) { /* keep indices and values */
694674ae819SStefano Zampini         temp_constraints++;
695674ae819SStefano Zampini         temp_indices[total_counts+1]=temp_indices[total_counts]+size_of_constraint;  /* store new starting point */
69651b0f6cfSStefano Zampini         change_basis[total_counts]=boolforchange;
697674ae819SStefano Zampini         total_counts++;
698674ae819SStefano Zampini       }
699674ae819SStefano Zampini     }
700674ae819SStefano Zampini     ierr = ISRestoreIndices(*used_IS,(const PetscInt**)&is_indices);CHKERRQ(ierr);
701674ae819SStefano Zampini     /* perform SVD on the constraint if use_nnsp_true has not be requested by the user */
702674ae819SStefano Zampini     if (!use_nnsp_true) {
703674ae819SStefano Zampini       ierr = PetscBLASIntCast(size_of_constraint,&Bs);CHKERRQ(ierr);
704674ae819SStefano Zampini       ierr = PetscBLASIntCast(temp_constraints,&Bt);CHKERRQ(ierr);
705674ae819SStefano Zampini 
706674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
707674ae819SStefano Zampini       ierr = PetscMemzero(correlation_mat,Bt*Bt*sizeof(PetscScalar));CHKERRQ(ierr);
708674ae819SStefano Zampini       /* Store upper triangular part of correlation matrix */
709674ae819SStefano Zampini       for (j=0;j<temp_constraints;j++) {
710674ae819SStefano Zampini         for (k=0;k<j+1;k++) {
7115b08dc53SStefano 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));
712674ae819SStefano Zampini 
713674ae819SStefano Zampini         }
714674ae819SStefano Zampini       }
715674ae819SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
716674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
717674ae819SStefano Zampini /*      LAPACKsyev_("V","U",&Bt,correlation_mat,&Bt,singular_vals,work,&lwork,&lierr); */
718f7c40c41SStefano 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));
719674ae819SStefano Zampini #else
720674ae819SStefano Zampini /*  LAPACK call is missing here! TODO */
721674ae819SStefano Zampini       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for complexes when PETSC_MISSING_GESVD = 1");
722674ae819SStefano Zampini #endif
723674ae819SStefano Zampini       if (lierr) {
724674ae819SStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYEVX Lapack routine %d",(int)lierr);
725674ae819SStefano Zampini       }
726674ae819SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
727674ae819SStefano Zampini       /* retain eigenvalues greater than tol: note that lapack SYEV gives eigs in ascending order */
728674ae819SStefano Zampini       j=0;
729674ae819SStefano Zampini       while (j < Bt && singular_vals[j] < tol) j++;
730674ae819SStefano Zampini       total_counts=total_counts-j;
731674ae819SStefano Zampini       if (j<temp_constraints) {
732674ae819SStefano Zampini         for (k=j;k<Bt;k++) {
733674ae819SStefano Zampini           singular_vals[k]=1.0/PetscSqrtReal(singular_vals[k]);
734674ae819SStefano Zampini         }
735674ae819SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
736f7c40c41SStefano 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));
737674ae819SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
738674ae819SStefano Zampini         /* copy POD basis into used quadrature memory */
739674ae819SStefano Zampini         for (k=0;k<Bt-j;k++) {
740674ae819SStefano Zampini           for (ii=0;ii<size_of_constraint;ii++) {
741674ae819SStefano 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];
742674ae819SStefano Zampini           }
743674ae819SStefano Zampini         }
744674ae819SStefano Zampini       }
745674ae819SStefano Zampini 
746674ae819SStefano Zampini #else  /* on missing GESVD */
747674ae819SStefano Zampini       PetscInt min_n = temp_constraints;
748674ae819SStefano Zampini       if (min_n > size_of_constraint) min_n = size_of_constraint;
749674ae819SStefano Zampini       dummy_int = Bs;
750674ae819SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
751674ae819SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
752f7c40c41SStefano 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));
753674ae819SStefano Zampini #else
754f7c40c41SStefano 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));
755674ae819SStefano Zampini #endif
756674ae819SStefano Zampini       if (lierr) {
757674ae819SStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr);
758674ae819SStefano Zampini       }
759674ae819SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
760674ae819SStefano Zampini       /* retain eigenvalues greater than tol: note that lapack SVD gives eigs in descending order */
761674ae819SStefano Zampini       j = 0;
762674ae819SStefano Zampini       while (j < min_n && singular_vals[min_n-j-1] < tol) j++;
763674ae819SStefano Zampini       total_counts = total_counts-(PetscInt)Bt+(min_n-j);
764674ae819SStefano Zampini #endif
765674ae819SStefano Zampini     }
766674ae819SStefano Zampini   }
767674ae819SStefano Zampini   /* free index sets of faces, edges and vertices */
768674ae819SStefano Zampini   for (i=0;i<n_ISForFaces;i++) {
769674ae819SStefano Zampini     ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
770674ae819SStefano Zampini   }
771674ae819SStefano Zampini   ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
772674ae819SStefano Zampini   for (i=0;i<n_ISForEdges;i++) {
773674ae819SStefano Zampini     ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
774674ae819SStefano Zampini   }
775674ae819SStefano Zampini   ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
776674ae819SStefano Zampini   ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
777674ae819SStefano Zampini 
778674ae819SStefano Zampini   /* set quantities in pcbddc data structure */
779674ae819SStefano Zampini   /* n_vertices defines the number of point primal dofs */
780674ae819SStefano Zampini   /* n_constraints defines the number of averages (they can be point primal dofs if change of basis is requested) */
781674ae819SStefano Zampini   local_primal_size = total_counts;
782674ae819SStefano Zampini   pcbddc->n_vertices = n_vertices;
783674ae819SStefano Zampini   pcbddc->n_constraints = total_counts-n_vertices;
784674ae819SStefano Zampini   pcbddc->local_primal_size = local_primal_size;
785674ae819SStefano Zampini 
786674ae819SStefano Zampini   /* Create constraint matrix */
787674ae819SStefano Zampini   /* The constraint matrix is used to compute the l2g map of primal dofs */
788674ae819SStefano Zampini   /* so we need to set it up properly either with or without change of basis */
789674ae819SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ConstraintMatrix);CHKERRQ(ierr);
790674ae819SStefano Zampini   ierr = MatSetType(pcbddc->ConstraintMatrix,impMatType);CHKERRQ(ierr);
791674ae819SStefano Zampini   ierr = MatSetSizes(pcbddc->ConstraintMatrix,local_primal_size,pcis->n,local_primal_size,pcis->n);CHKERRQ(ierr);
792674ae819SStefano Zampini   /* compute a local numbering of constraints : vertices first then constraints */
793674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
794674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
795674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_numbering);CHKERRQ(ierr);
796674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&aux_primal_permutation);CHKERRQ(ierr);
797674ae819SStefano Zampini   total_counts=0;
798674ae819SStefano Zampini   /* find vertices: subdomain corners plus dofs with basis changed */
799674ae819SStefano Zampini   for (i=0;i<local_primal_size;i++) {
800674ae819SStefano Zampini     size_of_constraint=temp_indices[i+1]-temp_indices[i];
801674ae819SStefano Zampini     if (change_basis[i] || size_of_constraint == 1) {
802674ae819SStefano Zampini       k=0;
803674ae819SStefano Zampini       while(k < size_of_constraint && array_vector[temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1]] != 0.0) {
804674ae819SStefano Zampini         k=k+1;
805674ae819SStefano Zampini       }
806674ae819SStefano Zampini       j=temp_indices_to_constraint[temp_indices[i]+size_of_constraint-k-1];
807674ae819SStefano Zampini       array_vector[j] = 1.0;
808674ae819SStefano Zampini       aux_primal_numbering[total_counts]=j;
809674ae819SStefano Zampini       aux_primal_permutation[total_counts]=total_counts;
810674ae819SStefano Zampini       total_counts++;
811674ae819SStefano Zampini     }
812674ae819SStefano Zampini   }
813674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array_vector);CHKERRQ(ierr);
814674ae819SStefano Zampini   /* permute indices in order to have a sorted set of vertices */
815674ae819SStefano Zampini   ierr = PetscSortIntWithPermutation(total_counts,aux_primal_numbering,aux_primal_permutation);
816674ae819SStefano Zampini   /* nonzero structure */
817674ae819SStefano Zampini   ierr = PetscMalloc(local_primal_size*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
818674ae819SStefano Zampini   for (i=0;i<total_counts;i++) {
819674ae819SStefano Zampini     nnz[i]=1;
820674ae819SStefano Zampini   }
821674ae819SStefano Zampini   j=total_counts;
822674ae819SStefano Zampini   for (i=n_vertices;i<local_primal_size;i++) {
823674ae819SStefano Zampini     if (!change_basis[i]) {
824674ae819SStefano Zampini       nnz[j]=temp_indices[i+1]-temp_indices[i];
825674ae819SStefano Zampini       j++;
826674ae819SStefano Zampini     }
827674ae819SStefano Zampini   }
828674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(pcbddc->ConstraintMatrix,0,nnz);CHKERRQ(ierr);
829674ae819SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
830674ae819SStefano Zampini   /* set values in constraint matrix */
831674ae819SStefano Zampini   for (i=0;i<total_counts;i++) {
832674ae819SStefano Zampini     j = aux_primal_permutation[i];
833674ae819SStefano Zampini     k = aux_primal_numbering[j];
834674ae819SStefano Zampini     ierr = MatSetValue(pcbddc->ConstraintMatrix,i,k,1.0,INSERT_VALUES);CHKERRQ(ierr);
835674ae819SStefano Zampini   }
836674ae819SStefano Zampini   for (i=n_vertices;i<local_primal_size;i++) {
837674ae819SStefano Zampini     if (!change_basis[i]) {
838674ae819SStefano Zampini       size_of_constraint=temp_indices[i+1]-temp_indices[i];
839674ae819SStefano 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);
840674ae819SStefano Zampini       total_counts++;
841674ae819SStefano Zampini     }
842674ae819SStefano Zampini   }
843674ae819SStefano Zampini   ierr = PetscFree(aux_primal_numbering);CHKERRQ(ierr);
844674ae819SStefano Zampini   ierr = PetscFree(aux_primal_permutation);CHKERRQ(ierr);
845674ae819SStefano Zampini   /* assembling */
846674ae819SStefano Zampini   ierr = MatAssemblyBegin(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
847674ae819SStefano Zampini   ierr = MatAssemblyEnd(pcbddc->ConstraintMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
848674ae819SStefano Zampini 
849674ae819SStefano Zampini   /* Create matrix for change of basis. We don't need it in case pcbddc->use_change_of_basis is FALSE */
850674ae819SStefano Zampini   if (pcbddc->use_change_of_basis) {
851674ae819SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr);
852674ae819SStefano Zampini     ierr = MatSetType(pcbddc->ChangeOfBasisMatrix,impMatType);CHKERRQ(ierr);
853674ae819SStefano Zampini     ierr = MatSetSizes(pcbddc->ChangeOfBasisMatrix,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
854674ae819SStefano Zampini     /* work arrays */
855674ae819SStefano Zampini     /* we need to reuse these arrays, so we free them */
856674ae819SStefano Zampini     ierr = PetscFree(temp_basis);CHKERRQ(ierr);
857674ae819SStefano Zampini     ierr = PetscFree(work);CHKERRQ(ierr);
858674ae819SStefano Zampini     ierr = PetscMalloc(pcis->n_B*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
859674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*(nnsp_addone+nnsp_size)*sizeof(PetscScalar),&temp_basis);CHKERRQ(ierr);
860674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscScalar),&work);CHKERRQ(ierr);
861674ae819SStefano Zampini     ierr = PetscMalloc((nnsp_addone+nnsp_size)*sizeof(PetscBLASInt),&ipiv);CHKERRQ(ierr);
862674ae819SStefano Zampini     for (i=0;i<pcis->n_B;i++) {
863674ae819SStefano Zampini       nnz[i]=1;
864674ae819SStefano Zampini     }
865674ae819SStefano Zampini     /* Overestimated nonzeros per row */
866674ae819SStefano Zampini     k=1;
867674ae819SStefano Zampini     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
868674ae819SStefano Zampini       if (change_basis[i]) {
869674ae819SStefano Zampini         size_of_constraint = temp_indices[i+1]-temp_indices[i];
870674ae819SStefano Zampini         if (k < size_of_constraint) {
871674ae819SStefano Zampini           k = size_of_constraint;
872674ae819SStefano Zampini         }
873674ae819SStefano Zampini         for (j=0;j<size_of_constraint;j++) {
874674ae819SStefano Zampini           nnz[temp_indices_to_constraint_B[temp_indices[i]+j]] = size_of_constraint;
875674ae819SStefano Zampini         }
876674ae819SStefano Zampini       }
877674ae819SStefano Zampini     }
878674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(pcbddc->ChangeOfBasisMatrix,0,nnz);CHKERRQ(ierr);
879674ae819SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
880674ae819SStefano Zampini     /* Temporary array to store indices */
881674ae819SStefano Zampini     ierr = PetscMalloc(k*sizeof(PetscInt),&is_indices);CHKERRQ(ierr);
882674ae819SStefano Zampini     /* Set initial identity in the matrix */
883674ae819SStefano Zampini     for (i=0;i<pcis->n_B;i++) {
884674ae819SStefano Zampini       ierr = MatSetValue(pcbddc->ChangeOfBasisMatrix,i,i,1.0,INSERT_VALUES);CHKERRQ(ierr);
885674ae819SStefano Zampini     }
886674ae819SStefano Zampini     /* Now we loop on the constraints which need a change of basis */
887674ae819SStefano Zampini     /* Change of basis matrix is evaluated as the FIRST APPROACH in */
888674ae819SStefano Zampini     /* Klawonn and Widlund, Dual-primal FETI-DP methods for linear elasticity, (6.2.1) */
889674ae819SStefano Zampini     temp_constraints = 0;
890674ae819SStefano Zampini     if (pcbddc->n_vertices < local_primal_size) {
891674ae819SStefano Zampini       temp_start_ptr = temp_indices_to_constraint_B[temp_indices[pcbddc->n_vertices]];
892674ae819SStefano Zampini     }
893674ae819SStefano Zampini     for (i=pcbddc->n_vertices;i<local_primal_size;i++) {
894674ae819SStefano Zampini       if (change_basis[i]) {
895674ae819SStefano Zampini         compute_submatrix = PETSC_FALSE;
896674ae819SStefano Zampini         useksp = PETSC_FALSE;
897674ae819SStefano Zampini         if (temp_start_ptr == temp_indices_to_constraint_B[temp_indices[i]]) {
898674ae819SStefano Zampini           temp_constraints++;
899674ae819SStefano Zampini           if (i == local_primal_size -1 ||  temp_start_ptr != temp_indices_to_constraint_B[temp_indices[i+1]]) {
900674ae819SStefano Zampini             compute_submatrix = PETSC_TRUE;
901674ae819SStefano Zampini           }
902674ae819SStefano Zampini         }
903674ae819SStefano Zampini         if (compute_submatrix) {
904674ae819SStefano Zampini           if (temp_constraints > 1 || pcbddc->use_nnsp_true) {
905674ae819SStefano Zampini             useksp = PETSC_TRUE;
906674ae819SStefano Zampini           }
907674ae819SStefano Zampini           size_of_constraint = temp_indices[i+1]-temp_indices[i];
908674ae819SStefano Zampini           if (useksp) { /* experimental TODO: reuse KSP and MAT instead of creating them each time */
909674ae819SStefano Zampini             ierr = MatCreate(PETSC_COMM_SELF,&temp_mat);CHKERRQ(ierr);
910674ae819SStefano Zampini             ierr = MatSetType(temp_mat,impMatType);CHKERRQ(ierr);
911674ae819SStefano Zampini             ierr = MatSetSizes(temp_mat,size_of_constraint,size_of_constraint,size_of_constraint,size_of_constraint);CHKERRQ(ierr);
912674ae819SStefano Zampini             ierr = MatSeqAIJSetPreallocation(temp_mat,size_of_constraint,NULL);CHKERRQ(ierr);
913674ae819SStefano Zampini           }
914674ae819SStefano Zampini           /* First _size_of_constraint-temp_constraints_ columns */
915674ae819SStefano Zampini           dual_dofs = size_of_constraint-temp_constraints;
916674ae819SStefano Zampini           start_constraint = i+1-temp_constraints;
917674ae819SStefano Zampini           for (s=0;s<dual_dofs;s++) {
918674ae819SStefano Zampini             is_indices[0] = s;
919674ae819SStefano Zampini             for (j=0;j<temp_constraints;j++) {
920674ae819SStefano Zampini               for (k=0;k<temp_constraints;k++) {
921674ae819SStefano Zampini                 temp_basis[j*temp_constraints+k]=temp_quadrature_constraint[temp_indices[start_constraint+k]+s+j+1];
922674ae819SStefano Zampini               }
923674ae819SStefano Zampini               work[j]=-temp_quadrature_constraint[temp_indices[start_constraint+j]+s];
924674ae819SStefano Zampini               is_indices[j+1]=s+j+1;
925674ae819SStefano Zampini             }
926674ae819SStefano Zampini             Bt = temp_constraints;
927674ae819SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
928f7c40c41SStefano Zampini             PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&Bt,&Bone,temp_basis,&Bt,ipiv,work,&Bt,&lierr));
929674ae819SStefano Zampini             if ( lierr ) {
930674ae819SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GESV Lapack routine %d",(int)lierr);
931674ae819SStefano Zampini             }
932674ae819SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
933674ae819SStefano Zampini             j = temp_indices_to_constraint_B[temp_indices[start_constraint]+s];
934674ae819SStefano 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);
935674ae819SStefano Zampini             if (useksp) {
936674ae819SStefano Zampini               /* temp mat with transposed rows and columns */
937674ae819SStefano Zampini               ierr = MatSetValues(temp_mat,1,&s,temp_constraints,&is_indices[1],work,INSERT_VALUES);CHKERRQ(ierr);
938674ae819SStefano Zampini               ierr = MatSetValue(temp_mat,is_indices[0],is_indices[0],1.0,INSERT_VALUES);CHKERRQ(ierr);
939674ae819SStefano Zampini             }
940674ae819SStefano Zampini           }
941674ae819SStefano Zampini           if (useksp) {
942674ae819SStefano Zampini             /* last rows of temp_mat */
943674ae819SStefano Zampini             for (j=0;j<size_of_constraint;j++) {
944674ae819SStefano Zampini               is_indices[j] = j;
945674ae819SStefano Zampini             }
946674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
947674ae819SStefano Zampini               k = s + dual_dofs;
948674ae819SStefano Zampini               ierr = MatSetValues(temp_mat,1,&k,size_of_constraint,is_indices,&temp_quadrature_constraint[temp_indices[start_constraint+s]],INSERT_VALUES);CHKERRQ(ierr);
949674ae819SStefano Zampini             }
950674ae819SStefano Zampini             ierr = MatAssemblyBegin(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
951674ae819SStefano Zampini             ierr = MatAssemblyEnd(temp_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
952674ae819SStefano Zampini             ierr = MatGetVecs(temp_mat,&temp_vec,NULL);CHKERRQ(ierr);
953674ae819SStefano Zampini             ierr = KSPCreate(PETSC_COMM_SELF,&temp_ksp);CHKERRQ(ierr);
954674ae819SStefano Zampini             ierr = KSPSetOperators(temp_ksp,temp_mat,temp_mat,SAME_PRECONDITIONER);CHKERRQ(ierr);
955674ae819SStefano Zampini             ierr = KSPSetType(temp_ksp,KSPPREONLY);CHKERRQ(ierr);
956674ae819SStefano Zampini             ierr = KSPGetPC(temp_ksp,&temp_pc);CHKERRQ(ierr);
957674ae819SStefano Zampini             ierr = PCSetType(temp_pc,PCLU);CHKERRQ(ierr);
958674ae819SStefano Zampini             ierr = KSPSetUp(temp_ksp);CHKERRQ(ierr);
959674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
960674ae819SStefano Zampini               ierr = VecSet(temp_vec,0.0);CHKERRQ(ierr);
961674ae819SStefano Zampini               ierr = VecSetValue(temp_vec,s+dual_dofs,1.0,INSERT_VALUES);CHKERRQ(ierr);
962674ae819SStefano Zampini               ierr = VecAssemblyBegin(temp_vec);CHKERRQ(ierr);
963674ae819SStefano Zampini               ierr = VecAssemblyEnd(temp_vec);CHKERRQ(ierr);
964674ae819SStefano Zampini               ierr = KSPSolve(temp_ksp,temp_vec,temp_vec);CHKERRQ(ierr);
965674ae819SStefano Zampini               ierr = VecGetArray(temp_vec,&array_vector);CHKERRQ(ierr);
966674ae819SStefano Zampini               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
967674ae819SStefano Zampini               /* last columns of change of basis matrix associated to new primal dofs */
968674ae819SStefano 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);
969674ae819SStefano Zampini               ierr = VecRestoreArray(temp_vec,&array_vector);CHKERRQ(ierr);
970674ae819SStefano Zampini             }
971674ae819SStefano Zampini             ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
972674ae819SStefano Zampini             ierr = KSPDestroy(&temp_ksp);CHKERRQ(ierr);
973674ae819SStefano Zampini             ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
974674ae819SStefano Zampini           } else {
975674ae819SStefano Zampini             /* last columns of change of basis matrix associated to new primal dofs */
976674ae819SStefano Zampini             for (s=0;s<temp_constraints;s++) {
977674ae819SStefano Zampini               j = temp_indices_to_constraint_B[temp_indices[start_constraint+s]+size_of_constraint-s-1];
978674ae819SStefano 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);
979674ae819SStefano Zampini             }
980674ae819SStefano Zampini           }
981674ae819SStefano Zampini           /* prepare for the next cycle */
982674ae819SStefano Zampini           temp_constraints = 0;
983674ae819SStefano Zampini           if (i != local_primal_size -1 ) {
984674ae819SStefano Zampini             temp_start_ptr = temp_indices_to_constraint_B[temp_indices[i+1]];
985674ae819SStefano Zampini           }
986674ae819SStefano Zampini         }
987674ae819SStefano Zampini       }
988674ae819SStefano Zampini     }
989674ae819SStefano Zampini     /* assembling */
990674ae819SStefano Zampini     ierr = MatAssemblyBegin(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
991674ae819SStefano Zampini     ierr = MatAssemblyEnd(pcbddc->ChangeOfBasisMatrix,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
992674ae819SStefano Zampini     ierr = PetscFree(ipiv);CHKERRQ(ierr);
993674ae819SStefano Zampini     ierr = PetscFree(is_indices);CHKERRQ(ierr);
994674ae819SStefano Zampini   }
995674ae819SStefano Zampini   /* free workspace no longer needed */
996674ae819SStefano Zampini   ierr = PetscFree(rwork);CHKERRQ(ierr);
997674ae819SStefano Zampini   ierr = PetscFree(work);CHKERRQ(ierr);
998674ae819SStefano Zampini   ierr = PetscFree(temp_basis);CHKERRQ(ierr);
999674ae819SStefano Zampini   ierr = PetscFree(singular_vals);CHKERRQ(ierr);
1000674ae819SStefano Zampini   ierr = PetscFree(correlation_mat);CHKERRQ(ierr);
1001674ae819SStefano Zampini   ierr = PetscFree(temp_indices);CHKERRQ(ierr);
1002674ae819SStefano Zampini   ierr = PetscFree(change_basis);CHKERRQ(ierr);
1003674ae819SStefano Zampini   ierr = PetscFree(temp_indices_to_constraint);CHKERRQ(ierr);
1004674ae819SStefano Zampini   ierr = PetscFree(temp_indices_to_constraint_B);CHKERRQ(ierr);
1005674ae819SStefano Zampini   ierr = PetscFree(local_to_B);CHKERRQ(ierr);
1006674ae819SStefano Zampini   ierr = PetscFree(temp_quadrature_constraint);CHKERRQ(ierr);
1007674ae819SStefano Zampini #if defined(PETSC_MISSING_LAPACK_GESVD)
1008674ae819SStefano Zampini   ierr = PetscFree(iwork);CHKERRQ(ierr);
1009674ae819SStefano Zampini   ierr = PetscFree(ifail);CHKERRQ(ierr);
1010674ae819SStefano Zampini   ierr = PetscFree(singular_vectors);CHKERRQ(ierr);
1011674ae819SStefano Zampini #endif
1012674ae819SStefano Zampini   for (k=0;k<nnsp_size;k++) {
1013674ae819SStefano Zampini     ierr = VecDestroy(&localnearnullsp[k]);CHKERRQ(ierr);
1014674ae819SStefano Zampini   }
1015674ae819SStefano Zampini   ierr = PetscFree(localnearnullsp);CHKERRQ(ierr);
1016674ae819SStefano Zampini   PetscFunctionReturn(0);
1017674ae819SStefano Zampini }
1018674ae819SStefano Zampini 
1019674ae819SStefano Zampini #undef __FUNCT__
1020674ae819SStefano Zampini #define __FUNCT__ "PCBDDCAnalyzeInterface"
1021674ae819SStefano Zampini PetscErrorCode PCBDDCAnalyzeInterface(PC pc)
1022674ae819SStefano Zampini {
1023674ae819SStefano Zampini   PC_BDDC     *pcbddc = (PC_BDDC*)pc->data;
1024674ae819SStefano Zampini   PC_IS       *pcis = (PC_IS*)pc->data;
1025674ae819SStefano Zampini   Mat_IS      *matis  = (Mat_IS*)pc->pmat->data;
1026674ae819SStefano Zampini   PetscInt    bs,ierr,i,vertex_size;
1027674ae819SStefano Zampini   PetscViewer viewer=pcbddc->dbg_viewer;
1028674ae819SStefano Zampini 
1029674ae819SStefano Zampini   PetscFunctionBegin;
1030674ae819SStefano Zampini   /* Init local Graph struct */
1031674ae819SStefano Zampini   ierr = PCBDDCGraphInit(pcbddc->mat_graph,matis->mapping);CHKERRQ(ierr);
1032674ae819SStefano Zampini 
1033575ad6abSStefano Zampini   /* Check validity of the csr graph passed in by the user */
1034575ad6abSStefano Zampini   if (pcbddc->mat_graph->nvtxs_csr != pcbddc->mat_graph->nvtxs) {
1035575ad6abSStefano Zampini     ierr = PCBDDCGraphResetCSR(pcbddc->mat_graph);CHKERRQ(ierr);
1036575ad6abSStefano Zampini   }
1037674ae819SStefano Zampini   /* Set default CSR adjacency of local dofs if not provided by the user with PCBDDCSetLocalAdjacencyGraph */
1038674ae819SStefano Zampini   if (!pcbddc->mat_graph->xadj || !pcbddc->mat_graph->adjncy) {
1039674ae819SStefano Zampini     Mat mat_adj;
1040674ae819SStefano Zampini     const PetscInt *xadj,*adjncy;
1041674ae819SStefano Zampini     PetscBool flg_row=PETSC_TRUE;
1042674ae819SStefano Zampini 
1043674ae819SStefano Zampini     ierr = MatConvert(matis->A,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
1044674ae819SStefano Zampini     ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1045674ae819SStefano Zampini     if (!flg_row) {
1046674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
1047674ae819SStefano Zampini     }
1048674ae819SStefano Zampini     ierr = PCBDDCSetLocalAdjacencyGraph(pc,i,xadj,adjncy,PETSC_COPY_VALUES);CHKERRQ(ierr);
1049674ae819SStefano Zampini     ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&i,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
1050674ae819SStefano Zampini     if (!flg_row) {
1051674ae819SStefano Zampini       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
1052674ae819SStefano Zampini     }
1053674ae819SStefano Zampini     ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
1054674ae819SStefano Zampini   }
1055674ae819SStefano Zampini 
1056674ae819SStefano Zampini   /* Set default dofs' splitting if no information has been provided by the user with PCBDDCSetDofsSplitting */
1057674ae819SStefano Zampini   vertex_size = 1;
1058674ae819SStefano Zampini   if (!pcbddc->n_ISForDofs) {
1059674ae819SStefano Zampini     IS *custom_ISForDofs;
1060674ae819SStefano Zampini 
1061674ae819SStefano Zampini     ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
1062674ae819SStefano Zampini     ierr = PetscMalloc(bs*sizeof(IS),&custom_ISForDofs);CHKERRQ(ierr);
1063674ae819SStefano Zampini     for (i=0;i<bs;i++) {
1064674ae819SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcis->n/bs,i,bs,&custom_ISForDofs[i]);CHKERRQ(ierr);
1065674ae819SStefano Zampini     }
1066674ae819SStefano Zampini     ierr = PCBDDCSetDofsSplitting(pc,bs,custom_ISForDofs);CHKERRQ(ierr);
1067674ae819SStefano Zampini     /* remove my references to IS objects */
1068674ae819SStefano Zampini     for (i=0;i<bs;i++) {
1069674ae819SStefano Zampini       ierr = ISDestroy(&custom_ISForDofs[i]);CHKERRQ(ierr);
1070674ae819SStefano Zampini     }
1071674ae819SStefano Zampini     ierr = PetscFree(custom_ISForDofs);CHKERRQ(ierr);
1072674ae819SStefano Zampini   } else { /* mat block size as vertex size (used for elasticity) */
1073674ae819SStefano Zampini     ierr = MatGetBlockSize(matis->A,&vertex_size);CHKERRQ(ierr);
1074674ae819SStefano Zampini   }
1075674ae819SStefano Zampini 
1076674ae819SStefano Zampini   /* Setup of Graph */
1077674ae819SStefano Zampini   ierr = PCBDDCGraphSetUp(pcbddc->mat_graph,vertex_size,pcbddc->NeumannBoundaries,pcbddc->DirichletBoundaries,pcbddc->n_ISForDofs,pcbddc->ISForDofs,pcbddc->user_primal_vertices);
1078674ae819SStefano Zampini 
1079674ae819SStefano Zampini   /* Graph's connected components analysis */
1080674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponents(pcbddc->mat_graph);CHKERRQ(ierr);
1081674ae819SStefano Zampini 
1082674ae819SStefano Zampini   /* print some info to stdout */
1083674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
1084e49050b4SStefano Zampini     ierr = PCBDDCGraphASCIIView(pcbddc->mat_graph,pcbddc->dbg_flag,viewer);
1085674ae819SStefano Zampini   }
1086674ae819SStefano Zampini   PetscFunctionReturn(0);
1087674ae819SStefano Zampini }
1088674ae819SStefano Zampini 
1089674ae819SStefano Zampini #undef __FUNCT__
1090674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalVerticesLocalIdx"
1091674ae819SStefano Zampini PetscErrorCode  PCBDDCGetPrimalVerticesLocalIdx(PC pc, PetscInt *n_vertices, PetscInt *vertices_idx[])
1092674ae819SStefano Zampini {
1093674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1094674ae819SStefano Zampini   PetscInt       *vertices,*row_cmat_indices,n,i,size_of_constraint,local_primal_size;
1095674ae819SStefano Zampini   PetscErrorCode ierr;
1096674ae819SStefano Zampini 
1097674ae819SStefano Zampini   PetscFunctionBegin;
1098674ae819SStefano Zampini   n = 0;
1099674ae819SStefano Zampini   vertices = 0;
1100674ae819SStefano Zampini   if (pcbddc->ConstraintMatrix) {
1101674ae819SStefano Zampini     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&i);CHKERRQ(ierr);
1102b120a5c6SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1103b120a5c6SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1104b120a5c6SStefano Zampini       if (size_of_constraint == 1) n++;
1105b120a5c6SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1106b120a5c6SStefano Zampini     }
1107b120a5c6SStefano Zampini     ierr = PetscMalloc(n*sizeof(PetscInt),&vertices);CHKERRQ(ierr);
1108b120a5c6SStefano Zampini     n = 0;
1109674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1110674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1111674ae819SStefano Zampini       if (size_of_constraint == 1) {
1112674ae819SStefano Zampini         vertices[n++]=row_cmat_indices[0];
1113674ae819SStefano Zampini       }
1114674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1115674ae819SStefano Zampini     }
1116674ae819SStefano Zampini   }
1117674ae819SStefano Zampini   *n_vertices = n;
1118674ae819SStefano Zampini   *vertices_idx = vertices;
1119674ae819SStefano Zampini   PetscFunctionReturn(0);
1120674ae819SStefano Zampini }
1121674ae819SStefano Zampini 
1122674ae819SStefano Zampini #undef __FUNCT__
1123674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGetPrimalConstraintsLocalIdx"
1124674ae819SStefano Zampini PetscErrorCode  PCBDDCGetPrimalConstraintsLocalIdx(PC pc, PetscInt *n_constraints, PetscInt *constraints_idx[])
1125674ae819SStefano Zampini {
1126674ae819SStefano Zampini   PC_BDDC        *pcbddc = (PC_BDDC*)(pc->data);
1127674ae819SStefano Zampini   PetscInt       *constraints_index,*row_cmat_indices,*row_cmat_global_indices;
1128674ae819SStefano Zampini   PetscInt       n,i,j,size_of_constraint,local_primal_size,local_size,max_size_of_constraint,min_index,min_loc;
1129674ae819SStefano Zampini   PetscBool      *touched;
1130674ae819SStefano Zampini   PetscErrorCode ierr;
1131674ae819SStefano Zampini 
1132674ae819SStefano Zampini   PetscFunctionBegin;
1133674ae819SStefano Zampini   n = 0;
1134674ae819SStefano Zampini   constraints_index = 0;
1135674ae819SStefano Zampini   if (pcbddc->ConstraintMatrix) {
1136674ae819SStefano Zampini     ierr = MatGetSize(pcbddc->ConstraintMatrix,&local_primal_size,&local_size);CHKERRQ(ierr);
1137674ae819SStefano Zampini     max_size_of_constraint = 0;
1138674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1139674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1140674ae819SStefano Zampini       if (size_of_constraint > 1) {
1141674ae819SStefano Zampini         n++;
1142674ae819SStefano Zampini       }
1143674ae819SStefano Zampini       max_size_of_constraint = PetscMax(size_of_constraint,max_size_of_constraint);
1144674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,NULL,NULL);CHKERRQ(ierr);
1145674ae819SStefano Zampini     }
1146674ae819SStefano Zampini     ierr = PetscMalloc(n*sizeof(PetscInt),&constraints_index);CHKERRQ(ierr);
1147674ae819SStefano Zampini     ierr = PetscMalloc(max_size_of_constraint*sizeof(PetscInt),&row_cmat_global_indices);CHKERRQ(ierr);
1148674ae819SStefano Zampini     ierr = PetscMalloc(local_size*sizeof(PetscBool),&touched);CHKERRQ(ierr);
1149674ae819SStefano Zampini     ierr = PetscMemzero(touched,local_size*sizeof(PetscBool));CHKERRQ(ierr);
1150674ae819SStefano Zampini     n = 0;
1151674ae819SStefano Zampini     for (i=0;i<local_primal_size;i++) {
1152674ae819SStefano Zampini       ierr = MatGetRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1153674ae819SStefano Zampini       if (size_of_constraint > 1) {
1154674ae819SStefano Zampini         ierr = ISLocalToGlobalMappingApply(pcbddc->mat_graph->l2gmap,size_of_constraint,row_cmat_indices,row_cmat_global_indices);CHKERRQ(ierr);
1155674ae819SStefano Zampini         min_index = row_cmat_global_indices[0];
1156674ae819SStefano Zampini         min_loc = 0;
1157674ae819SStefano Zampini         for (j=1;j<size_of_constraint;j++) {
1158674ae819SStefano Zampini           /* there can be more than one constraint on a single connected component */
1159674ae819SStefano Zampini           if (min_index > row_cmat_global_indices[j] && !touched[row_cmat_indices[j]]) {
1160674ae819SStefano Zampini             min_index = row_cmat_global_indices[j];
1161674ae819SStefano Zampini             min_loc = j;
1162674ae819SStefano Zampini           }
1163674ae819SStefano Zampini         }
1164674ae819SStefano Zampini         touched[row_cmat_indices[min_loc]] = PETSC_TRUE;
1165674ae819SStefano Zampini         constraints_index[n++] = row_cmat_indices[min_loc];
1166674ae819SStefano Zampini       }
1167674ae819SStefano Zampini       ierr = MatRestoreRow(pcbddc->ConstraintMatrix,i,&size_of_constraint,(const PetscInt**)&row_cmat_indices,NULL);CHKERRQ(ierr);
1168674ae819SStefano Zampini     }
1169674ae819SStefano Zampini   }
1170674ae819SStefano Zampini   ierr = PetscFree(touched);CHKERRQ(ierr);
1171674ae819SStefano Zampini   ierr = PetscFree(row_cmat_global_indices);CHKERRQ(ierr);
1172674ae819SStefano Zampini   *n_constraints = n;
1173674ae819SStefano Zampini   *constraints_idx = constraints_index;
1174674ae819SStefano Zampini   PetscFunctionReturn(0);
1175674ae819SStefano Zampini }
1176674ae819SStefano Zampini 
1177674ae819SStefano Zampini /* the next two functions has been adapted from pcis.c */
1178674ae819SStefano Zampini #undef __FUNCT__
1179674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchur"
1180674ae819SStefano Zampini PetscErrorCode  PCBDDCApplySchur(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1181674ae819SStefano Zampini {
1182674ae819SStefano Zampini   PetscErrorCode ierr;
1183674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
1184674ae819SStefano Zampini 
1185674ae819SStefano Zampini   PetscFunctionBegin;
1186674ae819SStefano Zampini   if (!vec2_B) { vec2_B = v; }
1187674ae819SStefano Zampini   ierr = MatMult(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1188674ae819SStefano Zampini   ierr = MatMult(pcis->A_IB,v,vec1_D);CHKERRQ(ierr);
1189674ae819SStefano Zampini   ierr = KSPSolve(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1190674ae819SStefano Zampini   ierr = MatMult(pcis->A_BI,vec2_D,vec2_B);CHKERRQ(ierr);
1191674ae819SStefano Zampini   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1192674ae819SStefano Zampini   PetscFunctionReturn(0);
1193674ae819SStefano Zampini }
1194674ae819SStefano Zampini 
1195674ae819SStefano Zampini #undef __FUNCT__
1196674ae819SStefano Zampini #define __FUNCT__ "PCBDDCApplySchurTranspose"
1197674ae819SStefano Zampini PetscErrorCode  PCBDDCApplySchurTranspose(PC pc, Vec v, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D)
1198674ae819SStefano Zampini {
1199674ae819SStefano Zampini   PetscErrorCode ierr;
1200674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)(pc->data);
1201674ae819SStefano Zampini 
1202674ae819SStefano Zampini   PetscFunctionBegin;
1203674ae819SStefano Zampini   if (!vec2_B) { vec2_B = v; }
1204674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_BB,v,vec1_B);CHKERRQ(ierr);
1205674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_BI,v,vec1_D);CHKERRQ(ierr);
1206674ae819SStefano Zampini   ierr = KSPSolveTranspose(pcis->ksp_D,vec1_D,vec2_D);CHKERRQ(ierr);
1207674ae819SStefano Zampini   ierr = MatMultTranspose(pcis->A_IB,vec2_D,vec2_B);CHKERRQ(ierr);
1208674ae819SStefano Zampini   ierr = VecAXPY(vec1_B,-1.0,vec2_B);CHKERRQ(ierr);
1209674ae819SStefano Zampini   PetscFunctionReturn(0);
1210674ae819SStefano Zampini }
1211674ae819SStefano Zampini 
1212674ae819SStefano Zampini #undef __FUNCT__
1213674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSubsetNumbering"
1214674ae819SStefano 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[])
1215674ae819SStefano Zampini {
1216674ae819SStefano Zampini   Vec            local_vec,global_vec;
1217674ae819SStefano Zampini   IS             seqis,paris;
1218674ae819SStefano Zampini   VecScatter     scatter_ctx;
1219674ae819SStefano Zampini   PetscScalar    *array;
1220674ae819SStefano Zampini   PetscInt       *temp_global_dofs;
1221674ae819SStefano Zampini   PetscScalar    globalsum;
1222674ae819SStefano Zampini   PetscInt       i,j,s;
1223674ae819SStefano Zampini   PetscInt       nlocals,first_index,old_index,max_local;
1224674ae819SStefano Zampini   PetscMPIInt    rank_prec_comm,size_prec_comm,max_global;
1225674ae819SStefano Zampini   PetscMPIInt    *dof_sizes,*dof_displs;
1226674ae819SStefano Zampini   PetscBool      first_found;
1227674ae819SStefano Zampini   PetscErrorCode ierr;
1228674ae819SStefano Zampini 
1229674ae819SStefano Zampini   PetscFunctionBegin;
1230674ae819SStefano Zampini   /* mpi buffers */
1231674ae819SStefano Zampini   MPI_Comm_size(comm,&size_prec_comm);
1232674ae819SStefano Zampini   MPI_Comm_rank(comm,&rank_prec_comm);
1233674ae819SStefano Zampini   j = ( !rank_prec_comm ? size_prec_comm : 0);
1234674ae819SStefano Zampini   ierr = PetscMalloc(j*sizeof(*dof_sizes),&dof_sizes);CHKERRQ(ierr);
1235674ae819SStefano Zampini   ierr = PetscMalloc(j*sizeof(*dof_displs),&dof_displs);CHKERRQ(ierr);
1236674ae819SStefano Zampini   /* get maximum size of subset */
1237674ae819SStefano Zampini   ierr = PetscMalloc(n_local_dofs*sizeof(PetscInt),&temp_global_dofs);CHKERRQ(ierr);
1238674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApply(l2gmap,n_local_dofs,local_dofs,temp_global_dofs);CHKERRQ(ierr);
1239674ae819SStefano Zampini   max_local = 0;
1240674ae819SStefano Zampini   if (n_local_dofs) {
1241674ae819SStefano Zampini     max_local = temp_global_dofs[0];
1242674ae819SStefano Zampini     for (i=1;i<n_local_dofs;i++) {
1243674ae819SStefano Zampini       if (max_local < temp_global_dofs[i] ) {
1244674ae819SStefano Zampini         max_local = temp_global_dofs[i];
1245674ae819SStefano Zampini       }
1246674ae819SStefano Zampini     }
1247674ae819SStefano Zampini   }
1248674ae819SStefano Zampini   ierr = MPI_Allreduce(&max_local,&max_global,1,MPIU_INT,MPI_MAX,comm);
1249674ae819SStefano Zampini   max_global++;
1250674ae819SStefano Zampini   max_local = 0;
1251674ae819SStefano Zampini   if (n_local_dofs) {
1252674ae819SStefano Zampini     max_local = local_dofs[0];
1253674ae819SStefano Zampini     for (i=1;i<n_local_dofs;i++) {
1254674ae819SStefano Zampini       if (max_local < local_dofs[i] ) {
1255674ae819SStefano Zampini         max_local = local_dofs[i];
1256674ae819SStefano Zampini       }
1257674ae819SStefano Zampini     }
1258674ae819SStefano Zampini   }
1259674ae819SStefano Zampini   max_local++;
1260674ae819SStefano Zampini   /* allocate workspace */
1261674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
1262674ae819SStefano Zampini   ierr = VecSetSizes(local_vec,PETSC_DECIDE,max_local);CHKERRQ(ierr);
1263674ae819SStefano Zampini   ierr = VecSetType(local_vec,VECSEQ);CHKERRQ(ierr);
1264674ae819SStefano Zampini   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
1265674ae819SStefano Zampini   ierr = VecSetSizes(global_vec,PETSC_DECIDE,max_global);CHKERRQ(ierr);
1266674ae819SStefano Zampini   ierr = VecSetType(global_vec,VECMPI);CHKERRQ(ierr);
1267674ae819SStefano Zampini   /* create scatter */
1268674ae819SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,n_local_dofs,local_dofs,PETSC_COPY_VALUES,&seqis);CHKERRQ(ierr);
1269674ae819SStefano Zampini   ierr = ISCreateGeneral(comm,n_local_dofs,temp_global_dofs,PETSC_COPY_VALUES,&paris);CHKERRQ(ierr);
1270674ae819SStefano Zampini   ierr = VecScatterCreate(local_vec,seqis,global_vec,paris,&scatter_ctx);CHKERRQ(ierr);
1271674ae819SStefano Zampini   ierr = ISDestroy(&seqis);CHKERRQ(ierr);
1272674ae819SStefano Zampini   ierr = ISDestroy(&paris);CHKERRQ(ierr);
1273674ae819SStefano Zampini   /* init array */
1274674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
1275674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1276674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1277674ae819SStefano Zampini   if (local_dofs_mult) {
1278674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1279674ae819SStefano Zampini       array[local_dofs[i]]=(PetscScalar)local_dofs_mult[i];
1280674ae819SStefano Zampini     }
1281674ae819SStefano Zampini   } else {
1282674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
1283674ae819SStefano Zampini       array[local_dofs[i]]=1.0;
1284674ae819SStefano Zampini     }
1285674ae819SStefano Zampini   }
1286674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1287674ae819SStefano Zampini   /* scatter into global vec and get total number of global dofs */
1288674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1289674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1290674ae819SStefano Zampini   ierr = VecSum(global_vec,&globalsum);CHKERRQ(ierr);
12915b08dc53SStefano Zampini   *n_global_subset = (PetscInt)PetscRealPart(globalsum);
1292674ae819SStefano Zampini   /* Fill global_vec with cumulative function for global numbering */
1293674ae819SStefano Zampini   ierr = VecGetArray(global_vec,&array);CHKERRQ(ierr);
1294674ae819SStefano Zampini   ierr = VecGetLocalSize(global_vec,&s);CHKERRQ(ierr);
1295674ae819SStefano Zampini   nlocals = 0;
1296674ae819SStefano Zampini   first_index = -1;
1297674ae819SStefano Zampini   first_found = PETSC_FALSE;
1298674ae819SStefano Zampini   for (i=0;i<s;i++) {
12995b08dc53SStefano Zampini     if (!first_found && PetscRealPart(array[i]) > 0.0) {
1300674ae819SStefano Zampini       first_found = PETSC_TRUE;
1301674ae819SStefano Zampini       first_index = i;
1302674ae819SStefano Zampini     }
13035b08dc53SStefano Zampini     nlocals += (PetscInt)PetscRealPart(array[i]);
1304674ae819SStefano Zampini   }
1305674ae819SStefano Zampini   ierr = MPI_Gather(&nlocals,1,MPIU_INT,dof_sizes,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1306674ae819SStefano Zampini   if (!rank_prec_comm) {
1307674ae819SStefano Zampini     dof_displs[0]=0;
1308674ae819SStefano Zampini     for (i=1;i<size_prec_comm;i++) {
1309674ae819SStefano Zampini       dof_displs[i] = dof_displs[i-1]+dof_sizes[i-1];
1310674ae819SStefano Zampini     }
1311674ae819SStefano Zampini   }
1312674ae819SStefano Zampini   ierr = MPI_Scatter(dof_displs,1,MPIU_INT,&nlocals,1,MPIU_INT,0,comm);CHKERRQ(ierr);
1313674ae819SStefano Zampini   if (first_found) {
1314674ae819SStefano Zampini     array[first_index] += (PetscScalar)nlocals;
1315674ae819SStefano Zampini     old_index = first_index;
1316674ae819SStefano Zampini     for (i=first_index+1;i<s;i++) {
13175b08dc53SStefano Zampini       if (PetscRealPart(array[i]) > 0.0) {
1318674ae819SStefano Zampini         array[i] += array[old_index];
1319674ae819SStefano Zampini         old_index = i;
1320674ae819SStefano Zampini       }
1321674ae819SStefano Zampini     }
1322674ae819SStefano Zampini   }
1323674ae819SStefano Zampini   ierr = VecRestoreArray(global_vec,&array);CHKERRQ(ierr);
1324674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
1325674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1326674ae819SStefano Zampini   ierr = VecScatterEnd  (scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1327674ae819SStefano Zampini   /* get global ordering of local dofs */
1328674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
1329674ae819SStefano Zampini   if (local_dofs_mult) {
1330674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
13315b08dc53SStefano Zampini       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-local_dofs_mult[i];
1332674ae819SStefano Zampini     }
1333674ae819SStefano Zampini   } else {
1334674ae819SStefano Zampini     for (i=0;i<n_local_dofs;i++) {
13355b08dc53SStefano Zampini       temp_global_dofs[i] = (PetscInt)PetscRealPart(array[local_dofs[i]])-1;
1336674ae819SStefano Zampini     }
1337674ae819SStefano Zampini   }
1338674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
1339674ae819SStefano Zampini   /* free workspace */
1340674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1341674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1342674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1343674ae819SStefano Zampini   ierr = PetscFree(dof_sizes);CHKERRQ(ierr);
1344674ae819SStefano Zampini   ierr = PetscFree(dof_displs);CHKERRQ(ierr);
1345674ae819SStefano Zampini   /* return pointer to global ordering of local dofs */
1346674ae819SStefano Zampini   *global_numbering_subset = temp_global_dofs;
1347674ae819SStefano Zampini   PetscFunctionReturn(0);
1348674ae819SStefano Zampini }
1349