xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision be8bcbd05a63c50268e4bcfeac28378736007e14)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
325d06dbeSstefano_zampini #include <petscblaslapack.h>
425d06dbeSstefano_zampini 
525d06dbeSstefano_zampini #undef __FUNCT__
625d06dbeSstefano_zampini #define __FUNCT__ "MatMult_BDdelta_deluxe_nonred"
725d06dbeSstefano_zampini static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
825d06dbeSstefano_zampini {
925d06dbeSstefano_zampini   BDdelta_DN     ctx;
1025d06dbeSstefano_zampini   PetscErrorCode ierr;
1125d06dbeSstefano_zampini 
1225d06dbeSstefano_zampini   PetscFunctionBegin;
1325d06dbeSstefano_zampini   ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
1425d06dbeSstefano_zampini   ierr = MatMultTranspose(ctx->BD,x,ctx->work);CHKERRQ(ierr);
1525d06dbeSstefano_zampini   ierr = KSPSolveTranspose(ctx->kBD,ctx->work,y);CHKERRQ(ierr);
1625d06dbeSstefano_zampini   PetscFunctionReturn(0);
1725d06dbeSstefano_zampini }
1825d06dbeSstefano_zampini 
1925d06dbeSstefano_zampini #undef __FUNCT__
2025d06dbeSstefano_zampini #define __FUNCT__ "MatMultTranspose_BDdelta_deluxe_nonred"
2125d06dbeSstefano_zampini static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
2225d06dbeSstefano_zampini {
2325d06dbeSstefano_zampini   BDdelta_DN     ctx;
2425d06dbeSstefano_zampini   PetscErrorCode ierr;
2525d06dbeSstefano_zampini 
2625d06dbeSstefano_zampini   PetscFunctionBegin;
2725d06dbeSstefano_zampini   ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
2825d06dbeSstefano_zampini   ierr = KSPSolve(ctx->kBD,x,ctx->work);CHKERRQ(ierr);
2925d06dbeSstefano_zampini   ierr = MatMult(ctx->BD,ctx->work,y);CHKERRQ(ierr);
3025d06dbeSstefano_zampini   PetscFunctionReturn(0);
3125d06dbeSstefano_zampini }
3225d06dbeSstefano_zampini 
3325d06dbeSstefano_zampini #undef __FUNCT__
3425d06dbeSstefano_zampini #define __FUNCT__ "MatDestroy_BDdelta_deluxe_nonred"
3525d06dbeSstefano_zampini static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
3625d06dbeSstefano_zampini {
3725d06dbeSstefano_zampini   BDdelta_DN     ctx;
3825d06dbeSstefano_zampini   PetscErrorCode ierr;
3925d06dbeSstefano_zampini 
4025d06dbeSstefano_zampini   PetscFunctionBegin;
4125d06dbeSstefano_zampini   ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
4225d06dbeSstefano_zampini   ierr = MatDestroy(&ctx->BD);CHKERRQ(ierr);
4325d06dbeSstefano_zampini   ierr = KSPDestroy(&ctx->kBD);CHKERRQ(ierr);
4425d06dbeSstefano_zampini   ierr = VecDestroy(&ctx->work);CHKERRQ(ierr);
4525d06dbeSstefano_zampini   ierr = PetscFree(ctx);CHKERRQ(ierr);
4625d06dbeSstefano_zampini   PetscFunctionReturn(0);
4725d06dbeSstefano_zampini }
4825d06dbeSstefano_zampini 
49674ae819SStefano Zampini 
50674ae819SStefano Zampini #undef __FUNCT__
51674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPMatContext"
52674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
53674ae819SStefano Zampini {
54674ae819SStefano Zampini   FETIDPMat_ctx  newctx;
55674ae819SStefano Zampini   PetscErrorCode ierr;
56674ae819SStefano Zampini 
57674ae819SStefano Zampini   PetscFunctionBegin;
58854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
59674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
60674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
61674ae819SStefano Zampini   newctx->pc              = pc;
62674ae819SStefano Zampini   *fetidpmat_ctx          = newctx;
63674ae819SStefano Zampini   PetscFunctionReturn(0);
64674ae819SStefano Zampini }
65674ae819SStefano Zampini 
66674ae819SStefano Zampini #undef __FUNCT__
67674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext"
68674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
69674ae819SStefano Zampini {
70674ae819SStefano Zampini   FETIDPPC_ctx   newctx;
71674ae819SStefano Zampini   PetscErrorCode ierr;
72674ae819SStefano Zampini 
73674ae819SStefano Zampini   PetscFunctionBegin;
74854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
75674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
76674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
77674ae819SStefano Zampini   newctx->pc              = pc;
78674ae819SStefano Zampini   *fetidppc_ctx           = newctx;
79674ae819SStefano Zampini   PetscFunctionReturn(0);
80674ae819SStefano Zampini }
81674ae819SStefano Zampini 
82674ae819SStefano Zampini #undef __FUNCT__
83674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat"
84674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
85674ae819SStefano Zampini {
86674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
87674ae819SStefano Zampini   PetscErrorCode ierr;
88674ae819SStefano Zampini 
89674ae819SStefano Zampini   PetscFunctionBegin;
90674ae819SStefano Zampini   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
91674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
92674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
93674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
94674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
95674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
96457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr);
97457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr);
98457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr);
99457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr);
100457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr);
1016cc1294bSstefano_zampini   ierr = VecDestroy(&mat_ctx->rhs_flip);CHKERRQ(ierr);
102457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr);
103457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr);
104457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr);
105674ae819SStefano Zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
106457d4a33Sstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr);
1079cc7774eSstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr);
108674ae819SStefano Zampini   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
109674ae819SStefano Zampini   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
110674ae819SStefano Zampini   PetscFunctionReturn(0);
111674ae819SStefano Zampini }
112674ae819SStefano Zampini 
113674ae819SStefano Zampini #undef __FUNCT__
114674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC"
115674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
116674ae819SStefano Zampini {
117674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
118674ae819SStefano Zampini   PetscErrorCode ierr;
119674ae819SStefano Zampini 
120674ae819SStefano Zampini   PetscFunctionBegin;
121674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
122674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
123674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
124674ae819SStefano Zampini   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
12568668abeSStefano Zampini   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
126674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
127457d4a33Sstefano_zampini   ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr);
128457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr);
129457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr);
130674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
131674ae819SStefano Zampini   PetscFunctionReturn(0);
132674ae819SStefano Zampini }
133674ae819SStefano Zampini 
134674ae819SStefano Zampini #undef __FUNCT__
135674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext"
136674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
137674ae819SStefano Zampini {
138674ae819SStefano Zampini   PetscErrorCode ierr;
139674ae819SStefano Zampini   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
140674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
141674ae819SStefano Zampini   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
142*be8bcbd0Sstefano_zampini #if defined(PETSC_USE_DEBUG)
143674ae819SStefano Zampini   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
144*be8bcbd0Sstefano_zampini #endif
145674ae819SStefano Zampini   MPI_Comm       comm;
14625d06dbeSstefano_zampini   Mat            ScalingMat,BD1,BD2;
147457d4a33Sstefano_zampini   Vec            fetidp_global;
148674ae819SStefano Zampini   IS             IS_l2g_lambda;
149a1c0d0daSstefano_zampini   IS             subset,subset_mult,subset_n,isvert;
150674ae819SStefano Zampini   PetscBool      skip_node,fully_redundant;
151674ae819SStefano Zampini   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
152dc456d91SStefano Zampini   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
15376ec1555SBarry Smith   PetscMPIInt    rank,size,buf_size,neigh;
154674ae819SStefano Zampini   PetscScalar    scalar_value;
155a1c0d0daSstefano_zampini   const PetscInt *vertex_indices;
156dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
157dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
158674ae819SStefano Zampini   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
159674ae819SStefano Zampini   PetscScalar    *array,*scaling_factors,*vals_B_delta;
16025d06dbeSstefano_zampini   PetscScalar    **all_factors;
161674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
162457d4a33Sstefano_zampini   PetscLayout    llay;
163a1c0d0daSstefano_zampini 
164457d4a33Sstefano_zampini   /* saddlepoint */
165457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
166457d4a33Sstefano_zampini   PetscLayout            play;
167457d4a33Sstefano_zampini   IS                     gP,pP;
168457d4a33Sstefano_zampini   PetscInt               nPl,nPg,nPgl;
169674ae819SStefano Zampini 
170674ae819SStefano Zampini   PetscFunctionBegin;
171674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
172674ae819SStefano Zampini   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
17376ec1555SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
174674ae819SStefano Zampini 
175457d4a33Sstefano_zampini   /* saddlepoint */
176457d4a33Sstefano_zampini   nPl      = 0;
177457d4a33Sstefano_zampini   nPg      = 0;
178457d4a33Sstefano_zampini   nPgl     = 0;
179457d4a33Sstefano_zampini   gP       = NULL;
180457d4a33Sstefano_zampini   pP       = NULL;
181457d4a33Sstefano_zampini   l2gmap_p = NULL;
182457d4a33Sstefano_zampini   play     = NULL;
183457d4a33Sstefano_zampini   ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr);
184022d8d2bSstefano_zampini   if (pP) { /* saddle point */
185457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
186457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr);
187457d4a33Sstefano_zampini     if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
188457d4a33Sstefano_zampini     ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr);
189457d4a33Sstefano_zampini     ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr);
190457d4a33Sstefano_zampini     ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr);
191457d4a33Sstefano_zampini     ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr);
192457d4a33Sstefano_zampini     ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr);
193457d4a33Sstefano_zampini 
194457d4a33Sstefano_zampini     /* interface pressure matrix */
195457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr);
196457d4a33Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */
197457d4a33Sstefano_zampini       IS Pg;
198457d4a33Sstefano_zampini 
199457d4a33Sstefano_zampini       ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr);
200457d4a33Sstefano_zampini       ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr);
201457d4a33Sstefano_zampini       ierr = ISDestroy(&Pg);CHKERRQ(ierr);
202457d4a33Sstefano_zampini       ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr);
203457d4a33Sstefano_zampini       ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr);
204457d4a33Sstefano_zampini       ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr);
205457d4a33Sstefano_zampini       ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr);
206457d4a33Sstefano_zampini       ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr);
207457d4a33Sstefano_zampini       ierr = PetscLayoutSetUp(play);CHKERRQ(ierr);
208457d4a33Sstefano_zampini     } else {
209457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr);
210457d4a33Sstefano_zampini       ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr);
211457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr);
212457d4a33Sstefano_zampini       ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr);
213457d4a33Sstefano_zampini       ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr);
214457d4a33Sstefano_zampini       ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr);
215457d4a33Sstefano_zampini       ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr);
216457d4a33Sstefano_zampini     }
217457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr);
218457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr);
219457d4a33Sstefano_zampini 
220457d4a33Sstefano_zampini     /* import matrices for interface pressures coupling */
221457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr);
222457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
223457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr);
224457d4a33Sstefano_zampini 
225457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr);
226457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
227457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr);
228457d4a33Sstefano_zampini 
229457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
230457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
231457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
232457d4a33Sstefano_zampini 
233457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
234457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
235457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
2366cc1294bSstefano_zampini 
2376cc1294bSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);CHKERRQ(ierr);
2386cc1294bSstefano_zampini     if (fetidpmat_ctx->rhs_flip) {
2396cc1294bSstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);CHKERRQ(ierr);
2406cc1294bSstefano_zampini     }
241457d4a33Sstefano_zampini   }
242457d4a33Sstefano_zampini 
243674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
244329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
245674ae819SStefano Zampini 
246674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
247674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
248674ae819SStefano Zampini   n_local_lambda = 0;
249674ae819SStefano Zampini   partial_sum = 0;
250674ae819SStefano Zampini   n_boundary_dofs = 0;
251674ae819SStefano Zampini   s = 0;
252a1c0d0daSstefano_zampini 
253674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
254a1c0d0daSstefano_zampini   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr);
255a1c0d0daSstefano_zampini   ierr = ISGetLocalSize(isvert,&n_vertices);CHKERRQ(ierr);
256a1c0d0daSstefano_zampini   ierr = ISGetIndices(isvert,&vertex_indices);CHKERRQ(ierr);
257a1c0d0daSstefano_zampini 
258674ae819SStefano Zampini   dual_size = pcis->n_B-n_vertices;
259785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
260785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
261785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
262674ae819SStefano Zampini 
263674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
264674ae819SStefano Zampini   for (i=0;i<pcis->n;i++){
265674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
2662d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
267674ae819SStefano Zampini     skip_node = PETSC_FALSE;
268674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
269674ae819SStefano Zampini       skip_node = PETSC_TRUE;
270674ae819SStefano Zampini       s++;
271674ae819SStefano Zampini     }
2722d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
2732d456bbaSstefano_zampini     if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
274674ae819SStefano Zampini     if (!skip_node) {
275674ae819SStefano Zampini       if (fully_redundant) {
276674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
277674ae819SStefano Zampini         n_lambda_for_dof = (j*(j+1))/2;
278674ae819SStefano Zampini       } else {
279674ae819SStefano Zampini         n_lambda_for_dof = j;
280674ae819SStefano Zampini       }
281674ae819SStefano Zampini       n_local_lambda += j;
282674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
283674ae819SStefano Zampini       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
284674ae819SStefano Zampini       /* store some data needed */
285674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
286674ae819SStefano Zampini       aux_local_numbering_1[partial_sum] = i;
287674ae819SStefano Zampini       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
288674ae819SStefano Zampini       partial_sum++;
289674ae819SStefano Zampini     }
290674ae819SStefano Zampini   }
291674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
292a1c0d0daSstefano_zampini   ierr = ISRestoreIndices(isvert,&vertex_indices);CHKERRQ(ierr);
293a1c0d0daSstefano_zampini   ierr = PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr);
2942d456bbaSstefano_zampini   dual_size = partial_sum;
295674ae819SStefano Zampini 
296674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
297dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
2983bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
299dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
300dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
3013d996552SStefano Zampini   ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr);
302dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
3033d996552SStefano Zampini 
3044fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG)
3054fe826edSStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
3064fe826edSStefano Zampini   ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
3074fe826edSStefano Zampini   ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
3084fe826edSStefano Zampini   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
3094fe826edSStefano Zampini   i = (PetscInt)PetscRealPart(scalar_value);
3106c4ed002SBarry 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);
3114fe826edSStefano Zampini #endif
312674ae819SStefano Zampini 
313674ae819SStefano Zampini   /* init data for scaling factors exchange */
31425d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
31525d06dbeSstefano_zampini     PetscInt    *ptrs_buffer,neigh_position;
31625d06dbeSstefano_zampini     PetscScalar *send_buffer,*recv_buffer;
31725d06dbeSstefano_zampini     MPI_Request *send_reqs,*recv_reqs;
31825d06dbeSstefano_zampini 
319674ae819SStefano Zampini     partial_sum = 0;
320785e854fSJed Brown     ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
3214b2aedd3SStefano Zampini     ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr);
3224b2aedd3SStefano Zampini     ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr);
32319c16490Sstefano_zampini     ierr = PetscMalloc1(pcis->n+1,&all_factors);CHKERRQ(ierr);
3244b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
325674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
326674ae819SStefano Zampini       partial_sum += pcis->n_shared[i];
327674ae819SStefano Zampini       ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
328674ae819SStefano Zampini     }
329785e854fSJed Brown     ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
330785e854fSJed Brown     ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
331785e854fSJed Brown     ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
332674ae819SStefano Zampini     for (i=0;i<pcis->n-1;i++) {
333674ae819SStefano Zampini       j = mat_graph->count[i];
334674ae819SStefano Zampini       all_factors[i+1]=all_factors[i]+j;
335674ae819SStefano Zampini     }
33625d06dbeSstefano_zampini 
337674ae819SStefano Zampini     /* scatter B scaling to N vec */
338674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
339674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
340674ae819SStefano Zampini     /* communications */
3412b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
342674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
343674ae819SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
344674ae819SStefano Zampini         send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
345674ae819SStefano Zampini       }
346674ae819SStefano Zampini       ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
347674ae819SStefano Zampini       ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
348674ae819SStefano Zampini       ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
349674ae819SStefano Zampini       ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
350674ae819SStefano Zampini     }
3512b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
3524b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) {
3534b2aedd3SStefano Zampini       ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3544b2aedd3SStefano Zampini     }
355674ae819SStefano Zampini     /* put values in correct places */
356674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
357674ae819SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
358674ae819SStefano Zampini         k = pcis->shared[i][j];
359674ae819SStefano Zampini         neigh_position = 0;
360674ae819SStefano Zampini         while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
361674ae819SStefano Zampini         all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
362674ae819SStefano Zampini       }
363674ae819SStefano Zampini     }
3644b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) {
3654b2aedd3SStefano Zampini       ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3664b2aedd3SStefano Zampini     }
367674ae819SStefano Zampini     ierr = PetscFree(send_reqs);CHKERRQ(ierr);
368674ae819SStefano Zampini     ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
369674ae819SStefano Zampini     ierr = PetscFree(send_buffer);CHKERRQ(ierr);
370674ae819SStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
371674ae819SStefano Zampini     ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
37225d06dbeSstefano_zampini   }
373674ae819SStefano Zampini 
374674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
375785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
376785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
377785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
378785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
37925d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
380785e854fSJed Brown     ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
38125d06dbeSstefano_zampini   } else {
38225d06dbeSstefano_zampini     scaling_factors = NULL;
38325d06dbeSstefano_zampini     all_factors     = NULL;
38425d06dbeSstefano_zampini   }
385dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
386674ae819SStefano Zampini   partial_sum=0;
387dc456d91SStefano Zampini   cum = 0;
388674ae819SStefano Zampini   for (i=0;i<dual_size;i++) {
389dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
390674ae819SStefano Zampini     j = mat_graph->count[aux_local_numbering_1[i]];
391674ae819SStefano Zampini     aux_sums[0]=0;
392674ae819SStefano Zampini     for (s=1;s<j;s++) {
393674ae819SStefano Zampini       aux_sums[s]=aux_sums[s-1]+j-s+1;
394674ae819SStefano Zampini     }
39525d06dbeSstefano_zampini     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
396674ae819SStefano Zampini     n_neg_values = 0;
3972a7da448SStefano Zampini     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
398674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
399674ae819SStefano Zampini     if (fully_redundant) {
400674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
401674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
402674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
403674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=-1.0;
40425d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s];
405674ae819SStefano Zampini       }
406674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
407674ae819SStefano Zampini         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
408674ae819SStefano Zampini         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
409674ae819SStefano Zampini         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
41025d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
411674ae819SStefano Zampini       }
412674ae819SStefano Zampini       partial_sum += j;
413674ae819SStefano Zampini     } else {
414674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
415674ae819SStefano Zampini       for (s=0;s<j;s++) {
416674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=n_global_lambda+s;
417674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
418674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=0.0;
419674ae819SStefano Zampini       }
420674ae819SStefano Zampini       /* B_delta */
421674ae819SStefano Zampini       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
422674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
423674ae819SStefano Zampini       }
424674ae819SStefano Zampini       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
425674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values]=1.0;
426674ae819SStefano Zampini       }
427674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999 */
42825d06dbeSstefano_zampini       if (!pcbddc->use_deluxe_scaling) {
429674ae819SStefano Zampini         for (s=0;s<n_neg_values;s++) {
430674ae819SStefano Zampini           scalar_value = 0.0;
43125d06dbeSstefano_zampini           for (k=0;k<s+1;k++) scalar_value += array[k];
432674ae819SStefano Zampini           scaling_factors[partial_sum+s] = -scalar_value;
433674ae819SStefano Zampini         }
434674ae819SStefano Zampini         for (s=0;s<n_pos_values;s++) {
435674ae819SStefano Zampini           scalar_value = 0.0;
43625d06dbeSstefano_zampini           for (k=s+n_neg_values;k<j;k++) scalar_value += array[k];
437674ae819SStefano Zampini           scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
438674ae819SStefano Zampini         }
43925d06dbeSstefano_zampini       }
440674ae819SStefano Zampini       partial_sum += j;
441674ae819SStefano Zampini     }
442dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
443674ae819SStefano Zampini   }
444dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
445dc456d91SStefano Zampini   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
446dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
447674ae819SStefano Zampini   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
448674ae819SStefano Zampini   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
449674ae819SStefano Zampini   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
45025d06dbeSstefano_zampini   if (all_factors) {
451674ae819SStefano Zampini     ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
452674ae819SStefano Zampini     ierr = PetscFree(all_factors);CHKERRQ(ierr);
45325d06dbeSstefano_zampini   }
454674ae819SStefano Zampini 
455674ae819SStefano Zampini   /* Create local part of B_delta */
456302440fdSBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
457674ae819SStefano Zampini   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
458674ae819SStefano Zampini   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
459674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
460674ae819SStefano Zampini   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
461674ae819SStefano Zampini   for (i=0;i<n_local_lambda;i++) {
462674ae819SStefano Zampini     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
463674ae819SStefano Zampini   }
464674ae819SStefano Zampini   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
465674ae819SStefano Zampini   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
466674ae819SStefano Zampini   ierr = MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
467674ae819SStefano Zampini 
46825d06dbeSstefano_zampini   BD1 = NULL;
46925d06dbeSstefano_zampini   BD2 = NULL;
470674ae819SStefano Zampini   if (fully_redundant) {
47125d06dbeSstefano_zampini     if (pcbddc->use_deluxe_scaling) SETERRQ(comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented");
472302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
473674ae819SStefano Zampini     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
474674ae819SStefano Zampini     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
475674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
476674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
477674ae819SStefano Zampini       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
478674ae819SStefano Zampini     }
479674ae819SStefano Zampini     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
480674ae819SStefano Zampini     ierr = MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481674ae819SStefano Zampini     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
482674ae819SStefano Zampini     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
483674ae819SStefano Zampini   } else {
484302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
485674ae819SStefano Zampini     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
48625d06dbeSstefano_zampini     if (!pcbddc->use_deluxe_scaling) {
487674ae819SStefano Zampini       ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
488674ae819SStefano Zampini       ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
489674ae819SStefano Zampini       for (i=0;i<n_local_lambda;i++) {
490674ae819SStefano Zampini         ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
491674ae819SStefano Zampini       }
492674ae819SStefano Zampini       ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
493674ae819SStefano Zampini       ierr = MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49425d06dbeSstefano_zampini     } else {
49525d06dbeSstefano_zampini       /* scaling as in Klawonn-Widlund 1999 */
49625d06dbeSstefano_zampini       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
49725d06dbeSstefano_zampini       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
49825d06dbeSstefano_zampini       Mat                 T;
49925d06dbeSstefano_zampini       PetscScalar         *W,lwork,*Bwork;
50025d06dbeSstefano_zampini       const PetscInt      *idxs;
50125d06dbeSstefano_zampini       PetscInt            cum,mss,*nnz;
50225d06dbeSstefano_zampini       PetscBLASInt        *pivots,B_lwork,B_N,B_ierr;
50325d06dbeSstefano_zampini 
50425d06dbeSstefano_zampini       if (!pcbddc->deluxe_singlemat) SETERRQ(comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
50525d06dbeSstefano_zampini       if (deluxe_ctx->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot compute B_Ddelta with deluxe scaling with active change context");
50625d06dbeSstefano_zampini 
50725d06dbeSstefano_zampini       mss  = 0;
50825d06dbeSstefano_zampini       ierr = PetscCalloc1(pcis->n_B,&nnz);CHKERRQ(ierr);
50925d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
51025d06dbeSstefano_zampini         ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
51125d06dbeSstefano_zampini         for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
51225d06dbeSstefano_zampini           PetscInt subset_size;
51325d06dbeSstefano_zampini 
51425d06dbeSstefano_zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
51525d06dbeSstefano_zampini           for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size;
51625d06dbeSstefano_zampini           mss  = PetscMax(mss,subset_size);
51725d06dbeSstefano_zampini           cum += subset_size;
51825d06dbeSstefano_zampini         }
51925d06dbeSstefano_zampini       }
52025d06dbeSstefano_zampini       ierr = MatCreate(PETSC_COMM_SELF,&T);CHKERRQ(ierr);
52125d06dbeSstefano_zampini       ierr = MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
52225d06dbeSstefano_zampini       ierr = MatSetType(T,MATSEQAIJ);CHKERRQ(ierr);
52325d06dbeSstefano_zampini       ierr = MatSeqAIJSetPreallocation(T,0,nnz);CHKERRQ(ierr);
52425d06dbeSstefano_zampini       ierr = PetscFree(nnz);CHKERRQ(ierr);
52525d06dbeSstefano_zampini 
52625d06dbeSstefano_zampini       /* workspace allocation */
5277f00194aSstefano_zampini       B_lwork = 0;
5287f00194aSstefano_zampini       if (mss) {
529*be8bcbd0Sstefano_zampini         PetscScalar dummy = 1;
530*be8bcbd0Sstefano_zampini 
53125d06dbeSstefano_zampini         B_lwork = -1;
53225d06dbeSstefano_zampini         ierr = PetscBLASIntCast(mss,&B_N);CHKERRQ(ierr);
53325d06dbeSstefano_zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
534*be8bcbd0Sstefano_zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr));
53525d06dbeSstefano_zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
53625d06dbeSstefano_zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
53725d06dbeSstefano_zampini         ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
5387f00194aSstefano_zampini       }
53925d06dbeSstefano_zampini       ierr = PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork);CHKERRQ(ierr);
54025d06dbeSstefano_zampini 
54125d06dbeSstefano_zampini       for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
54225d06dbeSstefano_zampini         PetscScalar  *M;
54325d06dbeSstefano_zampini         PetscInt     subset_size;
54425d06dbeSstefano_zampini 
54525d06dbeSstefano_zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
54625d06dbeSstefano_zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
54725d06dbeSstefano_zampini         ierr = MatDenseGetArray(deluxe_ctx->seq_mat[i],&M);CHKERRQ(ierr);
54825d06dbeSstefano_zampini         ierr = PetscMemcpy(W,M,subset_size*subset_size*sizeof(PetscScalar));CHKERRQ(ierr);
54925d06dbeSstefano_zampini         ierr = MatDenseRestoreArray(deluxe_ctx->seq_mat[i],&M);CHKERRQ(ierr);
55025d06dbeSstefano_zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
55125d06dbeSstefano_zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr));
55225d06dbeSstefano_zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
55325d06dbeSstefano_zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
55425d06dbeSstefano_zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
55525d06dbeSstefano_zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
55625d06dbeSstefano_zampini         ierr = MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES);CHKERRQ(ierr);
55725d06dbeSstefano_zampini         cum += subset_size;
55825d06dbeSstefano_zampini       }
55925d06dbeSstefano_zampini       ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56025d06dbeSstefano_zampini       ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56125d06dbeSstefano_zampini       ierr = MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1);CHKERRQ(ierr);
56225d06dbeSstefano_zampini       ierr = MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2);CHKERRQ(ierr);
56325d06dbeSstefano_zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
56425d06dbeSstefano_zampini       ierr = PetscFree3(W,pivots,Bwork);CHKERRQ(ierr);
56525d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
56625d06dbeSstefano_zampini         ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
56725d06dbeSstefano_zampini       }
56825d06dbeSstefano_zampini     }
569674ae819SStefano Zampini   }
570674ae819SStefano Zampini   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
571674ae819SStefano Zampini   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
572674ae819SStefano Zampini 
573457d4a33Sstefano_zampini   /* Layout of multipliers */
574457d4a33Sstefano_zampini   ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr);
575457d4a33Sstefano_zampini   ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr);
576457d4a33Sstefano_zampini   ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
577457d4a33Sstefano_zampini   ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr);
578457d4a33Sstefano_zampini   ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr);
579457d4a33Sstefano_zampini 
580457d4a33Sstefano_zampini   /* Local work vector of multipliers */
581457d4a33Sstefano_zampini   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
582457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
583457d4a33Sstefano_zampini   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
584457d4a33Sstefano_zampini 
58525d06dbeSstefano_zampini   if (BD2) {
58625d06dbeSstefano_zampini     ISLocalToGlobalMapping l2g;
58725d06dbeSstefano_zampini     Mat                    T,TA,*pT;
58825d06dbeSstefano_zampini     IS                     is;
58925d06dbeSstefano_zampini     PetscInt               nl,N;
59025d06dbeSstefano_zampini     BDdelta_DN             ctx;
59125d06dbeSstefano_zampini 
59225d06dbeSstefano_zampini     ierr = PetscLayoutGetLocalSize(llay,&nl);CHKERRQ(ierr);
59325d06dbeSstefano_zampini     ierr = PetscLayoutGetSize(llay,&N);CHKERRQ(ierr);
59425d06dbeSstefano_zampini     ierr = MatCreate(comm,&T);CHKERRQ(ierr);
59525d06dbeSstefano_zampini     ierr = MatSetSizes(T,nl,nl,N,N);CHKERRQ(ierr);
59625d06dbeSstefano_zampini     ierr = MatSetType(T,MATIS);CHKERRQ(ierr);
59725d06dbeSstefano_zampini     ierr = ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g);CHKERRQ(ierr);
59825d06dbeSstefano_zampini     ierr = MatSetLocalToGlobalMapping(T,l2g,l2g);CHKERRQ(ierr);
59925d06dbeSstefano_zampini     ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr);
60025d06dbeSstefano_zampini     ierr = MatISSetLocalMat(T,BD2);CHKERRQ(ierr);
60125d06dbeSstefano_zampini     ierr = MatDestroy(&BD2);CHKERRQ(ierr);
60225d06dbeSstefano_zampini     ierr = MatISGetMPIXAIJ(T,MAT_INITIAL_MATRIX,&TA);CHKERRQ(ierr);
60325d06dbeSstefano_zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
60425d06dbeSstefano_zampini     ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
60525d06dbeSstefano_zampini     ierr = MatGetSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT);CHKERRQ(ierr);
60625d06dbeSstefano_zampini     ierr = MatDestroy(&TA);CHKERRQ(ierr);
60725d06dbeSstefano_zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
60825d06dbeSstefano_zampini     BD2  = pT[0];
60925d06dbeSstefano_zampini     ierr = PetscFree(pT);CHKERRQ(ierr);
61025d06dbeSstefano_zampini 
61125d06dbeSstefano_zampini     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
61225d06dbeSstefano_zampini     ierr = PetscNew(&ctx);CHKERRQ(ierr);
61325d06dbeSstefano_zampini     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL);CHKERRQ(ierr);
61425d06dbeSstefano_zampini     ierr = MatShellSetContext(fetidpmat_ctx->B_Ddelta,(void *)ctx);CHKERRQ(ierr);
61525d06dbeSstefano_zampini     ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred);CHKERRQ(ierr);
61625d06dbeSstefano_zampini     ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred);CHKERRQ(ierr);
61725d06dbeSstefano_zampini     ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred);CHKERRQ(ierr);
61825d06dbeSstefano_zampini     ierr = MatSetUp(fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
61925d06dbeSstefano_zampini 
62025d06dbeSstefano_zampini     ierr = PetscObjectReference((PetscObject)BD1);CHKERRQ(ierr);
62125d06dbeSstefano_zampini     ctx->BD = BD1;
62225d06dbeSstefano_zampini     ierr = KSPCreate(PETSC_COMM_SELF,&ctx->kBD);CHKERRQ(ierr);
62325d06dbeSstefano_zampini     ierr = KSPSetOperators(ctx->kBD,BD2,BD2);CHKERRQ(ierr);
62425d06dbeSstefano_zampini     ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work);CHKERRQ(ierr);
62525d06dbeSstefano_zampini     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
62625d06dbeSstefano_zampini   }
62725d06dbeSstefano_zampini   ierr = MatDestroy(&BD1);CHKERRQ(ierr);
62825d06dbeSstefano_zampini   ierr = MatDestroy(&BD2);CHKERRQ(ierr);
62925d06dbeSstefano_zampini 
63025d06dbeSstefano_zampini   /* fetidpmat sizes */
63125d06dbeSstefano_zampini   fetidpmat_ctx->n += nPgl;
63225d06dbeSstefano_zampini   fetidpmat_ctx->N  = fetidpmat_ctx->n_lambda+nPg;
63325d06dbeSstefano_zampini 
634457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
635457d4a33Sstefano_zampini   ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr);
636457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr);
637457d4a33Sstefano_zampini   ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr);
638457d4a33Sstefano_zampini   ierr = VecSetUp(fetidp_global);CHKERRQ(ierr);
639457d4a33Sstefano_zampini 
6409eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
6419eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
642457d4a33Sstefano_zampini      of the FETI-DP linear system */
643457d4a33Sstefano_zampini   if (nPg) {
644af140850Sstefano_zampini     IS             IS_l2g_p,ais;
645457d4a33Sstefano_zampini     PetscLayout    alay;
646457d4a33Sstefano_zampini     const PetscInt *idxs,*pranges,*aranges,*lranges;
647af140850Sstefano_zampini     PetscInt       *l2g_indices_p,rst;
648457d4a33Sstefano_zampini 
649457d4a33Sstefano_zampini     ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr);
650457d4a33Sstefano_zampini     ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr);
651457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr);
652457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr);
653457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr);
654457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
655af140850Sstefano_zampini     /* shift local to global indices for pressure */
656457d4a33Sstefano_zampini     for (i=0;i<nPl;i++) {
657457d4a33Sstefano_zampini       PetscInt owner;
658457d4a33Sstefano_zampini 
659457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr);
660457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
661457d4a33Sstefano_zampini     }
662457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
663457d4a33Sstefano_zampini     ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr);
664af140850Sstefano_zampini 
665457d4a33Sstefano_zampini     /* local to global scatter for interface pressure */
666457d4a33Sstefano_zampini     ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr);
667457d4a33Sstefano_zampini     ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr);
668457d4a33Sstefano_zampini 
669af140850Sstefano_zampini     /* shift local to global indices for multipliers */
670457d4a33Sstefano_zampini     for (i=0;i<n_local_lambda;i++) {
671457d4a33Sstefano_zampini       PetscInt owner,ps;
672457d4a33Sstefano_zampini 
673457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr);
674457d4a33Sstefano_zampini       ps = pranges[owner+1]-pranges[owner];
675457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
676457d4a33Sstefano_zampini     }
677457d4a33Sstefano_zampini 
6789cc7774eSstefano_zampini     /* scatter from alldofs to interface pressures global fetidp vector */
6799cc7774eSstefano_zampini     ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr);
6809cc7774eSstefano_zampini     ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr);
681af140850Sstefano_zampini     ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr);
6829cc7774eSstefano_zampini     ierr = ISDestroy(&ais);CHKERRQ(ierr);
683457d4a33Sstefano_zampini   }
684457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr);
685457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr);
686457d4a33Sstefano_zampini   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
687a1c0d0daSstefano_zampini 
6889cc7774eSstefano_zampini   /* scatter from local to global multipliers */
689457d4a33Sstefano_zampini   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
690457d4a33Sstefano_zampini   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
691457d4a33Sstefano_zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr);
692a1c0d0daSstefano_zampini   ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr);
693457d4a33Sstefano_zampini 
694a1c0d0daSstefano_zampini   /* Create some work vectors needed by fetidp */
695674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
696674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
697674ae819SStefano Zampini   PetscFunctionReturn(0);
698674ae819SStefano Zampini }
699674ae819SStefano Zampini 
70025d06dbeSstefano_zampini 
701674ae819SStefano Zampini #undef __FUNCT__
702674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext"
703674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
704674ae819SStefano Zampini {
705674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
70625d06dbeSstefano_zampini   PC_BDDC        *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data;
70725d06dbeSstefano_zampini   PC_IS          *pcis = (PC_IS*)fetidppc_ctx->pc->data;
708f28b6018SStefano Zampini   PetscBool      lumped = PETSC_FALSE;
709674ae819SStefano Zampini   PetscErrorCode ierr;
710674ae819SStefano Zampini 
711674ae819SStefano Zampini   PetscFunctionBegin;
712674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
713674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
714674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
715674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
716674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
717674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
71825d06dbeSstefano_zampini   if (mat_ctx->deluxe_nonred) {
71925d06dbeSstefano_zampini     PC               pc,mpc;
72025d06dbeSstefano_zampini     BDdelta_DN       ctx;
72125d06dbeSstefano_zampini     MatSolverPackage solver;
72225d06dbeSstefano_zampini     const char       *prefix;
72325d06dbeSstefano_zampini 
72425d06dbeSstefano_zampini     ierr = MatShellGetContext(mat_ctx->B_Ddelta,&ctx);CHKERRQ(ierr);
72525d06dbeSstefano_zampini     ierr = KSPSetType(ctx->kBD,KSPPREONLY);CHKERRQ(ierr);
72625d06dbeSstefano_zampini     ierr = KSPGetPC(ctx->kBD,&mpc);CHKERRQ(ierr);
72725d06dbeSstefano_zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc);CHKERRQ(ierr);
72825d06dbeSstefano_zampini     ierr = PCSetType(mpc,PCLU);CHKERRQ(ierr);
72925d06dbeSstefano_zampini     ierr = PCFactorGetMatSolverPackage(pc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
73025d06dbeSstefano_zampini     if (solver) {
73125d06dbeSstefano_zampini       ierr = PCFactorSetMatSolverPackage(mpc,solver);CHKERRQ(ierr);
73225d06dbeSstefano_zampini     }
73325d06dbeSstefano_zampini     ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr);
73425d06dbeSstefano_zampini     ierr = KSPSetOptionsPrefix(ctx->kBD,prefix);CHKERRQ(ierr);
73525d06dbeSstefano_zampini     ierr = KSPAppendOptionsPrefix(ctx->kBD,"bddelta_");CHKERRQ(ierr);
73625d06dbeSstefano_zampini     ierr = KSPSetFromOptions(ctx->kBD);CHKERRQ(ierr);
73725d06dbeSstefano_zampini   }
73825d06dbeSstefano_zampini 
739674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
740674ae819SStefano Zampini   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
741f28b6018SStefano Zampini   /* Dirichlet preconditioner */
7429c2d02cdSstefano_zampini   ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL);CHKERRQ(ierr);
743f28b6018SStefano Zampini   if (!lumped) {
7449c2d02cdSstefano_zampini     IS        iP;
7459c2d02cdSstefano_zampini     PetscBool discrete_harmonic = PETSC_FALSE;
7469c2d02cdSstefano_zampini 
7479c2d02cdSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iP",(PetscObject*)&iP);CHKERRQ(ierr);
7489c2d02cdSstefano_zampini     if (iP) {
7499c2d02cdSstefano_zampini       ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL);CHKERRQ(ierr);
7509c2d02cdSstefano_zampini     }
7519c2d02cdSstefano_zampini     if (discrete_harmonic) {
7529c2d02cdSstefano_zampini       KSP        sksp;
7539c2d02cdSstefano_zampini       PC         pc;
7549c2d02cdSstefano_zampini       Mat        A_II,A_IB,A_BI;
7559c2d02cdSstefano_zampini       IS         aB;
7569c2d02cdSstefano_zampini       PetscInt   nb;
7575334bea6Sstefano_zampini       PetscBool  isshell;
7589c2d02cdSstefano_zampini       KSPType    ksptype;
7599c2d02cdSstefano_zampini       const char *prefix;
7609c2d02cdSstefano_zampini 
7619c2d02cdSstefano_zampini       /*
7629c2d02cdSstefano_zampini         We constructs a Schur complement for
7639c2d02cdSstefano_zampini 
7649c2d02cdSstefano_zampini         | A_II A_ID |
7659c2d02cdSstefano_zampini         | A_DI A_DD |
7669c2d02cdSstefano_zampini 
7679c2d02cdSstefano_zampini         instead of
7689c2d02cdSstefano_zampini 
7699c2d02cdSstefano_zampini         | A_II  B^t_II A_ID |
7709c2d02cdSstefano_zampini         | B_II -C_II   B_ID |
7719c2d02cdSstefano_zampini         | A_DI  B^t_ID A_DD |
7729c2d02cdSstefano_zampini 
7739c2d02cdSstefano_zampini       */
7749c2d02cdSstefano_zampini       ierr = ISGetLocalSize(pcis->is_B_local,&nb);CHKERRQ(ierr);
7759c2d02cdSstefano_zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB);CHKERRQ(ierr);
7769c2d02cdSstefano_zampini       ierr = MatGetSubMatrix(pcis->A_II,iP,iP,MAT_INITIAL_MATRIX,&A_II);CHKERRQ(ierr);
7779c2d02cdSstefano_zampini       ierr = MatGetSubMatrix(pcis->A_IB,iP,aB,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
7789c2d02cdSstefano_zampini       ierr = MatGetSubMatrix(pcis->A_BI,aB,iP,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
7799c2d02cdSstefano_zampini       ierr = MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
7809c2d02cdSstefano_zampini 
7819c2d02cdSstefano_zampini       /* propagate settings of solver */
7829c2d02cdSstefano_zampini       ierr = MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp);CHKERRQ(ierr);
7839c2d02cdSstefano_zampini       ierr = KSPGetType(pcis->ksp_D,&ksptype);CHKERRQ(ierr);
7849c2d02cdSstefano_zampini       ierr = KSPSetType(sksp,ksptype);CHKERRQ(ierr);
7859c2d02cdSstefano_zampini       ierr = KSPGetPC(pcis->ksp_D,&pc);CHKERRQ(ierr);
7865334bea6Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell);CHKERRQ(ierr);
7875334bea6Sstefano_zampini       if (!isshell) {
7885334bea6Sstefano_zampini         MatSolverPackage solver;
7895334bea6Sstefano_zampini         PCType           pctype;
7905334bea6Sstefano_zampini 
7919c2d02cdSstefano_zampini         ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
7929c2d02cdSstefano_zampini         ierr = PCFactorGetMatSolverPackage(pc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
7939c2d02cdSstefano_zampini         ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr);
7949c2d02cdSstefano_zampini         ierr = PCSetType(pc,pctype);CHKERRQ(ierr);
7959c2d02cdSstefano_zampini         if (solver) {
7969c2d02cdSstefano_zampini           ierr = PCFactorSetMatSolverPackage(pc,solver);CHKERRQ(ierr);
7979c2d02cdSstefano_zampini         }
7985334bea6Sstefano_zampini       } else {
7995334bea6Sstefano_zampini         ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr);
8005334bea6Sstefano_zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
8015334bea6Sstefano_zampini       }
8029c2d02cdSstefano_zampini       ierr = MatDestroy(&A_II);CHKERRQ(ierr);
8039c2d02cdSstefano_zampini       ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
8049c2d02cdSstefano_zampini       ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
8059c2d02cdSstefano_zampini       ierr = ISDestroy(&aB);CHKERRQ(ierr);
8069c2d02cdSstefano_zampini       ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr);
8079c2d02cdSstefano_zampini       ierr = KSPSetOptionsPrefix(sksp,prefix);CHKERRQ(ierr);
8089c2d02cdSstefano_zampini       ierr = KSPAppendOptionsPrefix(sksp,"harmonic_");CHKERRQ(ierr);
8093016320fSstefano_zampini       ierr = KSPSetFromOptions(sksp);CHKERRQ(ierr);
8109c2d02cdSstefano_zampini     } else { /* default Dirichlet preconditioner is pde-harmonic */
811ed6c3d69SStefano Zampini       ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
812ed6c3d69SStefano Zampini       ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
8139c2d02cdSstefano_zampini     }
814f28b6018SStefano Zampini   } else {
815f28b6018SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr);
816f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
817f28b6018SStefano Zampini   }
818af140850Sstefano_zampini   /* saddle-point */
819af140850Sstefano_zampini   if (mat_ctx->xPg) {
820af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr);
821af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
822af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr);
823af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
8246cc1294bSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PKSP",(PetscObject*)&fetidppc_ctx->kP);CHKERRQ(ierr);
8256cc1294bSstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidppc_ctx->kP);CHKERRQ(ierr);
826af140850Sstefano_zampini   }
827674ae819SStefano Zampini   PetscFunctionReturn(0);
828674ae819SStefano Zampini }
829674ae819SStefano Zampini 
830674ae819SStefano Zampini #undef __FUNCT__
831674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult"
832674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
833674ae819SStefano Zampini {
834674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
835617d11aeSStefano Zampini   PC_BDDC        *pcbddc;
836674ae819SStefano Zampini   PC_IS          *pcis;
837674ae819SStefano Zampini   PetscErrorCode ierr;
838674ae819SStefano Zampini 
839674ae819SStefano Zampini   PetscFunctionBegin;
840674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
841674ae819SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
842617d11aeSStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
843674ae819SStefano Zampini   /* Application of B_delta^T */
844af140850Sstefano_zampini   ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr);
845674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
846674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
847674ae819SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
848af140850Sstefano_zampini 
849af140850Sstefano_zampini   /* Add contribution from saddle point */
850af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
851af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
852af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
853af140850Sstefano_zampini     if (pcbddc->switch_static) {
854af140850Sstefano_zampini       ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr);
855af140850Sstefano_zampini     }
856af140850Sstefano_zampini     ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
857af140850Sstefano_zampini   } else {
858af140850Sstefano_zampini     if (pcbddc->switch_static) {
859674ae819SStefano Zampini       ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
860af140850Sstefano_zampini     }
861af140850Sstefano_zampini   }
862af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
863617d11aeSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
864dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
865c7ffc8ceSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
866af140850Sstefano_zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
867674ae819SStefano Zampini   /* Application of B_delta */
868674ae819SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
869af140850Sstefano_zampini   /* Contribution from boundary pressures */
870af140850Sstefano_zampini   if (mat_ctx->C) {
871af140850Sstefano_zampini     const PetscScalar *lx;
872af140850Sstefano_zampini     PetscScalar       *ly;
873af140850Sstefano_zampini 
874af140850Sstefano_zampini     /* pressures ordered first in x and y */
875af140850Sstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
876af140850Sstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
877af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr);
878af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr);
879af140850Sstefano_zampini     ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr);
880af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr);
881af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr);
882af140850Sstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
883af140850Sstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
884af140850Sstefano_zampini   }
885af140850Sstefano_zampini   /* Add contribution from saddle point */
886af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
887af140850Sstefano_zampini     ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
888af140850Sstefano_zampini     if (pcbddc->switch_static) {
889af140850Sstefano_zampini       ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
890af140850Sstefano_zampini     }
891af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
892af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
893af140850Sstefano_zampini   }
894674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
895674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
896674ae819SStefano Zampini   PetscFunctionReturn(0);
897674ae819SStefano Zampini }
898674ae819SStefano Zampini 
899674ae819SStefano Zampini #undef __FUNCT__
900edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose"
901edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
902edf7251bSStefano Zampini {
903edf7251bSStefano Zampini   FETIDPMat_ctx  mat_ctx;
904edf7251bSStefano Zampini   PC_IS          *pcis;
905edf7251bSStefano Zampini   PetscErrorCode ierr;
906edf7251bSStefano Zampini 
907edf7251bSStefano Zampini   PetscFunctionBegin;
908edf7251bSStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
909edf7251bSStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
910edf7251bSStefano Zampini   /* Application of B_delta^T */
911edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
912edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
913edf7251bSStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
914edf7251bSStefano Zampini   /* Application of \widetilde{S}^-1 */
915edf7251bSStefano Zampini   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
916edf7251bSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr);
917edf7251bSStefano Zampini   /* Application of B_delta */
918edf7251bSStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
919edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
920edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
921edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
922edf7251bSStefano Zampini   PetscFunctionReturn(0);
923edf7251bSStefano Zampini }
924edf7251bSStefano Zampini 
925edf7251bSStefano Zampini #undef __FUNCT__
926674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply"
927674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
928674ae819SStefano Zampini {
929674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
930674ae819SStefano Zampini   PC_IS          *pcis;
931674ae819SStefano Zampini   PetscErrorCode ierr;
932674ae819SStefano Zampini 
933674ae819SStefano Zampini   PetscFunctionBegin;
934302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
935674ae819SStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
936674ae819SStefano Zampini   /* Application of B_Ddelta^T */
937674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
938674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
939674ae819SStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
940674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
941ed6c3d69SStefano Zampini   /* Application of local Schur complement */
942ed6c3d69SStefano Zampini   ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
943edf7251bSStefano Zampini   /* Application of B_Ddelta */
944edf7251bSStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
945edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
946edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
947edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
948b519380fSstefano_zampini   /* interface pressure preconditioner */
949b519380fSstefano_zampini   if (pc_ctx->kP) {
950b519380fSstefano_zampini     const PetscScalar *lx;
951b519380fSstefano_zampini     PetscScalar *ly;
952b519380fSstefano_zampini 
953b519380fSstefano_zampini     /* pressures ordered first in x and y */
954b519380fSstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
955b519380fSstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
956b519380fSstefano_zampini     ierr = VecPlaceArray(pc_ctx->xPg,lx);CHKERRQ(ierr);
957b519380fSstefano_zampini     ierr = VecPlaceArray(pc_ctx->yPg,ly);CHKERRQ(ierr);
958b519380fSstefano_zampini     ierr = KSPSolve(pc_ctx->kP,pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr);
959b519380fSstefano_zampini     ierr = VecResetArray(pc_ctx->xPg);CHKERRQ(ierr);
960b519380fSstefano_zampini     ierr = VecResetArray(pc_ctx->yPg);CHKERRQ(ierr);
961b519380fSstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
962b519380fSstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
963b519380fSstefano_zampini   }
964edf7251bSStefano Zampini   PetscFunctionReturn(0);
965edf7251bSStefano Zampini }
966edf7251bSStefano Zampini 
967edf7251bSStefano Zampini #undef __FUNCT__
968edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose"
969edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
970edf7251bSStefano Zampini {
971edf7251bSStefano Zampini   FETIDPPC_ctx   pc_ctx;
972edf7251bSStefano Zampini   PC_IS          *pcis;
973edf7251bSStefano Zampini   PetscErrorCode ierr;
974edf7251bSStefano Zampini 
975edf7251bSStefano Zampini   PetscFunctionBegin;
976302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
977edf7251bSStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
978edf7251bSStefano Zampini   /* Application of B_Ddelta^T */
979edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
980edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
981edf7251bSStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
982edf7251bSStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
983ed6c3d69SStefano Zampini   /* Application of local Schur complement */
984ed6c3d69SStefano Zampini   ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
985674ae819SStefano Zampini   /* Application of B_Ddelta */
986674ae819SStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
987674ae819SStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
988674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
989674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
990674ae819SStefano Zampini   PetscFunctionReturn(0);
991674ae819SStefano Zampini }
992c45b8d2dSstefano_zampini 
993c45b8d2dSstefano_zampini #undef __FUNCT__
994c45b8d2dSstefano_zampini #define __FUNCT__ "FETIDPPCView"
995c45b8d2dSstefano_zampini PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
996c45b8d2dSstefano_zampini {
997c45b8d2dSstefano_zampini   FETIDPPC_ctx      pc_ctx;
998c45b8d2dSstefano_zampini   PetscBool         iascii;
999c45b8d2dSstefano_zampini   PetscViewer       sviewer;
1000c45b8d2dSstefano_zampini   PetscErrorCode    ierr;
1001c45b8d2dSstefano_zampini 
1002c45b8d2dSstefano_zampini   PetscFunctionBegin;
1003c45b8d2dSstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1004c45b8d2dSstefano_zampini   if (iascii) {
1005c45b8d2dSstefano_zampini     PetscMPIInt rank;
100625d06dbeSstefano_zampini     PetscBool   isschur,isshell;
1007c45b8d2dSstefano_zampini 
1008c45b8d2dSstefano_zampini     ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
1009c45b8d2dSstefano_zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
1010c45b8d2dSstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur);CHKERRQ(ierr);
1011c45b8d2dSstefano_zampini     if (isschur) {
1012c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP multipliers Dirichlet preconditioner (just from rank 0)\n");CHKERRQ(ierr);
1013c45b8d2dSstefano_zampini     } else {
1014c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP multipliers Lumped preconditioner (just from rank 0)\n");CHKERRQ(ierr);
1015c45b8d2dSstefano_zampini     }
1016c45b8d2dSstefano_zampini     ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
1017c45b8d2dSstefano_zampini     if (!rank) {
1018c45b8d2dSstefano_zampini       ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1019c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIAddTab(sviewer,2);CHKERRQ(ierr);
1020c45b8d2dSstefano_zampini       ierr = MatView(pc_ctx->S_j,sviewer);CHKERRQ(ierr);
1021c45b8d2dSstefano_zampini       ierr = PetscViewerASCIISubtractTab(sviewer,2);CHKERRQ(ierr);
1022c45b8d2dSstefano_zampini       ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr);
1023c45b8d2dSstefano_zampini     }
102425d06dbeSstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell);CHKERRQ(ierr);
102525d06dbeSstefano_zampini     if (isshell) {
102625d06dbeSstefano_zampini       BDdelta_DN ctx;
102725d06dbeSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n");CHKERRQ(ierr);
102825d06dbeSstefano_zampini       ierr = MatShellGetContext(pc_ctx->B_Ddelta,&ctx);CHKERRQ(ierr);
102925d06dbeSstefano_zampini       if (!rank) {
103025d06dbeSstefano_zampini         ierr = PetscViewerASCIIAddTab(sviewer,2);CHKERRQ(ierr);
103125d06dbeSstefano_zampini         ierr = KSPView(ctx->kBD,sviewer);CHKERRQ(ierr);
103225d06dbeSstefano_zampini         ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
103325d06dbeSstefano_zampini         ierr = MatView(ctx->BD,sviewer);CHKERRQ(ierr);
103425d06dbeSstefano_zampini         ierr = PetscViewerASCIISubtractTab(sviewer,2);CHKERRQ(ierr);
103525d06dbeSstefano_zampini         ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr);
103625d06dbeSstefano_zampini       }
103725d06dbeSstefano_zampini     }
1038c45b8d2dSstefano_zampini     ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
1039c45b8d2dSstefano_zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1040c45b8d2dSstefano_zampini     if (pc_ctx->kP) {
1041c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIPrintf(viewer,"  FETI-DP pressure preconditioner\n");CHKERRQ(ierr);
1042c45b8d2dSstefano_zampini       ierr = PetscViewerASCIIAddTab(viewer,2);CHKERRQ(ierr);
1043c45b8d2dSstefano_zampini       ierr = KSPView(pc_ctx->kP,viewer);CHKERRQ(ierr);
1044c45b8d2dSstefano_zampini       ierr = PetscViewerASCIISubtractTab(viewer,2);CHKERRQ(ierr);
1045c45b8d2dSstefano_zampini     }
1046c45b8d2dSstefano_zampini   }
1047c45b8d2dSstefano_zampini   PetscFunctionReturn(0);
1048c45b8d2dSstefano_zampini }
1049