xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 5334bea6d5df22191b4a24d5be0b5c9b77e92b8b)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3674ae819SStefano Zampini 
4674ae819SStefano Zampini #undef __FUNCT__
5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPMatContext"
6674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
7674ae819SStefano Zampini {
8674ae819SStefano Zampini   FETIDPMat_ctx  newctx;
9674ae819SStefano Zampini   PetscErrorCode ierr;
10674ae819SStefano Zampini 
11674ae819SStefano Zampini   PetscFunctionBegin;
12854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
13674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
14674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
15674ae819SStefano Zampini   newctx->pc              = pc;
16674ae819SStefano Zampini   *fetidpmat_ctx          = newctx;
17674ae819SStefano Zampini   PetscFunctionReturn(0);
18674ae819SStefano Zampini }
19674ae819SStefano Zampini 
20674ae819SStefano Zampini #undef __FUNCT__
21674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext"
22674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
23674ae819SStefano Zampini {
24674ae819SStefano Zampini   FETIDPPC_ctx   newctx;
25674ae819SStefano Zampini   PetscErrorCode ierr;
26674ae819SStefano Zampini 
27674ae819SStefano Zampini   PetscFunctionBegin;
28854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
29674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
30674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
31674ae819SStefano Zampini   newctx->pc              = pc;
32674ae819SStefano Zampini   *fetidppc_ctx           = newctx;
33674ae819SStefano Zampini   PetscFunctionReturn(0);
34674ae819SStefano Zampini }
35674ae819SStefano Zampini 
36674ae819SStefano Zampini #undef __FUNCT__
37674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat"
38674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
39674ae819SStefano Zampini {
40674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
41674ae819SStefano Zampini   PetscErrorCode ierr;
42674ae819SStefano Zampini 
43674ae819SStefano Zampini   PetscFunctionBegin;
44674ae819SStefano Zampini   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
45674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
46674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
47674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
48674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
49674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
50457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr);
51457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr);
52457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr);
53457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr);
54457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr);
556cc1294bSstefano_zampini   ierr = VecDestroy(&mat_ctx->rhs_flip);CHKERRQ(ierr);
56457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr);
57457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr);
58457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr);
59674ae819SStefano Zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
60457d4a33Sstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr);
619cc7774eSstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr);
62674ae819SStefano Zampini   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
63674ae819SStefano Zampini   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
64674ae819SStefano Zampini   PetscFunctionReturn(0);
65674ae819SStefano Zampini }
66674ae819SStefano Zampini 
67674ae819SStefano Zampini #undef __FUNCT__
68674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC"
69674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
70674ae819SStefano Zampini {
71674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
72674ae819SStefano Zampini   PetscErrorCode ierr;
73674ae819SStefano Zampini 
74674ae819SStefano Zampini   PetscFunctionBegin;
75674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
76674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
77674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
78674ae819SStefano Zampini   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
7968668abeSStefano Zampini   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
80674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
81457d4a33Sstefano_zampini   ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr);
82457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr);
83457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr);
84674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
85674ae819SStefano Zampini   PetscFunctionReturn(0);
86674ae819SStefano Zampini }
87674ae819SStefano Zampini 
88674ae819SStefano Zampini #undef __FUNCT__
89674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext"
90674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
91674ae819SStefano Zampini {
92674ae819SStefano Zampini   PetscErrorCode ierr;
93674ae819SStefano Zampini   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
94674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
95674ae819SStefano Zampini   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
96674ae819SStefano Zampini   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
97674ae819SStefano Zampini   MPI_Comm       comm;
98674ae819SStefano Zampini   Mat            ScalingMat;
99457d4a33Sstefano_zampini   Vec            fetidp_global;
100674ae819SStefano Zampini   IS             IS_l2g_lambda;
101a1c0d0daSstefano_zampini   IS             subset,subset_mult,subset_n,isvert;
102674ae819SStefano Zampini   PetscBool      skip_node,fully_redundant;
103674ae819SStefano Zampini   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
104dc456d91SStefano Zampini   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
10576ec1555SBarry Smith   PetscMPIInt    rank,size,buf_size,neigh;
106674ae819SStefano Zampini   PetscScalar    scalar_value;
107a1c0d0daSstefano_zampini   const PetscInt *vertex_indices;
108dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
109dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
110674ae819SStefano Zampini   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
111674ae819SStefano Zampini   PetscScalar    *array,*scaling_factors,*vals_B_delta;
112674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
113457d4a33Sstefano_zampini   PetscLayout    llay;
114a1c0d0daSstefano_zampini 
115674ae819SStefano Zampini   /* For communication of scaling factors */
116674ae819SStefano Zampini   PetscInt       *ptrs_buffer,neigh_position;
117674ae819SStefano Zampini   PetscScalar    **all_factors,*send_buffer,*recv_buffer;
118674ae819SStefano Zampini   MPI_Request    *send_reqs,*recv_reqs;
119a1c0d0daSstefano_zampini 
120457d4a33Sstefano_zampini   /* saddlepoint */
121457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
122457d4a33Sstefano_zampini   PetscLayout            play;
123457d4a33Sstefano_zampini   IS                     gP,pP;
124457d4a33Sstefano_zampini   PetscInt               nPl,nPg,nPgl;
125674ae819SStefano Zampini 
126674ae819SStefano Zampini   PetscFunctionBegin;
127674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
128674ae819SStefano Zampini   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12976ec1555SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
130674ae819SStefano Zampini 
131457d4a33Sstefano_zampini   /* saddlepoint */
132457d4a33Sstefano_zampini   nPl      = 0;
133457d4a33Sstefano_zampini   nPg      = 0;
134457d4a33Sstefano_zampini   nPgl     = 0;
135457d4a33Sstefano_zampini   gP       = NULL;
136457d4a33Sstefano_zampini   pP       = NULL;
137457d4a33Sstefano_zampini   l2gmap_p = NULL;
138457d4a33Sstefano_zampini   play     = NULL;
139457d4a33Sstefano_zampini   ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr);
140022d8d2bSstefano_zampini   if (pP) { /* saddle point */
141457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
142457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr);
143457d4a33Sstefano_zampini     if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
144457d4a33Sstefano_zampini     ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr);
145457d4a33Sstefano_zampini     ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr);
146457d4a33Sstefano_zampini     ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr);
147457d4a33Sstefano_zampini     ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr);
148457d4a33Sstefano_zampini     ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr);
149457d4a33Sstefano_zampini 
150457d4a33Sstefano_zampini     /* interface pressure matrix */
151457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr);
152457d4a33Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */
153457d4a33Sstefano_zampini       IS Pg;
154457d4a33Sstefano_zampini 
155457d4a33Sstefano_zampini       ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr);
156457d4a33Sstefano_zampini       ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr);
157457d4a33Sstefano_zampini       ierr = ISDestroy(&Pg);CHKERRQ(ierr);
158457d4a33Sstefano_zampini       ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr);
159457d4a33Sstefano_zampini       ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr);
160457d4a33Sstefano_zampini       ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr);
161457d4a33Sstefano_zampini       ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr);
162457d4a33Sstefano_zampini       ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr);
163457d4a33Sstefano_zampini       ierr = PetscLayoutSetUp(play);CHKERRQ(ierr);
164457d4a33Sstefano_zampini     } else {
165457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr);
166457d4a33Sstefano_zampini       ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr);
167457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr);
168457d4a33Sstefano_zampini       ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr);
169457d4a33Sstefano_zampini       ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr);
170457d4a33Sstefano_zampini       ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr);
171457d4a33Sstefano_zampini       ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr);
172457d4a33Sstefano_zampini     }
173457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr);
174457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr);
175457d4a33Sstefano_zampini 
176457d4a33Sstefano_zampini     /* import matrices for interface pressures coupling */
177457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr);
178457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
179457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr);
180457d4a33Sstefano_zampini 
181457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr);
182457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
183457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr);
184457d4a33Sstefano_zampini 
185457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
186457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
187457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
188457d4a33Sstefano_zampini 
189457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
190457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
191457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
1926cc1294bSstefano_zampini 
1936cc1294bSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);CHKERRQ(ierr);
1946cc1294bSstefano_zampini     if (fetidpmat_ctx->rhs_flip) {
1956cc1294bSstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);CHKERRQ(ierr);
1966cc1294bSstefano_zampini     }
197457d4a33Sstefano_zampini   }
198457d4a33Sstefano_zampini 
199674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
200329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
201674ae819SStefano Zampini 
202674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
203674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
204674ae819SStefano Zampini   n_local_lambda = 0;
205674ae819SStefano Zampini   partial_sum = 0;
206674ae819SStefano Zampini   n_boundary_dofs = 0;
207674ae819SStefano Zampini   s = 0;
208a1c0d0daSstefano_zampini 
209674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
210a1c0d0daSstefano_zampini   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr);
211a1c0d0daSstefano_zampini   ierr = ISGetLocalSize(isvert,&n_vertices);CHKERRQ(ierr);
212a1c0d0daSstefano_zampini   ierr = ISGetIndices(isvert,&vertex_indices);CHKERRQ(ierr);
213a1c0d0daSstefano_zampini 
214674ae819SStefano Zampini   dual_size = pcis->n_B-n_vertices;
215785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
216785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
217785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
218674ae819SStefano Zampini 
219674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
220674ae819SStefano Zampini   for (i=0;i<pcis->n;i++){
221674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
2222d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
223674ae819SStefano Zampini     skip_node = PETSC_FALSE;
224674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
225674ae819SStefano Zampini       skip_node = PETSC_TRUE;
226674ae819SStefano Zampini       s++;
227674ae819SStefano Zampini     }
2282d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
2292d456bbaSstefano_zampini     if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
230674ae819SStefano Zampini     if (!skip_node) {
231674ae819SStefano Zampini       if (fully_redundant) {
232674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
233674ae819SStefano Zampini         n_lambda_for_dof = (j*(j+1))/2;
234674ae819SStefano Zampini       } else {
235674ae819SStefano Zampini         n_lambda_for_dof = j;
236674ae819SStefano Zampini       }
237674ae819SStefano Zampini       n_local_lambda += j;
238674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
239674ae819SStefano Zampini       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
240674ae819SStefano Zampini       /* store some data needed */
241674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
242674ae819SStefano Zampini       aux_local_numbering_1[partial_sum] = i;
243674ae819SStefano Zampini       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
244674ae819SStefano Zampini       partial_sum++;
245674ae819SStefano Zampini     }
246674ae819SStefano Zampini   }
247674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
248a1c0d0daSstefano_zampini   ierr = ISRestoreIndices(isvert,&vertex_indices);CHKERRQ(ierr);
249a1c0d0daSstefano_zampini   ierr = PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr);
2502d456bbaSstefano_zampini   dual_size = partial_sum;
251674ae819SStefano Zampini 
252674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
253dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
2543bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
255dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
256dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
2573d996552SStefano Zampini   ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr);
258dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
2593d996552SStefano Zampini 
2604fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG)
2614fe826edSStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
2624fe826edSStefano Zampini   ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2634fe826edSStefano Zampini   ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2644fe826edSStefano Zampini   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
2654fe826edSStefano Zampini   i = (PetscInt)PetscRealPart(scalar_value);
2666c4ed002SBarry Smith   if (i != fetidpmat_ctx->n_lambda) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global number of multipliers mismatch! (%d!=%d)\n",fetidpmat_ctx->n_lambda,i);
2674fe826edSStefano Zampini #endif
268674ae819SStefano Zampini 
269674ae819SStefano Zampini   /* init data for scaling factors exchange */
270674ae819SStefano Zampini   partial_sum = 0;
271785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
2724b2aedd3SStefano Zampini   ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr);
2734b2aedd3SStefano Zampini   ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr);
27419c16490Sstefano_zampini   ierr = PetscMalloc1(pcis->n+1,&all_factors);CHKERRQ(ierr);
2754b2aedd3SStefano Zampini   if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
276674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
277674ae819SStefano Zampini     partial_sum += pcis->n_shared[i];
278674ae819SStefano Zampini     ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
279674ae819SStefano Zampini   }
280785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
281785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
282785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
283674ae819SStefano Zampini   for (i=0;i<pcis->n-1;i++) {
284674ae819SStefano Zampini     j = mat_graph->count[i];
285674ae819SStefano Zampini     all_factors[i+1]=all_factors[i]+j;
286674ae819SStefano Zampini   }
287674ae819SStefano Zampini   /* scatter B scaling to N vec */
288674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
289674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
290674ae819SStefano Zampini   /* communications */
2912b095fd8SStefano Zampini   ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
292674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
293674ae819SStefano Zampini     for (j=0;j<pcis->n_shared[i];j++) {
294674ae819SStefano Zampini       send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
295674ae819SStefano Zampini     }
296674ae819SStefano Zampini     ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
297674ae819SStefano Zampini     ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
298674ae819SStefano Zampini     ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
299674ae819SStefano Zampini     ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
300674ae819SStefano Zampini   }
3012b095fd8SStefano Zampini   ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
3024b2aedd3SStefano Zampini   if (pcis->n_neigh > 0) {
3034b2aedd3SStefano Zampini     ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3044b2aedd3SStefano Zampini   }
305674ae819SStefano Zampini   /* put values in correct places */
306674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
307674ae819SStefano Zampini     for (j=0;j<pcis->n_shared[i];j++) {
308674ae819SStefano Zampini       k = pcis->shared[i][j];
309674ae819SStefano Zampini       neigh_position = 0;
310674ae819SStefano Zampini       while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
311674ae819SStefano Zampini       all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
312674ae819SStefano Zampini     }
313674ae819SStefano Zampini   }
3144b2aedd3SStefano Zampini   if (pcis->n_neigh > 0) {
3154b2aedd3SStefano Zampini     ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3164b2aedd3SStefano Zampini   }
317674ae819SStefano Zampini   ierr = PetscFree(send_reqs);CHKERRQ(ierr);
318674ae819SStefano Zampini   ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
319674ae819SStefano Zampini   ierr = PetscFree(send_buffer);CHKERRQ(ierr);
320674ae819SStefano Zampini   ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
321674ae819SStefano Zampini   ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
322674ae819SStefano Zampini 
323674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
324785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
325785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
326785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
327785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
328785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
329dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
330674ae819SStefano Zampini   partial_sum=0;
331dc456d91SStefano Zampini   cum = 0;
332674ae819SStefano Zampini   for (i=0;i<dual_size;i++) {
333dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
334674ae819SStefano Zampini     j = mat_graph->count[aux_local_numbering_1[i]];
335674ae819SStefano Zampini     aux_sums[0]=0;
336674ae819SStefano Zampini     for (s=1;s<j;s++) {
337674ae819SStefano Zampini       aux_sums[s]=aux_sums[s-1]+j-s+1;
338674ae819SStefano Zampini     }
339674ae819SStefano Zampini     array = all_factors[aux_local_numbering_1[i]];
340674ae819SStefano Zampini     n_neg_values = 0;
3412a7da448SStefano Zampini     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
342674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
343674ae819SStefano Zampini     if (fully_redundant) {
344674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
345674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
346674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
347674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=-1.0;
348674ae819SStefano Zampini         scaling_factors[partial_sum+s]=array[s];
349674ae819SStefano Zampini       }
350674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
351674ae819SStefano Zampini         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
352674ae819SStefano Zampini         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
353674ae819SStefano Zampini         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
354674ae819SStefano Zampini         scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
355674ae819SStefano Zampini       }
356674ae819SStefano Zampini       partial_sum += j;
357674ae819SStefano Zampini     } else {
358674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
359674ae819SStefano Zampini       for (s=0;s<j;s++) {
360674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=n_global_lambda+s;
361674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
362674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=0.0;
363674ae819SStefano Zampini       }
364674ae819SStefano Zampini       /* B_delta */
365674ae819SStefano Zampini       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
366674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
367674ae819SStefano Zampini       }
368674ae819SStefano Zampini       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
369674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values]=1.0;
370674ae819SStefano Zampini       }
371674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999*/
372674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
373674ae819SStefano Zampini         scalar_value = 0.0;
374674ae819SStefano Zampini         for (k=0;k<s+1;k++) {
375674ae819SStefano Zampini           scalar_value += array[k];
376674ae819SStefano Zampini         }
377674ae819SStefano Zampini         scaling_factors[partial_sum+s] = -scalar_value;
378674ae819SStefano Zampini       }
379674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
380674ae819SStefano Zampini         scalar_value = 0.0;
381674ae819SStefano Zampini         for (k=s+n_neg_values;k<j;k++) {
382674ae819SStefano Zampini           scalar_value += array[k];
383674ae819SStefano Zampini         }
384674ae819SStefano Zampini         scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
385674ae819SStefano Zampini       }
386674ae819SStefano Zampini       partial_sum += j;
387674ae819SStefano Zampini     }
388dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
389674ae819SStefano Zampini   }
390dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
391dc456d91SStefano Zampini   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
392dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
393674ae819SStefano Zampini   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
394674ae819SStefano Zampini   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
395674ae819SStefano Zampini   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
396674ae819SStefano Zampini   ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
397674ae819SStefano Zampini   ierr = PetscFree(all_factors);CHKERRQ(ierr);
398674ae819SStefano Zampini 
399674ae819SStefano Zampini   /* Create local part of B_delta */
400302440fdSBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
401674ae819SStefano Zampini   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
402674ae819SStefano Zampini   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
403674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
404674ae819SStefano Zampini   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
405674ae819SStefano Zampini   for (i=0;i<n_local_lambda;i++) {
406674ae819SStefano Zampini     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
407674ae819SStefano Zampini   }
408674ae819SStefano Zampini   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
409674ae819SStefano Zampini   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
410674ae819SStefano Zampini   ierr = MatAssemblyEnd  (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411674ae819SStefano Zampini 
412674ae819SStefano Zampini   if (fully_redundant) {
413302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
414674ae819SStefano Zampini     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
415674ae819SStefano Zampini     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
416674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
417674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
418674ae819SStefano Zampini       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
419674ae819SStefano Zampini     }
420674ae819SStefano Zampini     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
421674ae819SStefano Zampini     ierr = MatAssemblyEnd  (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
422674ae819SStefano Zampini     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
423674ae819SStefano Zampini     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
424674ae819SStefano Zampini   } else {
425302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
426674ae819SStefano Zampini     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
427674ae819SStefano Zampini     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
428674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
429674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
430674ae819SStefano Zampini       ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
431674ae819SStefano Zampini     }
432674ae819SStefano Zampini     ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
433674ae819SStefano Zampini     ierr = MatAssemblyEnd  (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
434674ae819SStefano Zampini   }
435674ae819SStefano Zampini   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
436674ae819SStefano Zampini   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
437674ae819SStefano Zampini 
438457d4a33Sstefano_zampini   /* Layout of multipliers */
439457d4a33Sstefano_zampini   ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr);
440457d4a33Sstefano_zampini   ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr);
441457d4a33Sstefano_zampini   ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
442457d4a33Sstefano_zampini   ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr);
443457d4a33Sstefano_zampini   ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr);
444457d4a33Sstefano_zampini 
445457d4a33Sstefano_zampini   /* fetidpmat sizes */
446457d4a33Sstefano_zampini   fetidpmat_ctx->n += nPgl;
447457d4a33Sstefano_zampini   fetidpmat_ctx->N  = fetidpmat_ctx->n_lambda+nPg;
448457d4a33Sstefano_zampini 
449457d4a33Sstefano_zampini   /* Local work vector of multipliers */
450457d4a33Sstefano_zampini   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
451457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
452457d4a33Sstefano_zampini   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
453457d4a33Sstefano_zampini 
454457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
455457d4a33Sstefano_zampini   ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr);
456457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr);
457457d4a33Sstefano_zampini   ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr);
458457d4a33Sstefano_zampini   ierr = VecSetUp(fetidp_global);CHKERRQ(ierr);
459457d4a33Sstefano_zampini 
4609eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
4619eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
462457d4a33Sstefano_zampini      of the FETI-DP linear system */
463457d4a33Sstefano_zampini   if (nPg) {
464af140850Sstefano_zampini     IS             IS_l2g_p,ais;
465457d4a33Sstefano_zampini     PetscLayout    alay;
466457d4a33Sstefano_zampini     const PetscInt *idxs,*pranges,*aranges,*lranges;
467af140850Sstefano_zampini     PetscInt       *l2g_indices_p,rst;
468457d4a33Sstefano_zampini 
469457d4a33Sstefano_zampini     ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr);
470457d4a33Sstefano_zampini     ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr);
471457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr);
472457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr);
473457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr);
474457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
475af140850Sstefano_zampini     /* shift local to global indices for pressure */
476457d4a33Sstefano_zampini     for (i=0;i<nPl;i++) {
477457d4a33Sstefano_zampini       PetscInt owner;
478457d4a33Sstefano_zampini 
479457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr);
480457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
481457d4a33Sstefano_zampini     }
482457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
483457d4a33Sstefano_zampini     ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr);
484af140850Sstefano_zampini 
485457d4a33Sstefano_zampini     /* local to global scatter for interface pressure */
486457d4a33Sstefano_zampini     ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr);
487457d4a33Sstefano_zampini     ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr);
488457d4a33Sstefano_zampini 
489af140850Sstefano_zampini     /* shift local to global indices for multipliers */
490457d4a33Sstefano_zampini     for (i=0;i<n_local_lambda;i++) {
491457d4a33Sstefano_zampini       PetscInt owner,ps;
492457d4a33Sstefano_zampini 
493457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr);
494457d4a33Sstefano_zampini       ps = pranges[owner+1]-pranges[owner];
495457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
496457d4a33Sstefano_zampini     }
497457d4a33Sstefano_zampini 
4989cc7774eSstefano_zampini     /* scatter from alldofs to interface pressures global fetidp vector */
4999cc7774eSstefano_zampini     ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr);
5009cc7774eSstefano_zampini     ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr);
501af140850Sstefano_zampini     ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr);
5029cc7774eSstefano_zampini     ierr = ISDestroy(&ais);CHKERRQ(ierr);
503457d4a33Sstefano_zampini   }
504457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr);
505457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr);
506457d4a33Sstefano_zampini   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
507a1c0d0daSstefano_zampini 
5089cc7774eSstefano_zampini   /* scatter from local to global multipliers */
509457d4a33Sstefano_zampini   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
510457d4a33Sstefano_zampini   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
511457d4a33Sstefano_zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr);
512a1c0d0daSstefano_zampini   ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr);
513457d4a33Sstefano_zampini 
514a1c0d0daSstefano_zampini   /* Create some work vectors needed by fetidp */
515674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
516674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
517674ae819SStefano Zampini   PetscFunctionReturn(0);
518674ae819SStefano Zampini }
519674ae819SStefano Zampini 
520674ae819SStefano Zampini #undef __FUNCT__
521674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext"
522674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
523674ae819SStefano Zampini {
524674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
525ed6c3d69SStefano Zampini   PC_IS          *pcis;
526f28b6018SStefano Zampini   PetscBool      lumped = PETSC_FALSE;
527674ae819SStefano Zampini   PetscErrorCode ierr;
528674ae819SStefano Zampini 
529674ae819SStefano Zampini   PetscFunctionBegin;
530674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
531674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
532674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
533674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
534674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
535674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
536674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
537674ae819SStefano Zampini   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
538f28b6018SStefano Zampini   /* create preconditioner */
539ed6c3d69SStefano Zampini   pcis = (PC_IS*)fetidppc_ctx->pc->data;
540f28b6018SStefano Zampini   /* Dirichlet preconditioner */
5419c2d02cdSstefano_zampini   ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL);CHKERRQ(ierr);
542f28b6018SStefano Zampini   if (!lumped) {
5439c2d02cdSstefano_zampini     IS        iP;
5449c2d02cdSstefano_zampini     PetscBool discrete_harmonic = PETSC_FALSE;
5459c2d02cdSstefano_zampini 
5469c2d02cdSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iP",(PetscObject*)&iP);CHKERRQ(ierr);
5479c2d02cdSstefano_zampini     if (iP) {
5489c2d02cdSstefano_zampini       ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL);CHKERRQ(ierr);
5499c2d02cdSstefano_zampini     }
5509c2d02cdSstefano_zampini     if (discrete_harmonic) {
5519c2d02cdSstefano_zampini       KSP        sksp;
5529c2d02cdSstefano_zampini       PC         pc;
5539c2d02cdSstefano_zampini       Mat        A_II,A_IB,A_BI;
5549c2d02cdSstefano_zampini       IS         aB;
5559c2d02cdSstefano_zampini       PetscInt   nb;
556*5334bea6Sstefano_zampini       PetscBool  isshell;
5579c2d02cdSstefano_zampini       KSPType    ksptype;
5589c2d02cdSstefano_zampini       const char *prefix;
5599c2d02cdSstefano_zampini 
5609c2d02cdSstefano_zampini       /*
5619c2d02cdSstefano_zampini         We constructs a Schur complement for
5629c2d02cdSstefano_zampini 
5639c2d02cdSstefano_zampini         | A_II A_ID |
5649c2d02cdSstefano_zampini         | A_DI A_DD |
5659c2d02cdSstefano_zampini 
5669c2d02cdSstefano_zampini         instead of
5679c2d02cdSstefano_zampini 
5689c2d02cdSstefano_zampini         | A_II  B^t_II A_ID |
5699c2d02cdSstefano_zampini         | B_II -C_II   B_ID |
5709c2d02cdSstefano_zampini         | A_DI  B^t_ID A_DD |
5719c2d02cdSstefano_zampini 
5729c2d02cdSstefano_zampini       */
5739c2d02cdSstefano_zampini       ierr = ISGetLocalSize(pcis->is_B_local,&nb);CHKERRQ(ierr);
5749c2d02cdSstefano_zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB);CHKERRQ(ierr);
5759c2d02cdSstefano_zampini       ierr = MatGetSubMatrix(pcis->A_II,iP,iP,MAT_INITIAL_MATRIX,&A_II);CHKERRQ(ierr);
5769c2d02cdSstefano_zampini       ierr = MatGetSubMatrix(pcis->A_IB,iP,aB,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
5779c2d02cdSstefano_zampini       ierr = MatGetSubMatrix(pcis->A_BI,aB,iP,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
5789c2d02cdSstefano_zampini       ierr = MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
5799c2d02cdSstefano_zampini 
5809c2d02cdSstefano_zampini       /* propagate settings of solver */
5819c2d02cdSstefano_zampini       ierr = MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp);CHKERRQ(ierr);
5829c2d02cdSstefano_zampini       ierr = KSPGetType(pcis->ksp_D,&ksptype);CHKERRQ(ierr);
5839c2d02cdSstefano_zampini       ierr = KSPSetType(sksp,ksptype);CHKERRQ(ierr);
5849c2d02cdSstefano_zampini       ierr = KSPGetPC(pcis->ksp_D,&pc);CHKERRQ(ierr);
585*5334bea6Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell);CHKERRQ(ierr);
586*5334bea6Sstefano_zampini       if (!isshell) {
587*5334bea6Sstefano_zampini         MatSolverPackage solver;
588*5334bea6Sstefano_zampini         PCType           pctype;
589*5334bea6Sstefano_zampini 
5909c2d02cdSstefano_zampini         ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
5919c2d02cdSstefano_zampini         ierr = PCFactorGetMatSolverPackage(pc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
5929c2d02cdSstefano_zampini         ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr);
5939c2d02cdSstefano_zampini         ierr = PCSetType(pc,pctype);CHKERRQ(ierr);
5949c2d02cdSstefano_zampini         if (solver) {
5959c2d02cdSstefano_zampini           ierr = PCFactorSetMatSolverPackage(pc,solver);CHKERRQ(ierr);
5969c2d02cdSstefano_zampini         }
597*5334bea6Sstefano_zampini       } else {
598*5334bea6Sstefano_zampini         ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr);
599*5334bea6Sstefano_zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
600*5334bea6Sstefano_zampini       }
6019c2d02cdSstefano_zampini       ierr = MatDestroy(&A_II);CHKERRQ(ierr);
6029c2d02cdSstefano_zampini       ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
6039c2d02cdSstefano_zampini       ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
6049c2d02cdSstefano_zampini       ierr = ISDestroy(&aB);CHKERRQ(ierr);
6059c2d02cdSstefano_zampini       ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr);
6069c2d02cdSstefano_zampini       ierr = KSPSetOptionsPrefix(sksp,prefix);CHKERRQ(ierr);
6079c2d02cdSstefano_zampini       ierr = KSPAppendOptionsPrefix(sksp,"harmonic_");CHKERRQ(ierr);
6083016320fSstefano_zampini       ierr = KSPSetFromOptions(sksp);CHKERRQ(ierr);
6099c2d02cdSstefano_zampini     } else { /* default Dirichlet preconditioner is pde-harmonic */
610ed6c3d69SStefano Zampini       ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
611ed6c3d69SStefano Zampini       ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
6129c2d02cdSstefano_zampini     }
613f28b6018SStefano Zampini   } else {
614f28b6018SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr);
615f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
616f28b6018SStefano Zampini   }
617af140850Sstefano_zampini   /* saddle-point */
618af140850Sstefano_zampini   if (mat_ctx->xPg) {
619af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr);
620af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
621af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr);
622af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
6236cc1294bSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PKSP",(PetscObject*)&fetidppc_ctx->kP);CHKERRQ(ierr);
6246cc1294bSstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidppc_ctx->kP);CHKERRQ(ierr);
625af140850Sstefano_zampini   }
626674ae819SStefano Zampini   PetscFunctionReturn(0);
627674ae819SStefano Zampini }
628674ae819SStefano Zampini 
629674ae819SStefano Zampini #undef __FUNCT__
630674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult"
631674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
632674ae819SStefano Zampini {
633674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
634617d11aeSStefano Zampini   PC_BDDC        *pcbddc;
635674ae819SStefano Zampini   PC_IS          *pcis;
636674ae819SStefano Zampini   PetscErrorCode ierr;
637674ae819SStefano Zampini 
638674ae819SStefano Zampini   PetscFunctionBegin;
639674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
640674ae819SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
641617d11aeSStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
642674ae819SStefano Zampini   /* Application of B_delta^T */
643af140850Sstefano_zampini   ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr);
644674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
645674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
646674ae819SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
647af140850Sstefano_zampini 
648af140850Sstefano_zampini   /* Add contribution from saddle point */
649af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
650af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
651af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
652af140850Sstefano_zampini     if (pcbddc->switch_static) {
653af140850Sstefano_zampini       ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr);
654af140850Sstefano_zampini     }
655af140850Sstefano_zampini     ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
656af140850Sstefano_zampini   } else {
657af140850Sstefano_zampini     if (pcbddc->switch_static) {
658674ae819SStefano Zampini       ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
659af140850Sstefano_zampini     }
660af140850Sstefano_zampini   }
661af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
662617d11aeSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
663dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
664c7ffc8ceSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
665af140850Sstefano_zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
666674ae819SStefano Zampini   /* Application of B_delta */
667674ae819SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
668af140850Sstefano_zampini   /* Contribution from boundary pressures */
669af140850Sstefano_zampini   if (mat_ctx->C) {
670af140850Sstefano_zampini     const PetscScalar *lx;
671af140850Sstefano_zampini     PetscScalar       *ly;
672af140850Sstefano_zampini 
673af140850Sstefano_zampini     /* pressures ordered first in x and y */
674af140850Sstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
675af140850Sstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
676af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr);
677af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr);
678af140850Sstefano_zampini     ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr);
679af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr);
680af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr);
681af140850Sstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
682af140850Sstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
683af140850Sstefano_zampini   }
684af140850Sstefano_zampini   /* Add contribution from saddle point */
685af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
686af140850Sstefano_zampini     ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
687af140850Sstefano_zampini     if (pcbddc->switch_static) {
688af140850Sstefano_zampini       ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
689af140850Sstefano_zampini     }
690af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
691af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
692af140850Sstefano_zampini   }
693674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
694674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
695674ae819SStefano Zampini   PetscFunctionReturn(0);
696674ae819SStefano Zampini }
697674ae819SStefano Zampini 
698674ae819SStefano Zampini #undef __FUNCT__
699edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose"
700edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
701edf7251bSStefano Zampini {
702edf7251bSStefano Zampini   FETIDPMat_ctx  mat_ctx;
703edf7251bSStefano Zampini   PC_IS          *pcis;
704edf7251bSStefano Zampini   PetscErrorCode ierr;
705edf7251bSStefano Zampini 
706edf7251bSStefano Zampini   PetscFunctionBegin;
707edf7251bSStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
708edf7251bSStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
709edf7251bSStefano Zampini   /* Application of B_delta^T */
710edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
711edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
712edf7251bSStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
713edf7251bSStefano Zampini   /* Application of \widetilde{S}^-1 */
714edf7251bSStefano Zampini   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
715edf7251bSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr);
716edf7251bSStefano Zampini   /* Application of B_delta */
717edf7251bSStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
718edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
719edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
720edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
721edf7251bSStefano Zampini   PetscFunctionReturn(0);
722edf7251bSStefano Zampini }
723edf7251bSStefano Zampini 
724edf7251bSStefano Zampini #undef __FUNCT__
725674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply"
726674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
727674ae819SStefano Zampini {
728674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
729674ae819SStefano Zampini   PC_IS          *pcis;
730674ae819SStefano Zampini   PetscErrorCode ierr;
731674ae819SStefano Zampini 
732674ae819SStefano Zampini   PetscFunctionBegin;
733302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
734674ae819SStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
735674ae819SStefano Zampini   /* Application of B_Ddelta^T */
736674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
737674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
738674ae819SStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
739674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
740ed6c3d69SStefano Zampini   /* Application of local Schur complement */
741ed6c3d69SStefano Zampini   ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
742edf7251bSStefano Zampini   /* Application of B_Ddelta */
743edf7251bSStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
744edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
745edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
746edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
747b519380fSstefano_zampini   /* interface pressure preconditioner */
748b519380fSstefano_zampini   if (pc_ctx->kP) {
749b519380fSstefano_zampini     const PetscScalar *lx;
750b519380fSstefano_zampini     PetscScalar *ly;
751b519380fSstefano_zampini 
752b519380fSstefano_zampini     /* pressures ordered first in x and y */
753b519380fSstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
754b519380fSstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
755b519380fSstefano_zampini     ierr = VecPlaceArray(pc_ctx->xPg,lx);CHKERRQ(ierr);
756b519380fSstefano_zampini     ierr = VecPlaceArray(pc_ctx->yPg,ly);CHKERRQ(ierr);
757b519380fSstefano_zampini     ierr = KSPSolve(pc_ctx->kP,pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr);
758b519380fSstefano_zampini     ierr = VecResetArray(pc_ctx->xPg);CHKERRQ(ierr);
759b519380fSstefano_zampini     ierr = VecResetArray(pc_ctx->yPg);CHKERRQ(ierr);
760b519380fSstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
761b519380fSstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
762b519380fSstefano_zampini   }
763edf7251bSStefano Zampini   PetscFunctionReturn(0);
764edf7251bSStefano Zampini }
765edf7251bSStefano Zampini 
766edf7251bSStefano Zampini #undef __FUNCT__
767edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose"
768edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
769edf7251bSStefano Zampini {
770edf7251bSStefano Zampini   FETIDPPC_ctx   pc_ctx;
771edf7251bSStefano Zampini   PC_IS          *pcis;
772edf7251bSStefano Zampini   PetscErrorCode ierr;
773edf7251bSStefano Zampini 
774edf7251bSStefano Zampini   PetscFunctionBegin;
775302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
776edf7251bSStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
777edf7251bSStefano Zampini   /* Application of B_Ddelta^T */
778edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
779edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
780edf7251bSStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
781edf7251bSStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
782ed6c3d69SStefano Zampini   /* Application of local Schur complement */
783ed6c3d69SStefano Zampini   ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
784674ae819SStefano Zampini   /* Application of B_Ddelta */
785674ae819SStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
786674ae819SStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
787674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
788674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
789674ae819SStefano Zampini   PetscFunctionReturn(0);
790674ae819SStefano Zampini }
791c45b8d2dSstefano_zampini 
792c45b8d2dSstefano_zampini #undef __FUNCT__
793c45b8d2dSstefano_zampini #define __FUNCT__ "FETIDPPCView"
794c45b8d2dSstefano_zampini PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
795c45b8d2dSstefano_zampini {
796c45b8d2dSstefano_zampini   FETIDPPC_ctx      pc_ctx;
797c45b8d2dSstefano_zampini   PetscBool         iascii;
798c45b8d2dSstefano_zampini   PetscViewer       sviewer;
799c45b8d2dSstefano_zampini   PetscErrorCode    ierr;
800c45b8d2dSstefano_zampini 
801c45b8d2dSstefano_zampini   PetscFunctionBegin;
802c45b8d2dSstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
803c45b8d2dSstefano_zampini   if (iascii) {
804c45b8d2dSstefano_zampini     PetscMPIInt rank;
805c45b8d2dSstefano_zampini     PetscBool   isschur;
806c45b8d2dSstefano_zampini 
807c45b8d2dSstefano_zampini     ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
808c45b8d2dSstefano_zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
809c45b8d2dSstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur);CHKERRQ(ierr);
810c45b8d2dSstefano_zampini     if (isschur) {
811c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP multipliers Dirichlet preconditioner (just from rank 0)\n");CHKERRQ(ierr);
812c45b8d2dSstefano_zampini     } else {
813c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP multipliers Lumped preconditioner (just from rank 0)\n");CHKERRQ(ierr);
814c45b8d2dSstefano_zampini     }
815c45b8d2dSstefano_zampini     ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
816c45b8d2dSstefano_zampini     if (!rank) {
817c45b8d2dSstefano_zampini       ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
818c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIAddTab(sviewer,2);CHKERRQ(ierr);
819c45b8d2dSstefano_zampini       ierr = MatView(pc_ctx->S_j,sviewer);CHKERRQ(ierr);
820c45b8d2dSstefano_zampini       ierr = PetscViewerASCIISubtractTab(sviewer,2);CHKERRQ(ierr);
821c45b8d2dSstefano_zampini       ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr);
822c45b8d2dSstefano_zampini     }
823c45b8d2dSstefano_zampini     ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
824c45b8d2dSstefano_zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
825c45b8d2dSstefano_zampini     if (pc_ctx->kP) {
826c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP pressure preconditioner\n");CHKERRQ(ierr);
827c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIAddTab(viewer,2);CHKERRQ(ierr);
828c45b8d2dSstefano_zampini       ierr = KSPView(pc_ctx->kP,viewer);CHKERRQ(ierr);
829c45b8d2dSstefano_zampini       ierr = PetscViewerASCIISubtractTab(viewer,2);CHKERRQ(ierr);
830c45b8d2dSstefano_zampini     }
831c45b8d2dSstefano_zampini   }
832c45b8d2dSstefano_zampini   PetscFunctionReturn(0);
833c45b8d2dSstefano_zampini }
834