xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 2d456bbaf1b6b248fa65016de2a9c7278a2f8390)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3674ae819SStefano Zampini 
4674ae819SStefano Zampini #undef __FUNCT__
5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPMatContext"
6674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
7674ae819SStefano Zampini {
8674ae819SStefano Zampini   FETIDPMat_ctx  newctx;
9674ae819SStefano Zampini   PetscErrorCode ierr;
10674ae819SStefano Zampini 
11674ae819SStefano Zampini   PetscFunctionBegin;
12854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
13674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
14674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
15674ae819SStefano Zampini   newctx->pc              = pc;
16674ae819SStefano Zampini   *fetidpmat_ctx          = newctx;
17674ae819SStefano Zampini   PetscFunctionReturn(0);
18674ae819SStefano Zampini }
19674ae819SStefano Zampini 
20674ae819SStefano Zampini #undef __FUNCT__
21674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext"
22674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
23674ae819SStefano Zampini {
24674ae819SStefano Zampini   FETIDPPC_ctx   newctx;
25674ae819SStefano Zampini   PetscErrorCode ierr;
26674ae819SStefano Zampini 
27674ae819SStefano Zampini   PetscFunctionBegin;
28854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
29674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
30674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
31674ae819SStefano Zampini   newctx->pc              = pc;
32674ae819SStefano Zampini   *fetidppc_ctx           = newctx;
33674ae819SStefano Zampini   PetscFunctionReturn(0);
34674ae819SStefano Zampini }
35674ae819SStefano Zampini 
36674ae819SStefano Zampini #undef __FUNCT__
37674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat"
38674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
39674ae819SStefano Zampini {
40674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
41674ae819SStefano Zampini   PetscErrorCode ierr;
42674ae819SStefano Zampini 
43674ae819SStefano Zampini   PetscFunctionBegin;
44674ae819SStefano Zampini   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
45674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
46674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
47674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
48674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
49674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
50457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr);
51457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr);
52457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr);
53457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr);
54457d4a33Sstefano_zampini   ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr);
55457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr);
56457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr);
57457d4a33Sstefano_zampini   ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr);
58674ae819SStefano Zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
59457d4a33Sstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr);
609cc7774eSstefano_zampini   ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr);
61674ae819SStefano Zampini   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
62674ae819SStefano Zampini   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
63674ae819SStefano Zampini   PetscFunctionReturn(0);
64674ae819SStefano Zampini }
65674ae819SStefano Zampini 
66674ae819SStefano Zampini #undef __FUNCT__
67674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC"
68674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
69674ae819SStefano Zampini {
70674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
71674ae819SStefano Zampini   PetscErrorCode ierr;
72674ae819SStefano Zampini 
73674ae819SStefano Zampini   PetscFunctionBegin;
74674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
75674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
76674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
77674ae819SStefano Zampini   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
7868668abeSStefano Zampini   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
79674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
80457d4a33Sstefano_zampini   ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr);
81457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr);
82457d4a33Sstefano_zampini   ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr);
83674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
84674ae819SStefano Zampini   PetscFunctionReturn(0);
85674ae819SStefano Zampini }
86674ae819SStefano Zampini 
87674ae819SStefano Zampini #undef __FUNCT__
88674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext"
89674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
90674ae819SStefano Zampini {
91674ae819SStefano Zampini   PetscErrorCode ierr;
92674ae819SStefano Zampini   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
93674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
94674ae819SStefano Zampini   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
95674ae819SStefano Zampini   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
96674ae819SStefano Zampini   MPI_Comm       comm;
97674ae819SStefano Zampini   Mat            ScalingMat;
98457d4a33Sstefano_zampini   Vec            fetidp_global;
99674ae819SStefano Zampini   IS             IS_l2g_lambda;
100dc456d91SStefano Zampini   IS             subset,subset_mult,subset_n;
101674ae819SStefano Zampini   PetscBool      skip_node,fully_redundant;
102674ae819SStefano Zampini   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
103dc456d91SStefano Zampini   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
10476ec1555SBarry Smith   PetscMPIInt    rank,size,buf_size,neigh;
105674ae819SStefano Zampini   PetscScalar    scalar_value;
106674ae819SStefano Zampini   PetscInt       *vertex_indices;
107dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
108dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
109674ae819SStefano Zampini   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
110674ae819SStefano Zampini   PetscScalar    *array,*scaling_factors,*vals_B_delta;
111674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
112457d4a33Sstefano_zampini   PetscLayout    llay;
113674ae819SStefano Zampini   /* For communication of scaling factors */
114674ae819SStefano Zampini   PetscInt       *ptrs_buffer,neigh_position;
115674ae819SStefano Zampini   PetscScalar    **all_factors,*send_buffer,*recv_buffer;
116674ae819SStefano Zampini   MPI_Request    *send_reqs,*recv_reqs;
117674ae819SStefano Zampini   /* tests */
118674ae819SStefano Zampini   PetscBool      test_fetidp;
119674ae819SStefano Zampini   PetscViewer    viewer;
120457d4a33Sstefano_zampini   /* saddlepoint */
121457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
122457d4a33Sstefano_zampini   PetscLayout            play;
123457d4a33Sstefano_zampini   IS                     gP,pP;
124457d4a33Sstefano_zampini   PetscInt               nPl,nPg,nPgl;
125674ae819SStefano Zampini 
126674ae819SStefano Zampini   PetscFunctionBegin;
127674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
128674ae819SStefano Zampini   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12976ec1555SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
130674ae819SStefano Zampini 
131457d4a33Sstefano_zampini   /* saddlepoint */
132457d4a33Sstefano_zampini   nPl      = 0;
133457d4a33Sstefano_zampini   nPg      = 0;
134457d4a33Sstefano_zampini   nPgl     = 0;
135457d4a33Sstefano_zampini   gP       = NULL;
136457d4a33Sstefano_zampini   pP       = NULL;
137457d4a33Sstefano_zampini   l2gmap_p = NULL;
138457d4a33Sstefano_zampini   play     = NULL;
139457d4a33Sstefano_zampini   ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr);
140022d8d2bSstefano_zampini   if (pP) { /* saddle point */
141457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
142457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr);
143457d4a33Sstefano_zampini     if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
144457d4a33Sstefano_zampini     ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr);
145457d4a33Sstefano_zampini     ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr);
146457d4a33Sstefano_zampini     ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr);
147457d4a33Sstefano_zampini     ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr);
148457d4a33Sstefano_zampini     ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr);
149457d4a33Sstefano_zampini 
150457d4a33Sstefano_zampini     /* interface pressure matrix */
151457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr);
152457d4a33Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */
153457d4a33Sstefano_zampini       IS Pg;
154457d4a33Sstefano_zampini 
155457d4a33Sstefano_zampini       ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr);
156457d4a33Sstefano_zampini       ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr);
157457d4a33Sstefano_zampini       ierr = ISDestroy(&Pg);CHKERRQ(ierr);
158457d4a33Sstefano_zampini       ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr);
159457d4a33Sstefano_zampini       ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr);
160457d4a33Sstefano_zampini       ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr);
161457d4a33Sstefano_zampini       ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr);
162457d4a33Sstefano_zampini       ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr);
163457d4a33Sstefano_zampini       ierr = PetscLayoutSetUp(play);CHKERRQ(ierr);
164457d4a33Sstefano_zampini     } else {
165457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr);
166457d4a33Sstefano_zampini       ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr);
167457d4a33Sstefano_zampini       ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr);
168457d4a33Sstefano_zampini       ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr);
169457d4a33Sstefano_zampini       ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr);
170457d4a33Sstefano_zampini       ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr);
171457d4a33Sstefano_zampini       ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr);
172457d4a33Sstefano_zampini     }
173457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr);
174457d4a33Sstefano_zampini     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr);
175457d4a33Sstefano_zampini 
176457d4a33Sstefano_zampini     /* import matrices for interface pressures coupling */
177457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr);
178457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
179457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr);
180457d4a33Sstefano_zampini 
181457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr);
182457d4a33Sstefano_zampini     if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
183457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr);
184457d4a33Sstefano_zampini 
185457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
186457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
187457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
188457d4a33Sstefano_zampini 
189457d4a33Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
190457d4a33Sstefano_zampini     if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
191457d4a33Sstefano_zampini     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
192457d4a33Sstefano_zampini   }
193457d4a33Sstefano_zampini 
194674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
195329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
196674ae819SStefano Zampini 
197674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
198674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
199674ae819SStefano Zampini   n_local_lambda = 0;
200674ae819SStefano Zampini   partial_sum = 0;
201674ae819SStefano Zampini   n_boundary_dofs = 0;
202674ae819SStefano Zampini   s = 0;
203674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
204b371cd4fSStefano Zampini   n_vertices = pcbddc->n_vertices;
2050e6343abSStefano Zampini   vertex_indices = pcbddc->local_primal_ref_node;
206674ae819SStefano Zampini   dual_size = pcis->n_B-n_vertices;
207785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
208785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
209785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
210674ae819SStefano Zampini 
211674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
212674ae819SStefano Zampini   for (i=0;i<pcis->n;i++){
213674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
214*2d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
215674ae819SStefano Zampini     skip_node = PETSC_FALSE;
216674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
217674ae819SStefano Zampini       skip_node = PETSC_TRUE;
218674ae819SStefano Zampini       s++;
219674ae819SStefano Zampini     }
220*2d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
221*2d456bbaSstefano_zampini     if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
222674ae819SStefano Zampini     if (!skip_node) {
223674ae819SStefano Zampini       if (fully_redundant) {
224674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
225674ae819SStefano Zampini         n_lambda_for_dof = (j*(j+1))/2;
226674ae819SStefano Zampini       } else {
227674ae819SStefano Zampini         n_lambda_for_dof = j;
228674ae819SStefano Zampini       }
229674ae819SStefano Zampini       n_local_lambda += j;
230674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
231674ae819SStefano Zampini       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
232674ae819SStefano Zampini       /* store some data needed */
233674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
234674ae819SStefano Zampini       aux_local_numbering_1[partial_sum] = i;
235674ae819SStefano Zampini       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
236674ae819SStefano Zampini       partial_sum++;
237674ae819SStefano Zampini     }
238674ae819SStefano Zampini   }
239674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
240*2d456bbaSstefano_zampini   dual_size = partial_sum;
241674ae819SStefano Zampini 
242674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
243dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
2443bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
245dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
246dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
2473d996552SStefano Zampini   ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr);
248dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
2493d996552SStefano Zampini 
2504fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG)
2514fe826edSStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
2524fe826edSStefano Zampini   ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2534fe826edSStefano Zampini   ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2544fe826edSStefano Zampini   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
2554fe826edSStefano Zampini   i = (PetscInt)PetscRealPart(scalar_value);
2566c4ed002SBarry 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);
2574fe826edSStefano Zampini #endif
258674ae819SStefano Zampini 
259674ae819SStefano Zampini   /* init data for scaling factors exchange */
260674ae819SStefano Zampini   partial_sum = 0;
261785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
2624b2aedd3SStefano Zampini   ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr);
2634b2aedd3SStefano Zampini   ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr);
264785e854fSJed Brown   ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr);
2654b2aedd3SStefano Zampini   if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
266674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
267674ae819SStefano Zampini     partial_sum += pcis->n_shared[i];
268674ae819SStefano Zampini     ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
269674ae819SStefano Zampini   }
270785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
271785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
272785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
273674ae819SStefano Zampini   for (i=0;i<pcis->n-1;i++) {
274674ae819SStefano Zampini     j = mat_graph->count[i];
275674ae819SStefano Zampini     all_factors[i+1]=all_factors[i]+j;
276674ae819SStefano Zampini   }
277674ae819SStefano Zampini   /* scatter B scaling to N vec */
278674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
279674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
280674ae819SStefano Zampini   /* communications */
2812b095fd8SStefano Zampini   ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
282674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
283674ae819SStefano Zampini     for (j=0;j<pcis->n_shared[i];j++) {
284674ae819SStefano Zampini       send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
285674ae819SStefano Zampini     }
286674ae819SStefano Zampini     ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
287674ae819SStefano Zampini     ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
288674ae819SStefano Zampini     ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
289674ae819SStefano Zampini     ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
290674ae819SStefano Zampini   }
2912b095fd8SStefano Zampini   ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
2924b2aedd3SStefano Zampini   if (pcis->n_neigh > 0) {
2934b2aedd3SStefano Zampini     ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
2944b2aedd3SStefano Zampini   }
295674ae819SStefano Zampini   /* put values in correct places */
296674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
297674ae819SStefano Zampini     for (j=0;j<pcis->n_shared[i];j++) {
298674ae819SStefano Zampini       k = pcis->shared[i][j];
299674ae819SStefano Zampini       neigh_position = 0;
300674ae819SStefano Zampini       while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
301674ae819SStefano Zampini       all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
302674ae819SStefano Zampini     }
303674ae819SStefano Zampini   }
3044b2aedd3SStefano Zampini   if (pcis->n_neigh > 0) {
3054b2aedd3SStefano Zampini     ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
3064b2aedd3SStefano Zampini   }
307674ae819SStefano Zampini   ierr = PetscFree(send_reqs);CHKERRQ(ierr);
308674ae819SStefano Zampini   ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
309674ae819SStefano Zampini   ierr = PetscFree(send_buffer);CHKERRQ(ierr);
310674ae819SStefano Zampini   ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
311674ae819SStefano Zampini   ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
312674ae819SStefano Zampini 
313674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
314785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
315785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
316785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
317785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
318785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
319dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
320674ae819SStefano Zampini   partial_sum=0;
321dc456d91SStefano Zampini   cum = 0;
322674ae819SStefano Zampini   for (i=0;i<dual_size;i++) {
323dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
324674ae819SStefano Zampini     j = mat_graph->count[aux_local_numbering_1[i]];
325674ae819SStefano Zampini     aux_sums[0]=0;
326674ae819SStefano Zampini     for (s=1;s<j;s++) {
327674ae819SStefano Zampini       aux_sums[s]=aux_sums[s-1]+j-s+1;
328674ae819SStefano Zampini     }
329674ae819SStefano Zampini     array = all_factors[aux_local_numbering_1[i]];
330674ae819SStefano Zampini     n_neg_values = 0;
3312a7da448SStefano Zampini     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
332674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
333674ae819SStefano Zampini     if (fully_redundant) {
334674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
335674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
336674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
337674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=-1.0;
338674ae819SStefano Zampini         scaling_factors[partial_sum+s]=array[s];
339674ae819SStefano Zampini       }
340674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
341674ae819SStefano Zampini         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
342674ae819SStefano Zampini         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
343674ae819SStefano Zampini         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
344674ae819SStefano Zampini         scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
345674ae819SStefano Zampini       }
346674ae819SStefano Zampini       partial_sum += j;
347674ae819SStefano Zampini     } else {
348674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
349674ae819SStefano Zampini       for (s=0;s<j;s++) {
350674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=n_global_lambda+s;
351674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
352674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=0.0;
353674ae819SStefano Zampini       }
354674ae819SStefano Zampini       /* B_delta */
355674ae819SStefano Zampini       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
356674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
357674ae819SStefano Zampini       }
358674ae819SStefano Zampini       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
359674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values]=1.0;
360674ae819SStefano Zampini       }
361674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999*/
362674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
363674ae819SStefano Zampini         scalar_value = 0.0;
364674ae819SStefano Zampini         for (k=0;k<s+1;k++) {
365674ae819SStefano Zampini           scalar_value += array[k];
366674ae819SStefano Zampini         }
367674ae819SStefano Zampini         scaling_factors[partial_sum+s] = -scalar_value;
368674ae819SStefano Zampini       }
369674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
370674ae819SStefano Zampini         scalar_value = 0.0;
371674ae819SStefano Zampini         for (k=s+n_neg_values;k<j;k++) {
372674ae819SStefano Zampini           scalar_value += array[k];
373674ae819SStefano Zampini         }
374674ae819SStefano Zampini         scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
375674ae819SStefano Zampini       }
376674ae819SStefano Zampini       partial_sum += j;
377674ae819SStefano Zampini     }
378dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
379674ae819SStefano Zampini   }
380dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
381dc456d91SStefano Zampini   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
382dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
383674ae819SStefano Zampini   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
384674ae819SStefano Zampini   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
385674ae819SStefano Zampini   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
386674ae819SStefano Zampini   ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
387674ae819SStefano Zampini   ierr = PetscFree(all_factors);CHKERRQ(ierr);
388674ae819SStefano Zampini 
389674ae819SStefano Zampini   /* Create local part of B_delta */
390302440fdSBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
391674ae819SStefano Zampini   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
392674ae819SStefano Zampini   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
393674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
394674ae819SStefano Zampini   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
395674ae819SStefano Zampini   for (i=0;i<n_local_lambda;i++) {
396674ae819SStefano Zampini     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
397674ae819SStefano Zampini   }
398674ae819SStefano Zampini   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
399674ae819SStefano Zampini   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
400674ae819SStefano Zampini   ierr = MatAssemblyEnd  (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
401674ae819SStefano Zampini 
402674ae819SStefano Zampini   if (fully_redundant) {
403302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
404674ae819SStefano Zampini     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
405674ae819SStefano Zampini     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
406674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
407674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
408674ae819SStefano Zampini       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
409674ae819SStefano Zampini     }
410674ae819SStefano Zampini     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411674ae819SStefano Zampini     ierr = MatAssemblyEnd  (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
412674ae819SStefano Zampini     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
413674ae819SStefano Zampini     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
414674ae819SStefano Zampini   } else {
415302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
416674ae819SStefano Zampini     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
417674ae819SStefano Zampini     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
418674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
419674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
420674ae819SStefano Zampini       ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
421674ae819SStefano Zampini     }
422674ae819SStefano Zampini     ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
423674ae819SStefano Zampini     ierr = MatAssemblyEnd  (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
424674ae819SStefano Zampini   }
425674ae819SStefano Zampini   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
426674ae819SStefano Zampini   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
427674ae819SStefano Zampini 
428457d4a33Sstefano_zampini   /* Layout of multipliers */
429457d4a33Sstefano_zampini   ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr);
430457d4a33Sstefano_zampini   ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr);
431457d4a33Sstefano_zampini   ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
432457d4a33Sstefano_zampini   ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr);
433457d4a33Sstefano_zampini   ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr);
434457d4a33Sstefano_zampini 
435457d4a33Sstefano_zampini   /* fetidpmat sizes */
436457d4a33Sstefano_zampini   fetidpmat_ctx->n += nPgl;
437457d4a33Sstefano_zampini   fetidpmat_ctx->N  = fetidpmat_ctx->n_lambda+nPg;
438457d4a33Sstefano_zampini 
439457d4a33Sstefano_zampini   /* Local work vector of multipliers */
440457d4a33Sstefano_zampini   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
441457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
442457d4a33Sstefano_zampini   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
443457d4a33Sstefano_zampini 
444457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
445457d4a33Sstefano_zampini   ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr);
446457d4a33Sstefano_zampini   ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr);
447457d4a33Sstefano_zampini   ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr);
448457d4a33Sstefano_zampini   ierr = VecSetUp(fetidp_global);CHKERRQ(ierr);
449457d4a33Sstefano_zampini 
4509eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
4519eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
452457d4a33Sstefano_zampini      of the FETI-DP linear system */
453457d4a33Sstefano_zampini   if (nPg) {
454af140850Sstefano_zampini     IS             IS_l2g_p,ais;
455457d4a33Sstefano_zampini     PetscLayout    alay;
456457d4a33Sstefano_zampini     const PetscInt *idxs,*pranges,*aranges,*lranges;
457af140850Sstefano_zampini     PetscInt       *l2g_indices_p,rst;
458457d4a33Sstefano_zampini 
459457d4a33Sstefano_zampini     ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr);
460457d4a33Sstefano_zampini     ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr);
461457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr);
462457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr);
463457d4a33Sstefano_zampini     ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr);
464457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
465af140850Sstefano_zampini     /* shift local to global indices for pressure */
466457d4a33Sstefano_zampini     for (i=0;i<nPl;i++) {
467457d4a33Sstefano_zampini       PetscInt owner;
468457d4a33Sstefano_zampini 
469457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr);
470457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
471457d4a33Sstefano_zampini     }
472457d4a33Sstefano_zampini     ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
473457d4a33Sstefano_zampini     ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr);
474af140850Sstefano_zampini 
475457d4a33Sstefano_zampini     /* local to global scatter for interface pressure */
476457d4a33Sstefano_zampini     ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr);
477457d4a33Sstefano_zampini     ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr);
478457d4a33Sstefano_zampini 
479af140850Sstefano_zampini     /* shift local to global indices for multipliers */
480457d4a33Sstefano_zampini     for (i=0;i<n_local_lambda;i++) {
481457d4a33Sstefano_zampini       PetscInt owner,ps;
482457d4a33Sstefano_zampini 
483457d4a33Sstefano_zampini       ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr);
484457d4a33Sstefano_zampini       ps = pranges[owner+1]-pranges[owner];
485457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
486457d4a33Sstefano_zampini     }
487457d4a33Sstefano_zampini 
4889cc7774eSstefano_zampini     /* scatter from alldofs to interface pressures global fetidp vector */
4899cc7774eSstefano_zampini     ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr);
4909cc7774eSstefano_zampini     ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr);
491af140850Sstefano_zampini     ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr);
4929cc7774eSstefano_zampini     ierr = ISDestroy(&ais);CHKERRQ(ierr);
493457d4a33Sstefano_zampini   }
494457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr);
495457d4a33Sstefano_zampini   ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr);
496457d4a33Sstefano_zampini   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
4979cc7774eSstefano_zampini   /* scatter from local to global multipliers */
498457d4a33Sstefano_zampini   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
499457d4a33Sstefano_zampini   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
500457d4a33Sstefano_zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr);
501457d4a33Sstefano_zampini 
502674ae819SStefano Zampini   /* Create some vectors needed by fetidp */
503674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
504674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
505674ae819SStefano Zampini 
506674ae819SStefano Zampini   test_fetidp = PETSC_FALSE;
507c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr);
508674ae819SStefano Zampini 
509674ae819SStefano Zampini   if (test_fetidp && !pcbddc->use_deluxe_scaling) {
5109eec4de8Sstefano_zampini     Vec       test_vec,test_vec_p = NULL;
511*2d456bbaSstefano_zampini     IS        dirdofs;
5125b08dc53SStefano Zampini     PetscReal real_value;
5135b08dc53SStefano Zampini 
514674ae819SStefano Zampini     ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
5151575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
516457d4a33Sstefano_zampini     ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP SIZES--------------\n");CHKERRQ(ierr);
517457d4a33Sstefano_zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%04d]: Sizes %D %D %D %D %D\n",rank,fetidpmat_ctx->n,fetidpmat_ctx->N,nPgl,nPg,nPl);CHKERRQ(ierr);
518457d4a33Sstefano_zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
519674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr);
520674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr);
521674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr);
522674ae819SStefano Zampini     if (fully_redundant) {
523674ae819SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr);
524674ae819SStefano Zampini     } else {
525674ae819SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr);
526674ae819SStefano Zampini     }
527674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
528674ae819SStefano Zampini 
529674ae819SStefano Zampini     /******************************************************************/
5309eec4de8Sstefano_zampini     /* TEST A/B: Test numbering of global fetidp dofs             */
531674ae819SStefano Zampini     /******************************************************************/
532674ae819SStefano Zampini 
533674ae819SStefano Zampini     ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr);
534457d4a33Sstefano_zampini     ierr = VecSet(fetidp_global,1.0);CHKERRQ(ierr);
5359eec4de8Sstefano_zampini     ierr = VecSet(test_vec,1.);CHKERRQ(ierr);
536457d4a33Sstefano_zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
537457d4a33Sstefano_zampini     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5389eec4de8Sstefano_zampini     if (fetidpmat_ctx->l2g_p) {
5399eec4de8Sstefano_zampini       ierr = VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);CHKERRQ(ierr);
5409eec4de8Sstefano_zampini       ierr = VecSet(test_vec_p,1.);CHKERRQ(ierr);
5419eec4de8Sstefano_zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5429eec4de8Sstefano_zampini       ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5439eec4de8Sstefano_zampini     }
544674ae819SStefano Zampini     scalar_value = -1.0;
545674ae819SStefano Zampini     ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
5465b08dc53SStefano Zampini     ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr);
547674ae819SStefano Zampini     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
5485b08dc53SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
5499eec4de8Sstefano_zampini     if (fetidpmat_ctx->l2g_p) {
5509eec4de8Sstefano_zampini       ierr = VecAXPY(test_vec_p,scalar_value,fetidpmat_ctx->vP);CHKERRQ(ierr);
5519eec4de8Sstefano_zampini       ierr = VecNorm(test_vec_p,NORM_INFINITY,&real_value);CHKERRQ(ierr);
5529eec4de8Sstefano_zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc (p): % 1.14e\n",rank,real_value);CHKERRQ(ierr);
5539eec4de8Sstefano_zampini     }
554674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
555674ae819SStefano Zampini     if (fully_redundant) {
556457d4a33Sstefano_zampini       ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
557674ae819SStefano Zampini       ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr);
558457d4a33Sstefano_zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
559457d4a33Sstefano_zampini       ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
560457d4a33Sstefano_zampini       ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr);
5615b08dc53SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
562674ae819SStefano Zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
563674ae819SStefano Zampini     }
5649eec4de8Sstefano_zampini     if (fetidpmat_ctx->l2g_p) {
565af140850Sstefano_zampini       ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr);
566af140850Sstefano_zampini       ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
567af140850Sstefano_zampini       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
568af140850Sstefano_zampini       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
569af140850Sstefano_zampini 
5709eec4de8Sstefano_zampini       ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
571af140850Sstefano_zampini       ierr = VecSet(fetidpmat_ctx->vP,-1.0);CHKERRQ(ierr);
572af140850Sstefano_zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
573af140850Sstefano_zampini       ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
574af140850Sstefano_zampini       ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
575af140850Sstefano_zampini       ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
576af140850Sstefano_zampini       ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
577af140850Sstefano_zampini       ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5789eec4de8Sstefano_zampini       ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr);
579af140850Sstefano_zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob (p): % 1.14e\n",rank,PetscRealPart(scalar_value));CHKERRQ(ierr);
580*2d456bbaSstefano_zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5819eec4de8Sstefano_zampini     }
582674ae819SStefano Zampini     /******************************************************************/
583674ae819SStefano Zampini     /* TEST C: It should holds B_delta*w=0, w\in\widehat{W}           */
584674ae819SStefano Zampini     /* This is the meaning of the B matrix                            */
585674ae819SStefano Zampini     /******************************************************************/
586674ae819SStefano Zampini 
587674ae819SStefano Zampini     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
588674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
589e176bc59SStefano Zampini     ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
590e176bc59SStefano Zampini     ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
591e176bc59SStefano Zampini     ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
592e176bc59SStefano Zampini     ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
593674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
594674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
595674ae819SStefano Zampini     /* Action of B_delta */
596674ae819SStefano Zampini     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
597457d4a33Sstefano_zampini     ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
598457d4a33Sstefano_zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
599457d4a33Sstefano_zampini     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
600457d4a33Sstefano_zampini     ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
6015b08dc53SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr);
602674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
603674ae819SStefano Zampini 
604674ae819SStefano Zampini     /******************************************************************/
605674ae819SStefano Zampini     /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W}     */
606674ae819SStefano Zampini     /* E_D = R_D^TR                                                   */
607674ae819SStefano Zampini     /* P_D = B_{D,delta}^T B_{delta}                                  */
608674ae819SStefano Zampini     /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
609674ae819SStefano Zampini     /******************************************************************/
610674ae819SStefano Zampini 
611674ae819SStefano Zampini     /* compute a random vector in \widetilde{W} */
612674ae819SStefano Zampini     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
613*2d456bbaSstefano_zampini     /* set zero at vertices and essential dofs */
614*2d456bbaSstefano_zampini     scalar_value = 0.0;
615674ae819SStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
616674ae819SStefano Zampini     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
617*2d456bbaSstefano_zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirdofs);CHKERRQ(ierr);
618*2d456bbaSstefano_zampini     if (dirdofs) {
619*2d456bbaSstefano_zampini       const PetscInt *idxs;
620*2d456bbaSstefano_zampini       PetscInt       ndir;
621*2d456bbaSstefano_zampini 
622*2d456bbaSstefano_zampini       ierr = ISGetLocalSize(dirdofs,&ndir);CHKERRQ(ierr);
623*2d456bbaSstefano_zampini       ierr = ISGetIndices(dirdofs,&idxs);CHKERRQ(ierr);
624*2d456bbaSstefano_zampini       for (i=0;i<ndir;i++) { array[idxs[i]]=scalar_value; }
625*2d456bbaSstefano_zampini       ierr = ISRestoreIndices(dirdofs,&idxs);CHKERRQ(ierr);
626*2d456bbaSstefano_zampini     }
627674ae819SStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
628674ae819SStefano Zampini     /* store w for final comparison */
629674ae819SStefano Zampini     ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr);
630674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
631674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632674ae819SStefano Zampini 
633674ae819SStefano Zampini     /* Jump operator P_D : results stored in pcis->vec1_B */
634674ae819SStefano Zampini 
635674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
636674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
637674ae819SStefano Zampini     /* Action of B_delta */
638674ae819SStefano Zampini     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
639457d4a33Sstefano_zampini     ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
640457d4a33Sstefano_zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
641457d4a33Sstefano_zampini     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
642674ae819SStefano Zampini     /* Action of B_Ddelta^T */
643457d4a33Sstefano_zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
644457d4a33Sstefano_zampini     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
645674ae819SStefano Zampini     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
646674ae819SStefano Zampini 
647674ae819SStefano Zampini     /* Average operator E_D : results stored in pcis->vec2_B */
648674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
649674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
650674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr);
651674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
652674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
653674ae819SStefano Zampini 
654674ae819SStefano Zampini     /* test E_D=I-P_D */
655674ae819SStefano Zampini     scalar_value = 1.0;
656674ae819SStefano Zampini     ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr);
657674ae819SStefano Zampini     scalar_value = -1.0;
658674ae819SStefano Zampini     ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr);
6595b08dc53SStefano Zampini     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr);
660674ae819SStefano Zampini     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
6615b08dc53SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
662674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
663674ae819SStefano Zampini 
664674ae819SStefano Zampini     /******************************************************************/
665*2d456bbaSstefano_zampini     /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W}           */
666674ae819SStefano Zampini     /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
667674ae819SStefano Zampini     /******************************************************************/
668674ae819SStefano Zampini 
669674ae819SStefano Zampini     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
670*2d456bbaSstefano_zampini     /* set zero at vertices and essential dofs */
671*2d456bbaSstefano_zampini     scalar_value = 0.0;
672674ae819SStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
673674ae819SStefano Zampini     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
674*2d456bbaSstefano_zampini     if (dirdofs) {
675*2d456bbaSstefano_zampini       const PetscInt *idxs;
676*2d456bbaSstefano_zampini       PetscInt       ndir;
677*2d456bbaSstefano_zampini 
678*2d456bbaSstefano_zampini       ierr = ISGetLocalSize(dirdofs,&ndir);CHKERRQ(ierr);
679*2d456bbaSstefano_zampini       ierr = ISGetIndices(dirdofs,&idxs);CHKERRQ(ierr);
680*2d456bbaSstefano_zampini       for (i=0;i<ndir;i++) { array[idxs[i]]=scalar_value; }
681*2d456bbaSstefano_zampini       ierr = ISRestoreIndices(dirdofs,&idxs);CHKERRQ(ierr);
682*2d456bbaSstefano_zampini     }
683674ae819SStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
684674ae819SStefano Zampini 
685674ae819SStefano Zampini     /* Jump operator P_D : results stored in pcis->vec1_B */
686674ae819SStefano Zampini 
687674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
688674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
689674ae819SStefano Zampini     /* Action of B_delta */
690674ae819SStefano Zampini     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
691457d4a33Sstefano_zampini     ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
692457d4a33Sstefano_zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
693457d4a33Sstefano_zampini     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
694674ae819SStefano Zampini     /* Action of B_Ddelta^T */
695457d4a33Sstefano_zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
696457d4a33Sstefano_zampini     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
697674ae819SStefano Zampini     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
698674ae819SStefano Zampini     /* scaling */
699674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
7005b08dc53SStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
7015b08dc53SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr);
702674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
703674ae819SStefano Zampini 
704674ae819SStefano Zampini     if (!fully_redundant) {
705674ae819SStefano Zampini       /******************************************************************/
706674ae819SStefano Zampini       /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
707674ae819SStefano Zampini       /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
708674ae819SStefano Zampini       /******************************************************************/
709457d4a33Sstefano_zampini       ierr = VecDuplicate(fetidp_global,&test_vec);CHKERRQ(ierr);
710457d4a33Sstefano_zampini       ierr = VecSetRandom(fetidp_global,NULL);CHKERRQ(ierr);
711457d4a33Sstefano_zampini       if (fetidpmat_ctx->l2g_p) {
712457d4a33Sstefano_zampini         ierr = VecSet(fetidpmat_ctx->vP,0.);CHKERRQ(ierr);
713457d4a33Sstefano_zampini         ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
714457d4a33Sstefano_zampini         ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
715457d4a33Sstefano_zampini       }
716674ae819SStefano Zampini       /* Action of B_Ddelta^T */
717457d4a33Sstefano_zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
718457d4a33Sstefano_zampini       ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
719674ae819SStefano Zampini       ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
720674ae819SStefano Zampini       /* Action of B_delta */
721674ae819SStefano Zampini       ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
722674ae819SStefano Zampini       ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
723674ae819SStefano Zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
724674ae819SStefano Zampini       ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
725674ae819SStefano Zampini       scalar_value = -1.0;
726457d4a33Sstefano_zampini       ierr = VecAXPY(fetidp_global,scalar_value,test_vec);CHKERRQ(ierr);
727457d4a33Sstefano_zampini       ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
7285b08dc53SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr);
729674ae819SStefano Zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
730674ae819SStefano Zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
731674ae819SStefano Zampini       ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
732674ae819SStefano Zampini     }
7339eec4de8Sstefano_zampini     ierr = VecDestroy(&test_vec_p);CHKERRQ(ierr);
734*2d456bbaSstefano_zampini     ierr = ISDestroy(&dirdofs);CHKERRQ(ierr);
735674ae819SStefano Zampini   }
736674ae819SStefano Zampini   /* final cleanup */
737457d4a33Sstefano_zampini   ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr);
738674ae819SStefano Zampini   PetscFunctionReturn(0);
739674ae819SStefano Zampini }
740674ae819SStefano Zampini 
741674ae819SStefano Zampini #undef __FUNCT__
742674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext"
743674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
744674ae819SStefano Zampini {
745674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
746ed6c3d69SStefano Zampini   PC_IS          *pcis;
747f28b6018SStefano Zampini   PetscBool      lumped = PETSC_FALSE;
748674ae819SStefano Zampini   PetscErrorCode ierr;
749674ae819SStefano Zampini 
750674ae819SStefano Zampini   PetscFunctionBegin;
751674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
752674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
753674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
754674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
755674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
756674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
757674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
758674ae819SStefano Zampini   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
759f28b6018SStefano Zampini   /* create preconditioner */
760ed6c3d69SStefano Zampini   pcis = (PC_IS*)fetidppc_ctx->pc->data;
761f28b6018SStefano Zampini   /* Dirichlet preconditioner */
762f28b6018SStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_lumped",&lumped,NULL);CHKERRQ(ierr);
763f28b6018SStefano Zampini   if (!lumped) {
764ed6c3d69SStefano Zampini     ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
765ed6c3d69SStefano Zampini     ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
766f28b6018SStefano Zampini   } else {
767f28b6018SStefano Zampini     ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr);
768f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
769f28b6018SStefano Zampini   }
770af140850Sstefano_zampini   /* saddle-point */
771af140850Sstefano_zampini   if (mat_ctx->xPg) {
772af140850Sstefano_zampini     Mat PAmat = NULL,PPmat = NULL;
773af140850Sstefano_zampini 
774af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr);
775af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
776af140850Sstefano_zampini     ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr);
777af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
778af140850Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);CHKERRQ(ierr);
779af140850Sstefano_zampini     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PAmat",(PetscObject*)&PAmat);CHKERRQ(ierr);
780af140850Sstefano_zampini     ierr = KSPCreate(PetscObjectComm((PetscObject)fetidppc_ctx->pc),&fetidppc_ctx->kP);CHKERRQ(ierr);
781af140850Sstefano_zampini     ierr = KSPSetOperators(fetidppc_ctx->kP,PAmat,PPmat);CHKERRQ(ierr);
782af140850Sstefano_zampini     ierr = KSPSetOptionsPrefix(fetidppc_ctx->kP,"fetidp_pmass_");CHKERRQ(ierr);
783af140850Sstefano_zampini     ierr = KSPSetFromOptions(fetidppc_ctx->kP);CHKERRQ(ierr);
784af140850Sstefano_zampini   }
785674ae819SStefano Zampini   PetscFunctionReturn(0);
786674ae819SStefano Zampini }
787674ae819SStefano Zampini 
788674ae819SStefano Zampini #undef __FUNCT__
789674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult"
790674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
791674ae819SStefano Zampini {
792674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
793617d11aeSStefano Zampini   PC_BDDC        *pcbddc;
794674ae819SStefano Zampini   PC_IS          *pcis;
795674ae819SStefano Zampini   PetscErrorCode ierr;
796674ae819SStefano Zampini 
797674ae819SStefano Zampini   PetscFunctionBegin;
798674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
799674ae819SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
800617d11aeSStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
801674ae819SStefano Zampini   /* Application of B_delta^T */
802af140850Sstefano_zampini   ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr);
803674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
804674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
805674ae819SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
806af140850Sstefano_zampini 
807af140850Sstefano_zampini   /* Add contribution from saddle point */
808af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
809af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
810af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
811af140850Sstefano_zampini     if (pcbddc->switch_static) {
812af140850Sstefano_zampini       ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr);
813af140850Sstefano_zampini     }
814af140850Sstefano_zampini     ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
815af140850Sstefano_zampini   } else {
816af140850Sstefano_zampini     if (pcbddc->switch_static) {
817674ae819SStefano Zampini       ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
818af140850Sstefano_zampini     }
819af140850Sstefano_zampini   }
820af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
821617d11aeSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
822dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
823c7ffc8ceSStefano Zampini   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
824af140850Sstefano_zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
825674ae819SStefano Zampini   /* Application of B_delta */
826674ae819SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
827af140850Sstefano_zampini   /* Contribution from boundary pressures */
828af140850Sstefano_zampini   if (mat_ctx->C) {
829af140850Sstefano_zampini     const PetscScalar *lx;
830af140850Sstefano_zampini     PetscScalar       *ly;
831af140850Sstefano_zampini 
832af140850Sstefano_zampini     /* pressures ordered first in x and y */
833af140850Sstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
834af140850Sstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
835af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr);
836af140850Sstefano_zampini     ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr);
837af140850Sstefano_zampini     ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr);
838af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr);
839af140850Sstefano_zampini     ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr);
840af140850Sstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
841af140850Sstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
842af140850Sstefano_zampini   }
843af140850Sstefano_zampini   /* Add contribution from saddle point */
844af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
845af140850Sstefano_zampini     ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
846af140850Sstefano_zampini     if (pcbddc->switch_static) {
847af140850Sstefano_zampini       ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
848af140850Sstefano_zampini     }
849af140850Sstefano_zampini     ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
850af140850Sstefano_zampini     ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
851af140850Sstefano_zampini   }
852674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
853674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
854674ae819SStefano Zampini   PetscFunctionReturn(0);
855674ae819SStefano Zampini }
856674ae819SStefano Zampini 
857674ae819SStefano Zampini #undef __FUNCT__
858edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose"
859edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
860edf7251bSStefano Zampini {
861edf7251bSStefano Zampini   FETIDPMat_ctx  mat_ctx;
862edf7251bSStefano Zampini   PC_IS          *pcis;
863edf7251bSStefano Zampini   PetscErrorCode ierr;
864edf7251bSStefano Zampini 
865edf7251bSStefano Zampini   PetscFunctionBegin;
866edf7251bSStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
867edf7251bSStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
868edf7251bSStefano Zampini   /* Application of B_delta^T */
869edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
870edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
871edf7251bSStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
872edf7251bSStefano Zampini   /* Application of \widetilde{S}^-1 */
873edf7251bSStefano Zampini   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
874edf7251bSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr);
875edf7251bSStefano Zampini   /* Application of B_delta */
876edf7251bSStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
877edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
878edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
879edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
880edf7251bSStefano Zampini   PetscFunctionReturn(0);
881edf7251bSStefano Zampini }
882edf7251bSStefano Zampini 
883edf7251bSStefano Zampini #undef __FUNCT__
884674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply"
885674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
886674ae819SStefano Zampini {
887674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
888674ae819SStefano Zampini   PC_IS          *pcis;
889674ae819SStefano Zampini   PetscErrorCode ierr;
890674ae819SStefano Zampini 
891674ae819SStefano Zampini   PetscFunctionBegin;
892302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
893674ae819SStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
894674ae819SStefano Zampini   /* Application of B_Ddelta^T */
895674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
896674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
897674ae819SStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
898674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
899ed6c3d69SStefano Zampini   /* Application of local Schur complement */
900ed6c3d69SStefano Zampini   ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
901edf7251bSStefano Zampini   /* Application of B_Ddelta */
902edf7251bSStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
903edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
904edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
905edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
906b519380fSstefano_zampini   /* interface pressure preconditioner */
907b519380fSstefano_zampini   if (pc_ctx->kP) {
908b519380fSstefano_zampini     const PetscScalar *lx;
909b519380fSstefano_zampini     PetscScalar *ly;
910b519380fSstefano_zampini 
911b519380fSstefano_zampini     /* pressures ordered first in x and y */
912b519380fSstefano_zampini     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
913b519380fSstefano_zampini     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
914b519380fSstefano_zampini     ierr = VecPlaceArray(pc_ctx->xPg,lx);CHKERRQ(ierr);
915b519380fSstefano_zampini     ierr = VecPlaceArray(pc_ctx->yPg,ly);CHKERRQ(ierr);
916b519380fSstefano_zampini     ierr = KSPSolve(pc_ctx->kP,pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr);
917b519380fSstefano_zampini     ierr = VecCopy(pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr);
918b519380fSstefano_zampini     ierr = VecResetArray(pc_ctx->xPg);CHKERRQ(ierr);
919b519380fSstefano_zampini     ierr = VecResetArray(pc_ctx->yPg);CHKERRQ(ierr);
920b519380fSstefano_zampini     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
921b519380fSstefano_zampini     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
922b519380fSstefano_zampini   }
923edf7251bSStefano Zampini   PetscFunctionReturn(0);
924edf7251bSStefano Zampini }
925edf7251bSStefano Zampini 
926edf7251bSStefano Zampini #undef __FUNCT__
927edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose"
928edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
929edf7251bSStefano Zampini {
930edf7251bSStefano Zampini   FETIDPPC_ctx   pc_ctx;
931edf7251bSStefano Zampini   PC_IS          *pcis;
932edf7251bSStefano Zampini   PetscErrorCode ierr;
933edf7251bSStefano Zampini 
934edf7251bSStefano Zampini   PetscFunctionBegin;
935302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
936edf7251bSStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
937edf7251bSStefano Zampini   /* Application of B_Ddelta^T */
938edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
939edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
940edf7251bSStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
941edf7251bSStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
942ed6c3d69SStefano Zampini   /* Application of local Schur complement */
943ed6c3d69SStefano Zampini   ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
944674ae819SStefano Zampini   /* Application of B_Ddelta */
945674ae819SStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
946674ae819SStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
947674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
948674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
949674ae819SStefano Zampini   PetscFunctionReturn(0);
950674ae819SStefano Zampini }
951