xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision e5afcf2835ad2df3c79a70d4d9a0fbb86e97247e)
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 static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
625d06dbeSstefano_zampini {
725d06dbeSstefano_zampini   BDdelta_DN     ctx;
825d06dbeSstefano_zampini   PetscErrorCode ierr;
925d06dbeSstefano_zampini 
1025d06dbeSstefano_zampini   PetscFunctionBegin;
1125d06dbeSstefano_zampini   ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
1225d06dbeSstefano_zampini   ierr = MatMultTranspose(ctx->BD,x,ctx->work);CHKERRQ(ierr);
1325d06dbeSstefano_zampini   ierr = KSPSolveTranspose(ctx->kBD,ctx->work,y);CHKERRQ(ierr);
14c0decd05SBarry Smith   /* No PC so cannot propagate up failure in KSPSolveTranspose() */
1525d06dbeSstefano_zampini   PetscFunctionReturn(0);
1625d06dbeSstefano_zampini }
1725d06dbeSstefano_zampini 
1825d06dbeSstefano_zampini static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
1925d06dbeSstefano_zampini {
2025d06dbeSstefano_zampini   BDdelta_DN     ctx;
2125d06dbeSstefano_zampini   PetscErrorCode ierr;
2225d06dbeSstefano_zampini 
2325d06dbeSstefano_zampini   PetscFunctionBegin;
2425d06dbeSstefano_zampini   ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
2525d06dbeSstefano_zampini   ierr = KSPSolve(ctx->kBD,x,ctx->work);CHKERRQ(ierr);
26c0decd05SBarry Smith   /* No PC so cannot propagate up failure in KSPSolve() */
2725d06dbeSstefano_zampini   ierr = MatMult(ctx->BD,ctx->work,y);CHKERRQ(ierr);
2825d06dbeSstefano_zampini   PetscFunctionReturn(0);
2925d06dbeSstefano_zampini }
3025d06dbeSstefano_zampini 
3125d06dbeSstefano_zampini static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
3225d06dbeSstefano_zampini {
3325d06dbeSstefano_zampini   BDdelta_DN     ctx;
3425d06dbeSstefano_zampini   PetscErrorCode ierr;
3525d06dbeSstefano_zampini 
3625d06dbeSstefano_zampini   PetscFunctionBegin;
3725d06dbeSstefano_zampini   ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
3825d06dbeSstefano_zampini   ierr = MatDestroy(&ctx->BD);CHKERRQ(ierr);
3925d06dbeSstefano_zampini   ierr = KSPDestroy(&ctx->kBD);CHKERRQ(ierr);
4025d06dbeSstefano_zampini   ierr = VecDestroy(&ctx->work);CHKERRQ(ierr);
4125d06dbeSstefano_zampini   ierr = PetscFree(ctx);CHKERRQ(ierr);
4225d06dbeSstefano_zampini   PetscFunctionReturn(0);
4325d06dbeSstefano_zampini }
4425d06dbeSstefano_zampini 
45674ae819SStefano Zampini 
46674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
47674ae819SStefano Zampini {
48674ae819SStefano Zampini   FETIDPMat_ctx  newctx;
49674ae819SStefano Zampini   PetscErrorCode ierr;
50674ae819SStefano Zampini 
51674ae819SStefano Zampini   PetscFunctionBegin;
52854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
53674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
54674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
55674ae819SStefano Zampini   newctx->pc              = pc;
56674ae819SStefano Zampini   *fetidpmat_ctx          = newctx;
57674ae819SStefano Zampini   PetscFunctionReturn(0);
58674ae819SStefano Zampini }
59674ae819SStefano Zampini 
60674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
61674ae819SStefano Zampini {
62674ae819SStefano Zampini   FETIDPPC_ctx   newctx;
63674ae819SStefano Zampini   PetscErrorCode ierr;
64674ae819SStefano Zampini 
65674ae819SStefano Zampini   PetscFunctionBegin;
66854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
67674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
68674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
69674ae819SStefano Zampini   newctx->pc              = pc;
70674ae819SStefano Zampini   *fetidppc_ctx           = newctx;
71674ae819SStefano Zampini   PetscFunctionReturn(0);
72674ae819SStefano Zampini }
73674ae819SStefano Zampini 
74674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
75674ae819SStefano Zampini {
76674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
77674ae819SStefano Zampini   PetscErrorCode ierr;
78674ae819SStefano Zampini 
79674ae819SStefano Zampini   PetscFunctionBegin;
80674ae819SStefano Zampini   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
81674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
82674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
83674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
84674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
85674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
86457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr);
87457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr);
88457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr);
89457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr);
90457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr);
916cc1294bSstefano_zampini   ierr = VecDestroy(&mat_ctx->rhs_flip);CHKERRQ(ierr);
92457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr);
93457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr);
94457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr);
95674ae819SStefano Zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
96e1214c54Sstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda_only);CHKERRQ(ierr);
97457d4a33Sstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr);
989cc7774eSstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr);
99674ae819SStefano Zampini   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
100e1214c54Sstefano_zampini   ierr = ISDestroy(&mat_ctx->pressure);CHKERRQ(ierr);
101e1214c54Sstefano_zampini   ierr = ISDestroy(&mat_ctx->lagrange);CHKERRQ(ierr);
102674ae819SStefano Zampini   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
103674ae819SStefano Zampini   PetscFunctionReturn(0);
104674ae819SStefano Zampini }
105674ae819SStefano Zampini 
106674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
107674ae819SStefano Zampini {
108674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
109674ae819SStefano Zampini   PetscErrorCode ierr;
110674ae819SStefano Zampini 
111674ae819SStefano Zampini   PetscFunctionBegin;
112674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
113674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
114674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
115674ae819SStefano Zampini   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
11668668abeSStefano Zampini   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
117674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
118457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr);
119457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr);
120674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
121674ae819SStefano Zampini   PetscFunctionReturn(0);
122674ae819SStefano Zampini }
123674ae819SStefano Zampini 
124674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
125674ae819SStefano Zampini {
126674ae819SStefano Zampini   PetscErrorCode ierr;
127674ae819SStefano Zampini   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
128674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
129674ae819SStefano Zampini   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
130be8bcbd0Sstefano_zampini #if defined(PETSC_USE_DEBUG)
131674ae819SStefano Zampini   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
132be8bcbd0Sstefano_zampini #endif
133674ae819SStefano Zampini   MPI_Comm       comm;
13425d06dbeSstefano_zampini   Mat            ScalingMat,BD1,BD2;
135457d4a33Sstefano_zampini   Vec            fetidp_global;
136674ae819SStefano Zampini   IS             IS_l2g_lambda;
137a1c0d0daSstefano_zampini   IS             subset,subset_mult,subset_n,isvert;
138674ae819SStefano Zampini   PetscBool      skip_node,fully_redundant;
139674ae819SStefano Zampini   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
140dc456d91SStefano Zampini   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
14176ec1555SBarry Smith   PetscMPIInt    rank,size,buf_size,neigh;
142674ae819SStefano Zampini   PetscScalar    scalar_value;
143a1c0d0daSstefano_zampini   const PetscInt *vertex_indices;
144dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
145dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
146674ae819SStefano Zampini   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
147674ae819SStefano Zampini   PetscScalar    *array,*scaling_factors,*vals_B_delta;
14825d06dbeSstefano_zampini   PetscScalar    **all_factors;
149674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
150457d4a33Sstefano_zampini   PetscLayout    llay;
151a1c0d0daSstefano_zampini 
152457d4a33Sstefano_zampini   /* saddlepoint */
153457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
154457d4a33Sstefano_zampini   PetscLayout            play;
155457d4a33Sstefano_zampini   IS                     gP,pP;
156457d4a33Sstefano_zampini   PetscInt               nPl,nPg,nPgl;
157674ae819SStefano Zampini 
158674ae819SStefano Zampini   PetscFunctionBegin;
159674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
160674ae819SStefano Zampini   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
16176ec1555SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
162674ae819SStefano Zampini 
163457d4a33Sstefano_zampini   /* saddlepoint */
164457d4a33Sstefano_zampini   nPl      = 0;
165457d4a33Sstefano_zampini   nPg      = 0;
166457d4a33Sstefano_zampini   nPgl     = 0;
167457d4a33Sstefano_zampini   gP       = NULL;
168457d4a33Sstefano_zampini   pP       = NULL;
169457d4a33Sstefano_zampini   l2gmap_p = NULL;
170457d4a33Sstefano_zampini   play     = NULL;
171457d4a33Sstefano_zampini   ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr);
172022d8d2bSstefano_zampini   if (pP) { /* saddle point */
173457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
174457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr);
175457d4a33Sstefano_zampini     if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
176457d4a33Sstefano_zampini     ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr);
177457d4a33Sstefano_zampini     ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr);
178457d4a33Sstefano_zampini     ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr);
179457d4a33Sstefano_zampini     ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr);
180457d4a33Sstefano_zampini     ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr);
181457d4a33Sstefano_zampini 
182e1214c54Sstefano_zampini     /* pressure matrix */
183457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr);
184e1214c54Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
185457d4a33Sstefano_zampini       IS Pg;
186457d4a33Sstefano_zampini 
187457d4a33Sstefano_zampini       ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr);
188457d4a33Sstefano_zampini       ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr);
189457d4a33Sstefano_zampini       ierr = ISDestroy(&Pg);CHKERRQ(ierr);
190457d4a33Sstefano_zampini       ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr);
191457d4a33Sstefano_zampini       ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr);
192457d4a33Sstefano_zampini       ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr);
193457d4a33Sstefano_zampini       ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr);
194457d4a33Sstefano_zampini       ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr);
195457d4a33Sstefano_zampini       ierr = PetscLayoutSetUp(play);CHKERRQ(ierr);
196457d4a33Sstefano_zampini     } else {
197457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr);
198457d4a33Sstefano_zampini       ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr);
199457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr);
200457d4a33Sstefano_zampini       ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr);
201457d4a33Sstefano_zampini       ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr);
202457d4a33Sstefano_zampini       ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr);
203457d4a33Sstefano_zampini       ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr);
204457d4a33Sstefano_zampini     }
205457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr);
206457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr);
207457d4a33Sstefano_zampini 
208e1214c54Sstefano_zampini     /* import matrices for pressures coupling */
209457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr);
210457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
211457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr);
212457d4a33Sstefano_zampini 
213457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr);
214457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
215457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr);
216457d4a33Sstefano_zampini 
217457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
218457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
219457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
220457d4a33Sstefano_zampini 
221457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
222457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
223457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
2246cc1294bSstefano_zampini 
2256cc1294bSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);CHKERRQ(ierr);
2266cc1294bSstefano_zampini     if (fetidpmat_ctx->rhs_flip) {
2276cc1294bSstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);CHKERRQ(ierr);
2286cc1294bSstefano_zampini     }
229457d4a33Sstefano_zampini   }
230457d4a33Sstefano_zampini 
231674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
232329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
233674ae819SStefano Zampini 
234674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
235674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
236674ae819SStefano Zampini   n_local_lambda = 0;
237674ae819SStefano Zampini   partial_sum = 0;
238674ae819SStefano Zampini   n_boundary_dofs = 0;
239674ae819SStefano Zampini   s = 0;
240a1c0d0daSstefano_zampini 
241674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
242a1c0d0daSstefano_zampini   ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr);
243a1c0d0daSstefano_zampini   ierr = ISGetLocalSize(isvert,&n_vertices);CHKERRQ(ierr);
244a1c0d0daSstefano_zampini   ierr = ISGetIndices(isvert,&vertex_indices);CHKERRQ(ierr);
245a1c0d0daSstefano_zampini 
246674ae819SStefano Zampini   dual_size = pcis->n_B-n_vertices;
247785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
248785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
249785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
250674ae819SStefano Zampini 
251674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
252674ae819SStefano Zampini   for (i=0;i<pcis->n;i++){
253674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
2542d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
255674ae819SStefano Zampini     skip_node = PETSC_FALSE;
256674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
257674ae819SStefano Zampini       skip_node = PETSC_TRUE;
258674ae819SStefano Zampini       s++;
259674ae819SStefano Zampini     }
2602d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
2612d456bbaSstefano_zampini     if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
262674ae819SStefano Zampini     if (!skip_node) {
263674ae819SStefano Zampini       if (fully_redundant) {
264674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
265674ae819SStefano Zampini         n_lambda_for_dof = (j*(j+1))/2;
266674ae819SStefano Zampini       } else {
267674ae819SStefano Zampini         n_lambda_for_dof = j;
268674ae819SStefano Zampini       }
269674ae819SStefano Zampini       n_local_lambda += j;
270674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
271674ae819SStefano Zampini       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
272674ae819SStefano Zampini       /* store some data needed */
273674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
274674ae819SStefano Zampini       aux_local_numbering_1[partial_sum] = i;
275674ae819SStefano Zampini       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
276674ae819SStefano Zampini       partial_sum++;
277674ae819SStefano Zampini     }
278674ae819SStefano Zampini   }
279674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
280a1c0d0daSstefano_zampini   ierr = ISRestoreIndices(isvert,&vertex_indices);CHKERRQ(ierr);
281a1c0d0daSstefano_zampini   ierr = PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr);
2822d456bbaSstefano_zampini   dual_size = partial_sum;
283674ae819SStefano Zampini 
284674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
285dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
2863bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
287dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
288dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
2893d996552SStefano Zampini   ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr);
290dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
2913d996552SStefano Zampini 
2924fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG)
2934fe826edSStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
2944fe826edSStefano Zampini   ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2954fe826edSStefano Zampini   ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2964fe826edSStefano Zampini   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
2974fe826edSStefano Zampini   i = (PetscInt)PetscRealPart(scalar_value);
2986080607fSStefano Zampini   if (i != fetidpmat_ctx->n_lambda) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global number of multipliers mismatch! (%D != %D)",fetidpmat_ctx->n_lambda,i);
2994fe826edSStefano Zampini #endif
300674ae819SStefano Zampini 
301674ae819SStefano Zampini   /* init data for scaling factors exchange */
30225d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
30325d06dbeSstefano_zampini     PetscInt    *ptrs_buffer,neigh_position;
30425d06dbeSstefano_zampini     PetscScalar *send_buffer,*recv_buffer;
30525d06dbeSstefano_zampini     MPI_Request *send_reqs,*recv_reqs;
30625d06dbeSstefano_zampini 
307674ae819SStefano Zampini     partial_sum = 0;
308785e854fSJed Brown     ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
3094b2aedd3SStefano Zampini     ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr);
3104b2aedd3SStefano Zampini     ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr);
31119c16490Sstefano_zampini     ierr = PetscMalloc1(pcis->n+1,&all_factors);CHKERRQ(ierr);
3124b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
313674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
314674ae819SStefano Zampini       partial_sum += pcis->n_shared[i];
315674ae819SStefano Zampini       ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
316674ae819SStefano Zampini     }
317785e854fSJed Brown     ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
318785e854fSJed Brown     ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
319785e854fSJed Brown     ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
320674ae819SStefano Zampini     for (i=0;i<pcis->n-1;i++) {
321674ae819SStefano Zampini       j = mat_graph->count[i];
322674ae819SStefano Zampini       all_factors[i+1]=all_factors[i]+j;
323674ae819SStefano Zampini     }
32425d06dbeSstefano_zampini 
325674ae819SStefano Zampini     /* scatter B scaling to N vec */
326674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
327674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
328674ae819SStefano Zampini     /* communications */
3292b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
330674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
331674ae819SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
332674ae819SStefano Zampini         send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
333674ae819SStefano Zampini       }
334674ae819SStefano Zampini       ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
335674ae819SStefano Zampini       ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
336674ae819SStefano Zampini       ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
337674ae819SStefano Zampini       ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
338674ae819SStefano Zampini     }
3392b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
3404b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) {
3414b2aedd3SStefano Zampini       ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3424b2aedd3SStefano Zampini     }
343674ae819SStefano Zampini     /* put values in correct places */
344674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
345674ae819SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
346674ae819SStefano Zampini         k = pcis->shared[i][j];
347674ae819SStefano Zampini         neigh_position = 0;
348674ae819SStefano Zampini         while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
349674ae819SStefano Zampini         all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
350674ae819SStefano Zampini       }
351674ae819SStefano Zampini     }
3524b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) {
3534b2aedd3SStefano Zampini       ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3544b2aedd3SStefano Zampini     }
355674ae819SStefano Zampini     ierr = PetscFree(send_reqs);CHKERRQ(ierr);
356674ae819SStefano Zampini     ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
357674ae819SStefano Zampini     ierr = PetscFree(send_buffer);CHKERRQ(ierr);
358674ae819SStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
359674ae819SStefano Zampini     ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
36025d06dbeSstefano_zampini   }
361674ae819SStefano Zampini 
362674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
363785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
364785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
365785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
366785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
36725d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
368785e854fSJed Brown     ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
36925d06dbeSstefano_zampini   } else {
37025d06dbeSstefano_zampini     scaling_factors = NULL;
37125d06dbeSstefano_zampini     all_factors     = NULL;
37225d06dbeSstefano_zampini   }
373dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
374674ae819SStefano Zampini   partial_sum=0;
375dc456d91SStefano Zampini   cum = 0;
376674ae819SStefano Zampini   for (i=0;i<dual_size;i++) {
377dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
378674ae819SStefano Zampini     j = mat_graph->count[aux_local_numbering_1[i]];
379674ae819SStefano Zampini     aux_sums[0]=0;
380674ae819SStefano Zampini     for (s=1;s<j;s++) {
381674ae819SStefano Zampini       aux_sums[s]=aux_sums[s-1]+j-s+1;
382674ae819SStefano Zampini     }
38325d06dbeSstefano_zampini     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
384674ae819SStefano Zampini     n_neg_values = 0;
3852a7da448SStefano Zampini     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
386674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
387674ae819SStefano Zampini     if (fully_redundant) {
388674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
389674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
390674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
391674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=-1.0;
39225d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s];
393674ae819SStefano Zampini       }
394674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
395674ae819SStefano Zampini         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
396674ae819SStefano Zampini         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
397674ae819SStefano Zampini         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
39825d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
399674ae819SStefano Zampini       }
400674ae819SStefano Zampini       partial_sum += j;
401674ae819SStefano Zampini     } else {
402674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
403674ae819SStefano Zampini       for (s=0;s<j;s++) {
404674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=n_global_lambda+s;
405674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
406674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=0.0;
407674ae819SStefano Zampini       }
408674ae819SStefano Zampini       /* B_delta */
409674ae819SStefano Zampini       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
410674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
411674ae819SStefano Zampini       }
412674ae819SStefano Zampini       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
413674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values]=1.0;
414674ae819SStefano Zampini       }
415674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999 */
41625d06dbeSstefano_zampini       if (!pcbddc->use_deluxe_scaling) {
417674ae819SStefano Zampini         for (s=0;s<n_neg_values;s++) {
418674ae819SStefano Zampini           scalar_value = 0.0;
41925d06dbeSstefano_zampini           for (k=0;k<s+1;k++) scalar_value += array[k];
420674ae819SStefano Zampini           scaling_factors[partial_sum+s] = -scalar_value;
421674ae819SStefano Zampini         }
422674ae819SStefano Zampini         for (s=0;s<n_pos_values;s++) {
423674ae819SStefano Zampini           scalar_value = 0.0;
42425d06dbeSstefano_zampini           for (k=s+n_neg_values;k<j;k++) scalar_value += array[k];
425674ae819SStefano Zampini           scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
426674ae819SStefano Zampini         }
42725d06dbeSstefano_zampini       }
428674ae819SStefano Zampini       partial_sum += j;
429674ae819SStefano Zampini     }
430dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
431674ae819SStefano Zampini   }
432dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
433dc456d91SStefano Zampini   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
434dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
435674ae819SStefano Zampini   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
436674ae819SStefano Zampini   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
437674ae819SStefano Zampini   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
43825d06dbeSstefano_zampini   if (all_factors) {
439674ae819SStefano Zampini     ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
440674ae819SStefano Zampini     ierr = PetscFree(all_factors);CHKERRQ(ierr);
44125d06dbeSstefano_zampini   }
442674ae819SStefano Zampini 
443674ae819SStefano Zampini   /* Create local part of B_delta */
444302440fdSBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
445674ae819SStefano Zampini   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
446674ae819SStefano Zampini   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
447674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
448674ae819SStefano Zampini   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
449674ae819SStefano Zampini   for (i=0;i<n_local_lambda;i++) {
450674ae819SStefano Zampini     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
451674ae819SStefano Zampini   }
452674ae819SStefano Zampini   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
453674ae819SStefano Zampini   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
454674ae819SStefano Zampini   ierr = MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
455674ae819SStefano Zampini 
45625d06dbeSstefano_zampini   BD1 = NULL;
45725d06dbeSstefano_zampini   BD2 = NULL;
458674ae819SStefano Zampini   if (fully_redundant) {
45925d06dbeSstefano_zampini     if (pcbddc->use_deluxe_scaling) SETERRQ(comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented");
460302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
461674ae819SStefano Zampini     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
462674ae819SStefano Zampini     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
463674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
464674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
465674ae819SStefano Zampini       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
466674ae819SStefano Zampini     }
467674ae819SStefano Zampini     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
468674ae819SStefano Zampini     ierr = MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
469674ae819SStefano Zampini     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
470674ae819SStefano Zampini     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
471674ae819SStefano Zampini   } else {
472302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
473674ae819SStefano Zampini     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
474524221abSStefano Zampini     if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
475674ae819SStefano Zampini       ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
476674ae819SStefano Zampini       ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
477674ae819SStefano Zampini       for (i=0;i<n_local_lambda;i++) {
478674ae819SStefano Zampini         ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
479674ae819SStefano Zampini       }
480674ae819SStefano Zampini       ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481674ae819SStefano Zampini       ierr = MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48225d06dbeSstefano_zampini     } else {
48325d06dbeSstefano_zampini       /* scaling as in Klawonn-Widlund 1999 */
48425d06dbeSstefano_zampini       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
48525d06dbeSstefano_zampini       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
48625d06dbeSstefano_zampini       Mat                 T;
48725d06dbeSstefano_zampini       PetscScalar         *W,lwork,*Bwork;
488451a39c7SStefano Zampini       const PetscInt      *idxs = NULL;
48925d06dbeSstefano_zampini       PetscInt            cum,mss,*nnz;
49025d06dbeSstefano_zampini       PetscBLASInt        *pivots,B_lwork,B_N,B_ierr;
49125d06dbeSstefano_zampini 
49225d06dbeSstefano_zampini       if (!pcbddc->deluxe_singlemat) SETERRQ(comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
49325d06dbeSstefano_zampini       mss  = 0;
49425d06dbeSstefano_zampini       ierr = PetscCalloc1(pcis->n_B,&nnz);CHKERRQ(ierr);
49525d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
49625d06dbeSstefano_zampini         ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
49725d06dbeSstefano_zampini         for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
49825d06dbeSstefano_zampini           PetscInt subset_size;
49925d06dbeSstefano_zampini 
50025d06dbeSstefano_zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
50125d06dbeSstefano_zampini           for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size;
50225d06dbeSstefano_zampini           mss  = PetscMax(mss,subset_size);
50325d06dbeSstefano_zampini           cum += subset_size;
50425d06dbeSstefano_zampini         }
50525d06dbeSstefano_zampini       }
50625d06dbeSstefano_zampini       ierr = MatCreate(PETSC_COMM_SELF,&T);CHKERRQ(ierr);
50725d06dbeSstefano_zampini       ierr = MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr);
50825d06dbeSstefano_zampini       ierr = MatSetType(T,MATSEQAIJ);CHKERRQ(ierr);
50925d06dbeSstefano_zampini       ierr = MatSeqAIJSetPreallocation(T,0,nnz);CHKERRQ(ierr);
51025d06dbeSstefano_zampini       ierr = PetscFree(nnz);CHKERRQ(ierr);
51125d06dbeSstefano_zampini 
51225d06dbeSstefano_zampini       /* workspace allocation */
5137f00194aSstefano_zampini       B_lwork = 0;
5147f00194aSstefano_zampini       if (mss) {
515be8bcbd0Sstefano_zampini         PetscScalar dummy = 1;
516be8bcbd0Sstefano_zampini 
51725d06dbeSstefano_zampini         B_lwork = -1;
51825d06dbeSstefano_zampini         ierr = PetscBLASIntCast(mss,&B_N);CHKERRQ(ierr);
51925d06dbeSstefano_zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
520be8bcbd0Sstefano_zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr));
52125d06dbeSstefano_zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
52225d06dbeSstefano_zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
52325d06dbeSstefano_zampini         ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
5247f00194aSstefano_zampini       }
52525d06dbeSstefano_zampini       ierr = PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork);CHKERRQ(ierr);
52625d06dbeSstefano_zampini 
52725d06dbeSstefano_zampini       for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
5281683a169SBarry Smith         const PetscScalar *M;
52925d06dbeSstefano_zampini         PetscInt          subset_size;
53025d06dbeSstefano_zampini 
53125d06dbeSstefano_zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
53225d06dbeSstefano_zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
5331683a169SBarry Smith         ierr = MatDenseGetArrayRead(deluxe_ctx->seq_mat[i],&M);CHKERRQ(ierr);
53425d06dbeSstefano_zampini         ierr = PetscMemcpy(W,M,subset_size*subset_size*sizeof(PetscScalar));CHKERRQ(ierr);
5351683a169SBarry Smith         ierr = MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i],&M);CHKERRQ(ierr);
53625d06dbeSstefano_zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
53725d06dbeSstefano_zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr));
53825d06dbeSstefano_zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
53925d06dbeSstefano_zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
54025d06dbeSstefano_zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
54125d06dbeSstefano_zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
542451a39c7SStefano Zampini         /* silent static analyzer */
543451a39c7SStefano Zampini         if (!idxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"IDXS not present");
54425d06dbeSstefano_zampini         ierr = MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES);CHKERRQ(ierr);
54525d06dbeSstefano_zampini         cum += subset_size;
54625d06dbeSstefano_zampini       }
54725d06dbeSstefano_zampini       ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
54825d06dbeSstefano_zampini       ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
54925d06dbeSstefano_zampini       ierr = MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1);CHKERRQ(ierr);
55025d06dbeSstefano_zampini       ierr = MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2);CHKERRQ(ierr);
55125d06dbeSstefano_zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
55225d06dbeSstefano_zampini       ierr = PetscFree3(W,pivots,Bwork);CHKERRQ(ierr);
55325d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
55425d06dbeSstefano_zampini         ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
55525d06dbeSstefano_zampini       }
55625d06dbeSstefano_zampini     }
557674ae819SStefano Zampini   }
558674ae819SStefano Zampini   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
559674ae819SStefano Zampini   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
560674ae819SStefano Zampini 
561457d4a33Sstefano_zampini   /* Layout of multipliers */
562457d4a33Sstefano_zampini   ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr);
563457d4a33Sstefano_zampini   ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr);
564457d4a33Sstefano_zampini   ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
565457d4a33Sstefano_zampini   ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr);
566457d4a33Sstefano_zampini   ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr);
567457d4a33Sstefano_zampini 
568457d4a33Sstefano_zampini   /* Local work vector of multipliers */
569457d4a33Sstefano_zampini   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
570457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
571457d4a33Sstefano_zampini   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
572457d4a33Sstefano_zampini 
57325d06dbeSstefano_zampini   if (BD2) {
57425d06dbeSstefano_zampini     ISLocalToGlobalMapping l2g;
57525d06dbeSstefano_zampini     Mat                    T,TA,*pT;
57625d06dbeSstefano_zampini     IS                     is;
57725d06dbeSstefano_zampini     PetscInt               nl,N;
57825d06dbeSstefano_zampini     BDdelta_DN             ctx;
57925d06dbeSstefano_zampini 
58025d06dbeSstefano_zampini     ierr = PetscLayoutGetLocalSize(llay,&nl);CHKERRQ(ierr);
58125d06dbeSstefano_zampini     ierr = PetscLayoutGetSize(llay,&N);CHKERRQ(ierr);
58225d06dbeSstefano_zampini     ierr = MatCreate(comm,&T);CHKERRQ(ierr);
58325d06dbeSstefano_zampini     ierr = MatSetSizes(T,nl,nl,N,N);CHKERRQ(ierr);
58425d06dbeSstefano_zampini     ierr = MatSetType(T,MATIS);CHKERRQ(ierr);
58525d06dbeSstefano_zampini     ierr = ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g);CHKERRQ(ierr);
58625d06dbeSstefano_zampini     ierr = MatSetLocalToGlobalMapping(T,l2g,l2g);CHKERRQ(ierr);
58725d06dbeSstefano_zampini     ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr);
58825d06dbeSstefano_zampini     ierr = MatISSetLocalMat(T,BD2);CHKERRQ(ierr);
589487b449aSStefano Zampini     ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
590487b449aSStefano Zampini     ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59125d06dbeSstefano_zampini     ierr = MatDestroy(&BD2);CHKERRQ(ierr);
592487b449aSStefano Zampini     ierr = MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&TA);CHKERRQ(ierr);
59325d06dbeSstefano_zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
59425d06dbeSstefano_zampini     ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
5957dae84e0SHong Zhang     ierr = MatCreateSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT);CHKERRQ(ierr);
59625d06dbeSstefano_zampini     ierr = MatDestroy(&TA);CHKERRQ(ierr);
59725d06dbeSstefano_zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
59825d06dbeSstefano_zampini     BD2  = pT[0];
59925d06dbeSstefano_zampini     ierr = PetscFree(pT);CHKERRQ(ierr);
60025d06dbeSstefano_zampini 
60125d06dbeSstefano_zampini     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
60225d06dbeSstefano_zampini     ierr = PetscNew(&ctx);CHKERRQ(ierr);
60325d06dbeSstefano_zampini     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL);CHKERRQ(ierr);
60425d06dbeSstefano_zampini     ierr = MatShellSetContext(fetidpmat_ctx->B_Ddelta,(void *)ctx);CHKERRQ(ierr);
60525d06dbeSstefano_zampini     ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred);CHKERRQ(ierr);
60625d06dbeSstefano_zampini     ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred);CHKERRQ(ierr);
60725d06dbeSstefano_zampini     ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred);CHKERRQ(ierr);
60825d06dbeSstefano_zampini     ierr = MatSetUp(fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
60925d06dbeSstefano_zampini 
61025d06dbeSstefano_zampini     ierr = PetscObjectReference((PetscObject)BD1);CHKERRQ(ierr);
61125d06dbeSstefano_zampini     ctx->BD = BD1;
61225d06dbeSstefano_zampini     ierr = KSPCreate(PETSC_COMM_SELF,&ctx->kBD);CHKERRQ(ierr);
61325d06dbeSstefano_zampini     ierr = KSPSetOperators(ctx->kBD,BD2,BD2);CHKERRQ(ierr);
61425d06dbeSstefano_zampini     ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work);CHKERRQ(ierr);
61525d06dbeSstefano_zampini     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
61625d06dbeSstefano_zampini   }
61725d06dbeSstefano_zampini   ierr = MatDestroy(&BD1);CHKERRQ(ierr);
61825d06dbeSstefano_zampini   ierr = MatDestroy(&BD2);CHKERRQ(ierr);
61925d06dbeSstefano_zampini 
62025d06dbeSstefano_zampini   /* fetidpmat sizes */
62125d06dbeSstefano_zampini   fetidpmat_ctx->n += nPgl;
62225d06dbeSstefano_zampini   fetidpmat_ctx->N  = fetidpmat_ctx->n_lambda+nPg;
62325d06dbeSstefano_zampini 
624457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
625457d4a33Sstefano_zampini   ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr);
626457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr);
627457d4a33Sstefano_zampini   ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr);
628457d4a33Sstefano_zampini   ierr = VecSetUp(fetidp_global);CHKERRQ(ierr);
629457d4a33Sstefano_zampini 
6309eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
6319eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
632457d4a33Sstefano_zampini      of the FETI-DP linear system */
633457d4a33Sstefano_zampini   if (nPg) {
634e1214c54Sstefano_zampini     Vec            v;
635af140850Sstefano_zampini     IS             IS_l2g_p,ais;
636457d4a33Sstefano_zampini     PetscLayout    alay;
637457d4a33Sstefano_zampini     const PetscInt *idxs,*pranges,*aranges,*lranges;
638af140850Sstefano_zampini     PetscInt       *l2g_indices_p,rst;
639e1214c54Sstefano_zampini     PetscMPIInt    rank;
640457d4a33Sstefano_zampini 
641457d4a33Sstefano_zampini     ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr);
642457d4a33Sstefano_zampini     ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr);
643457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr);
644457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr);
645457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr);
646e1214c54Sstefano_zampini 
647e1214c54Sstefano_zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global),&rank);CHKERRQ(ierr);
648e1214c54Sstefano_zampini     ierr = ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),pranges[rank+1]-pranges[rank],aranges[rank],1,&fetidpmat_ctx->pressure);CHKERRQ(ierr);
64915579a77SStefano Zampini     ierr = PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure,"F_P");CHKERRQ(ierr);
650e1214c54Sstefano_zampini     ierr = ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),lranges[rank+1]-lranges[rank],aranges[rank]+pranges[rank+1]-pranges[rank],1,&fetidpmat_ctx->lagrange);CHKERRQ(ierr);
65115579a77SStefano Zampini     ierr = PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange,"F_L");CHKERRQ(ierr);
652457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
653af140850Sstefano_zampini     /* shift local to global indices for pressure */
654457d4a33Sstefano_zampini     for (i=0;i<nPl;i++) {
655457d4a33Sstefano_zampini       PetscInt owner;
656457d4a33Sstefano_zampini 
657457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr);
658457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
659457d4a33Sstefano_zampini     }
660457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
661457d4a33Sstefano_zampini     ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr);
662af140850Sstefano_zampini 
663e1214c54Sstefano_zampini     /* local to global scatter for pressure */
6649448b7f1SJunchao Zhang     ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr);
665457d4a33Sstefano_zampini     ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr);
666457d4a33Sstefano_zampini 
667e1214c54Sstefano_zampini     /* scatter for lagrange multipliers only */
668e1214c54Sstefano_zampini     ierr = VecCreate(comm,&v);CHKERRQ(ierr);
669e1214c54Sstefano_zampini     ierr = VecSetType(v,VECSTANDARD);CHKERRQ(ierr);
670e1214c54Sstefano_zampini     ierr = VecSetLayout(v,llay);CHKERRQ(ierr);
671e1214c54Sstefano_zampini     ierr = VecSetUp(v);CHKERRQ(ierr);
672e1214c54Sstefano_zampini     ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&ais);CHKERRQ(ierr);
6739448b7f1SJunchao Zhang     ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,v,ais,&fetidpmat_ctx->l2g_lambda_only);CHKERRQ(ierr);
674e1214c54Sstefano_zampini     ierr = ISDestroy(&ais);CHKERRQ(ierr);
675e1214c54Sstefano_zampini     ierr = VecDestroy(&v);CHKERRQ(ierr);
676e1214c54Sstefano_zampini 
677af140850Sstefano_zampini     /* shift local to global indices for multipliers */
678457d4a33Sstefano_zampini     for (i=0;i<n_local_lambda;i++) {
679457d4a33Sstefano_zampini       PetscInt owner,ps;
680457d4a33Sstefano_zampini 
681457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr);
682457d4a33Sstefano_zampini       ps = pranges[owner+1]-pranges[owner];
683457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
684457d4a33Sstefano_zampini     }
685457d4a33Sstefano_zampini 
686e1214c54Sstefano_zampini     /* scatter from alldofs to pressures global fetidp vector */
6879cc7774eSstefano_zampini     ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr);
6889cc7774eSstefano_zampini     ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr);
6899448b7f1SJunchao Zhang     ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr);
6909cc7774eSstefano_zampini     ierr = ISDestroy(&ais);CHKERRQ(ierr);
691457d4a33Sstefano_zampini   }
692457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr);
693457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr);
694457d4a33Sstefano_zampini   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
695a1c0d0daSstefano_zampini 
6969cc7774eSstefano_zampini   /* scatter from local to global multipliers */
6979448b7f1SJunchao Zhang   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
698457d4a33Sstefano_zampini   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
699457d4a33Sstefano_zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr);
700a1c0d0daSstefano_zampini   ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr);
701457d4a33Sstefano_zampini 
702a1c0d0daSstefano_zampini   /* Create some work vectors needed by fetidp */
703674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
704674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
705674ae819SStefano Zampini   PetscFunctionReturn(0);
706674ae819SStefano Zampini }
707674ae819SStefano Zampini 
70825d06dbeSstefano_zampini 
709674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
710674ae819SStefano Zampini {
711674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
71225d06dbeSstefano_zampini   PC_BDDC        *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data;
71325d06dbeSstefano_zampini   PC_IS          *pcis = (PC_IS*)fetidppc_ctx->pc->data;
714f28b6018SStefano Zampini   PetscBool      lumped = PETSC_FALSE;
715674ae819SStefano Zampini   PetscErrorCode ierr;
716674ae819SStefano Zampini 
717674ae819SStefano Zampini   PetscFunctionBegin;
718674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
719674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
720674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
721674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
722674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
723674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
72425d06dbeSstefano_zampini   if (mat_ctx->deluxe_nonred) {
72525d06dbeSstefano_zampini     PC               pc,mpc;
72625d06dbeSstefano_zampini     BDdelta_DN       ctx;
7273ca39a21SBarry Smith     MatSolverType    solver;
72825d06dbeSstefano_zampini     const char       *prefix;
72925d06dbeSstefano_zampini 
73025d06dbeSstefano_zampini     ierr = MatShellGetContext(mat_ctx->B_Ddelta,&ctx);CHKERRQ(ierr);
73125d06dbeSstefano_zampini     ierr = KSPSetType(ctx->kBD,KSPPREONLY);CHKERRQ(ierr);
73225d06dbeSstefano_zampini     ierr = KSPGetPC(ctx->kBD,&mpc);CHKERRQ(ierr);
73325d06dbeSstefano_zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc);CHKERRQ(ierr);
73425d06dbeSstefano_zampini     ierr = PCSetType(mpc,PCLU);CHKERRQ(ierr);
735ea799195SBarry Smith     ierr = PCFactorGetMatSolverType(pc,(MatSolverType*)&solver);CHKERRQ(ierr);
73625d06dbeSstefano_zampini     if (solver) {
7373ca39a21SBarry Smith       ierr = PCFactorSetMatSolverType(mpc,solver);CHKERRQ(ierr);
73825d06dbeSstefano_zampini     }
73925d06dbeSstefano_zampini     ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr);
74025d06dbeSstefano_zampini     ierr = KSPSetOptionsPrefix(ctx->kBD,prefix);CHKERRQ(ierr);
74125d06dbeSstefano_zampini     ierr = KSPAppendOptionsPrefix(ctx->kBD,"bddelta_");CHKERRQ(ierr);
74225d06dbeSstefano_zampini     ierr = KSPSetFromOptions(ctx->kBD);CHKERRQ(ierr);
74325d06dbeSstefano_zampini   }
74425d06dbeSstefano_zampini 
745e1214c54Sstefano_zampini   if (mat_ctx->l2g_lambda_only) {
746e1214c54Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only);CHKERRQ(ierr);
747e1214c54Sstefano_zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
748e1214c54Sstefano_zampini   } else {
749674ae819SStefano Zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
750674ae819SStefano Zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
751e1214c54Sstefano_zampini   }
752f28b6018SStefano Zampini   /* Dirichlet preconditioner */
7539c2d02cdSstefano_zampini   ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL);CHKERRQ(ierr);
754f28b6018SStefano Zampini   if (!lumped) {
7559feb908dSstefano_zampini     IS        iV;
7569c2d02cdSstefano_zampini     PetscBool discrete_harmonic = PETSC_FALSE;
7579c2d02cdSstefano_zampini 
7589feb908dSstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iV",(PetscObject*)&iV);CHKERRQ(ierr);
7599feb908dSstefano_zampini     if (iV) {
7609c2d02cdSstefano_zampini       ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL);CHKERRQ(ierr);
7619c2d02cdSstefano_zampini     }
7629c2d02cdSstefano_zampini     if (discrete_harmonic) {
7639c2d02cdSstefano_zampini       KSP             sksp;
7649c2d02cdSstefano_zampini       PC              pc;
765111315fdSstefano_zampini       PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
7669c2d02cdSstefano_zampini       Mat             A_II,A_IB,A_BI;
767111315fdSstefano_zampini       IS              iP = NULL;
768111315fdSstefano_zampini       PetscBool       isshell,reuse = PETSC_FALSE;
7699c2d02cdSstefano_zampini       KSPType         ksptype;
7709c2d02cdSstefano_zampini       const char      *prefix;
7719c2d02cdSstefano_zampini 
7729c2d02cdSstefano_zampini       /*
7739c2d02cdSstefano_zampini         We constructs a Schur complement for
7749c2d02cdSstefano_zampini 
7759c2d02cdSstefano_zampini         | A_II A_ID |
7769c2d02cdSstefano_zampini         | A_DI A_DD |
7779c2d02cdSstefano_zampini 
7789c2d02cdSstefano_zampini         instead of
7799c2d02cdSstefano_zampini 
7809c2d02cdSstefano_zampini         | A_II  B^t_II A_ID |
7819c2d02cdSstefano_zampini         | B_II -C_II   B_ID |
7829c2d02cdSstefano_zampini         | A_DI  B^t_ID A_DD |
7839c2d02cdSstefano_zampini 
7849c2d02cdSstefano_zampini       */
785111315fdSstefano_zampini       if (sub_schurs && sub_schurs->reuse_solver) {
786111315fdSstefano_zampini         ierr = PetscObjectQuery((PetscObject)sub_schurs->A,"__KSPFETIDP_iP",(PetscObject*)&iP);CHKERRQ(ierr);
787111315fdSstefano_zampini         if (iP) reuse = PETSC_TRUE;
788111315fdSstefano_zampini       }
789111315fdSstefano_zampini       if (!reuse) {
790111315fdSstefano_zampini         IS       aB;
791111315fdSstefano_zampini         PetscInt nb;
7929c2d02cdSstefano_zampini         ierr = ISGetLocalSize(pcis->is_B_local,&nb);CHKERRQ(ierr);
7939c2d02cdSstefano_zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB);CHKERRQ(ierr);
7949feb908dSstefano_zampini         ierr = MatCreateSubMatrix(pcis->A_II,iV,iV,MAT_INITIAL_MATRIX,&A_II);CHKERRQ(ierr);
7959feb908dSstefano_zampini         ierr = MatCreateSubMatrix(pcis->A_IB,iV,aB,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
7969feb908dSstefano_zampini         ierr = MatCreateSubMatrix(pcis->A_BI,aB,iV,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
797111315fdSstefano_zampini         ierr = ISDestroy(&aB);CHKERRQ(ierr);
798111315fdSstefano_zampini       } else {
799111315fdSstefano_zampini         ierr = MatCreateSubMatrix(sub_schurs->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
800111315fdSstefano_zampini         ierr = MatCreateSubMatrix(sub_schurs->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
801111315fdSstefano_zampini         ierr = PetscObjectReference((PetscObject)pcis->A_II);CHKERRQ(ierr);
802111315fdSstefano_zampini         A_II = pcis->A_II;
803111315fdSstefano_zampini       }
8049c2d02cdSstefano_zampini       ierr = MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
8059c2d02cdSstefano_zampini 
8069c2d02cdSstefano_zampini       /* propagate settings of solver */
8079c2d02cdSstefano_zampini       ierr = MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp);CHKERRQ(ierr);
8089c2d02cdSstefano_zampini       ierr = KSPGetType(pcis->ksp_D,&ksptype);CHKERRQ(ierr);
8099c2d02cdSstefano_zampini       ierr = KSPSetType(sksp,ksptype);CHKERRQ(ierr);
8109c2d02cdSstefano_zampini       ierr = KSPGetPC(pcis->ksp_D,&pc);CHKERRQ(ierr);
8115334bea6Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell);CHKERRQ(ierr);
8125334bea6Sstefano_zampini       if (!isshell) {
8133ca39a21SBarry Smith         MatSolverType    solver;
8145334bea6Sstefano_zampini         PCType           pctype;
8155334bea6Sstefano_zampini 
8169c2d02cdSstefano_zampini         ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
817ea799195SBarry Smith         ierr = PCFactorGetMatSolverType(pc,(MatSolverType*)&solver);CHKERRQ(ierr);
8189c2d02cdSstefano_zampini         ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr);
8199c2d02cdSstefano_zampini         ierr = PCSetType(pc,pctype);CHKERRQ(ierr);
8209c2d02cdSstefano_zampini         if (solver) {
8213ca39a21SBarry Smith           ierr = PCFactorSetMatSolverType(pc,solver);CHKERRQ(ierr);
8229c2d02cdSstefano_zampini         }
8235334bea6Sstefano_zampini       } else {
8245334bea6Sstefano_zampini         ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr);
8255334bea6Sstefano_zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
8265334bea6Sstefano_zampini       }
8279c2d02cdSstefano_zampini       ierr = MatDestroy(&A_II);CHKERRQ(ierr);
8289c2d02cdSstefano_zampini       ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
8299c2d02cdSstefano_zampini       ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
8309c2d02cdSstefano_zampini       ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr);
8319c2d02cdSstefano_zampini       ierr = KSPSetOptionsPrefix(sksp,prefix);CHKERRQ(ierr);
8329c2d02cdSstefano_zampini       ierr = KSPAppendOptionsPrefix(sksp,"harmonic_");CHKERRQ(ierr);
8333016320fSstefano_zampini       ierr = KSPSetFromOptions(sksp);CHKERRQ(ierr);
834111315fdSstefano_zampini       if (reuse) {
835111315fdSstefano_zampini         ierr = KSPSetPC(sksp,sub_schurs->reuse_solver->interior_solver);CHKERRQ(ierr);
83615579a77SStefano Zampini         ierr = PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver,(PetscObject)sksp,0);CHKERRQ(ierr);
837111315fdSstefano_zampini       }
8389c2d02cdSstefano_zampini     } else { /* default Dirichlet preconditioner is pde-harmonic */
839ed6c3d69SStefano Zampini       ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
840ed6c3d69SStefano Zampini       ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
8419c2d02cdSstefano_zampini     }
842f28b6018SStefano Zampini   } else {
843f28b6018SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr);
844f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
845f28b6018SStefano Zampini   }
846af140850Sstefano_zampini   /* saddle-point */
847af140850Sstefano_zampini   if (mat_ctx->xPg) {
848af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr);
849af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
850af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr);
851af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
852af140850Sstefano_zampini   }
853674ae819SStefano Zampini   PetscFunctionReturn(0);
854674ae819SStefano Zampini }
855674ae819SStefano Zampini 
856ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
857674ae819SStefano Zampini {
858674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
859617d11aeSStefano Zampini   PC_BDDC        *pcbddc;
860674ae819SStefano Zampini   PC_IS          *pcis;
861674ae819SStefano Zampini   PetscErrorCode ierr;
862674ae819SStefano Zampini 
863674ae819SStefano Zampini   PetscFunctionBegin;
864674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
865674ae819SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
866617d11aeSStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
867674ae819SStefano Zampini   /* Application of B_delta^T */
868af140850Sstefano_zampini   ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr);
869674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
870674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
871674ae819SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
872af140850Sstefano_zampini 
873af140850Sstefano_zampini   /* Add contribution from saddle point */
874af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
875af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
876af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
877af140850Sstefano_zampini     if (pcbddc->switch_static) {
878ef425aeaSstefano_zampini       if (trans) {
879ef425aeaSstefano_zampini         ierr = MatMultTranspose(mat_ctx->B_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr);
880ef425aeaSstefano_zampini       } else {
881af140850Sstefano_zampini         ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr);
882af140850Sstefano_zampini       }
883ef425aeaSstefano_zampini     }
884ef425aeaSstefano_zampini     if (trans) {
885ef425aeaSstefano_zampini       ierr = MatMultTransposeAdd(mat_ctx->B_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
886ef425aeaSstefano_zampini     } else {
887af140850Sstefano_zampini       ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
888ef425aeaSstefano_zampini     }
889af140850Sstefano_zampini   } else {
890af140850Sstefano_zampini     if (pcbddc->switch_static) {
891674ae819SStefano Zampini       ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
892af140850Sstefano_zampini     }
893af140850Sstefano_zampini   }
894af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
895617d11aeSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
896ef425aeaSstefano_zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,trans);CHKERRQ(ierr);
897c7ffc8ceSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
898af140850Sstefano_zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
899674ae819SStefano Zampini   /* Application of B_delta */
900674ae819SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
901af140850Sstefano_zampini   /* Contribution from boundary pressures */
902af140850Sstefano_zampini   if (mat_ctx->C) {
903af140850Sstefano_zampini     const PetscScalar *lx;
904af140850Sstefano_zampini     PetscScalar       *ly;
905af140850Sstefano_zampini 
90615579a77SStefano Zampini     /* pressure ordered first in the local part of x and y */
907af140850Sstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
908af140850Sstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
909af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr);
910af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr);
911ef425aeaSstefano_zampini     if (trans) {
912ef425aeaSstefano_zampini       ierr = MatMultTranspose(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr);
913ef425aeaSstefano_zampini     } else {
914af140850Sstefano_zampini       ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr);
915ef425aeaSstefano_zampini     }
916af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr);
917af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr);
918af140850Sstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
919af140850Sstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
920af140850Sstefano_zampini   }
921af140850Sstefano_zampini   /* Add contribution from saddle point */
922af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
923ef425aeaSstefano_zampini     if (trans) {
924ef425aeaSstefano_zampini       ierr = MatMultTranspose(mat_ctx->Bt_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
925ef425aeaSstefano_zampini     } else {
926af140850Sstefano_zampini       ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
927ef425aeaSstefano_zampini     }
928af140850Sstefano_zampini     if (pcbddc->switch_static) {
929ef425aeaSstefano_zampini       if (trans) {
930ef425aeaSstefano_zampini         ierr = MatMultTransposeAdd(mat_ctx->Bt_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
931ef425aeaSstefano_zampini       } else {
932af140850Sstefano_zampini         ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
933af140850Sstefano_zampini       }
934ef425aeaSstefano_zampini     }
935af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
936af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
937af140850Sstefano_zampini   }
938674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
939674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
940674ae819SStefano Zampini   PetscFunctionReturn(0);
941674ae819SStefano Zampini }
942674ae819SStefano Zampini 
943ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
944edf7251bSStefano Zampini {
945edf7251bSStefano Zampini   PetscErrorCode ierr;
946edf7251bSStefano Zampini 
947edf7251bSStefano Zampini   PetscFunctionBegin;
948ef425aeaSstefano_zampini   ierr = FETIDPMatMult_Kernel(fetimat,x,y,PETSC_FALSE);CHKERRQ(ierr);
949ef425aeaSstefano_zampini   PetscFunctionReturn(0);
950ef425aeaSstefano_zampini }
951ef425aeaSstefano_zampini 
952ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
953ef425aeaSstefano_zampini {
954ef425aeaSstefano_zampini   PetscErrorCode ierr;
955ef425aeaSstefano_zampini 
956ef425aeaSstefano_zampini   PetscFunctionBegin;
957ef425aeaSstefano_zampini   ierr = FETIDPMatMult_Kernel(fetimat,x,y,PETSC_TRUE);CHKERRQ(ierr);
958edf7251bSStefano Zampini   PetscFunctionReturn(0);
959edf7251bSStefano Zampini }
960edf7251bSStefano Zampini 
96116597391Sstefano_zampini PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
962674ae819SStefano Zampini {
963674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
964674ae819SStefano Zampini   PC_IS          *pcis;
965674ae819SStefano Zampini   PetscErrorCode ierr;
966674ae819SStefano Zampini 
967674ae819SStefano Zampini   PetscFunctionBegin;
968302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
969674ae819SStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
970674ae819SStefano Zampini   /* Application of B_Ddelta^T */
971674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
972674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
973674ae819SStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
974674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
975ed6c3d69SStefano Zampini   /* Application of local Schur complement */
97616597391Sstefano_zampini   if (trans) {
97716597391Sstefano_zampini     ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
97816597391Sstefano_zampini   } else {
979ed6c3d69SStefano Zampini     ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
98016597391Sstefano_zampini   }
981edf7251bSStefano Zampini   /* Application of B_Ddelta */
982edf7251bSStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
983edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
984edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
985edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
986edf7251bSStefano Zampini   PetscFunctionReturn(0);
987edf7251bSStefano Zampini }
988edf7251bSStefano Zampini 
98916597391Sstefano_zampini PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
990edf7251bSStefano Zampini {
991edf7251bSStefano Zampini   PetscErrorCode ierr;
992edf7251bSStefano Zampini 
993edf7251bSStefano Zampini   PetscFunctionBegin;
99416597391Sstefano_zampini   ierr = FETIDPPCApply_Kernel(pc,x,y,PETSC_FALSE);CHKERRQ(ierr);
99516597391Sstefano_zampini   PetscFunctionReturn(0);
99616597391Sstefano_zampini }
99716597391Sstefano_zampini 
99816597391Sstefano_zampini PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
99916597391Sstefano_zampini {
100016597391Sstefano_zampini   PetscErrorCode ierr;
100116597391Sstefano_zampini 
100216597391Sstefano_zampini   PetscFunctionBegin;
100316597391Sstefano_zampini   ierr = FETIDPPCApply_Kernel(pc,x,y,PETSC_TRUE);CHKERRQ(ierr);
1004674ae819SStefano Zampini   PetscFunctionReturn(0);
1005674ae819SStefano Zampini }
1006c45b8d2dSstefano_zampini 
1007c45b8d2dSstefano_zampini PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
1008c45b8d2dSstefano_zampini {
1009c45b8d2dSstefano_zampini   FETIDPPC_ctx      pc_ctx;
1010c45b8d2dSstefano_zampini   PetscBool         iascii;
1011c45b8d2dSstefano_zampini   PetscViewer       sviewer;
1012c45b8d2dSstefano_zampini   PetscErrorCode    ierr;
1013c45b8d2dSstefano_zampini 
1014c45b8d2dSstefano_zampini   PetscFunctionBegin;
1015c45b8d2dSstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1016c45b8d2dSstefano_zampini   if (iascii) {
1017c45b8d2dSstefano_zampini     PetscMPIInt rank;
101825d06dbeSstefano_zampini     PetscBool   isschur,isshell;
1019c45b8d2dSstefano_zampini 
1020c45b8d2dSstefano_zampini     ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
1021c45b8d2dSstefano_zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
1022c45b8d2dSstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur);CHKERRQ(ierr);
1023c45b8d2dSstefano_zampini     if (isschur) {
102415579a77SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"  Dirichlet preconditioner (just from rank 0)\n");CHKERRQ(ierr);
1025c45b8d2dSstefano_zampini     } else {
102615579a77SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"  Lumped preconditioner (just from rank 0)\n");CHKERRQ(ierr);
1027c45b8d2dSstefano_zampini     }
1028c45b8d2dSstefano_zampini     ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
1029c45b8d2dSstefano_zampini     if (!rank) {
1030c45b8d2dSstefano_zampini       ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
103115579a77SStefano Zampini       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
1032c45b8d2dSstefano_zampini       ierr = MatView(pc_ctx->S_j,sviewer);CHKERRQ(ierr);
103315579a77SStefano Zampini       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
1034c45b8d2dSstefano_zampini       ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr);
1035c45b8d2dSstefano_zampini     }
1036*e5afcf28SBarry Smith     ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
103725d06dbeSstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell);CHKERRQ(ierr);
103825d06dbeSstefano_zampini     if (isshell) {
103925d06dbeSstefano_zampini       BDdelta_DN ctx;
104025d06dbeSstefano_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);
104125d06dbeSstefano_zampini       ierr = MatShellGetContext(pc_ctx->B_Ddelta,&ctx);CHKERRQ(ierr);
1042*e5afcf28SBarry Smith       ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
104325d06dbeSstefano_zampini       if (!rank) {
104415579a77SStefano Zampini         PetscInt tl;
104515579a77SStefano Zampini 
104615579a77SStefano Zampini         ierr = PetscViewerASCIIGetTab(sviewer,&tl);CHKERRQ(ierr);
104715579a77SStefano Zampini         ierr = PetscObjectSetTabLevel((PetscObject)ctx->kBD,tl);CHKERRQ(ierr);
104825d06dbeSstefano_zampini         ierr = KSPView(ctx->kBD,sviewer);CHKERRQ(ierr);
104925d06dbeSstefano_zampini         ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
105025d06dbeSstefano_zampini         ierr = MatView(ctx->BD,sviewer);CHKERRQ(ierr);
105125d06dbeSstefano_zampini         ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr);
105225d06dbeSstefano_zampini       }
1053c45b8d2dSstefano_zampini       ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr);
1054*e5afcf28SBarry Smith     }
1055c45b8d2dSstefano_zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1056c45b8d2dSstefano_zampini   }
1057c45b8d2dSstefano_zampini   PetscFunctionReturn(0);
1058c45b8d2dSstefano_zampini }
1059