xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision af1408500dea18abb78742575d10d03b134ffd5c)
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PCBDDCCreateFETIDPMatContext"
6 PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
7 {
8   FETIDPMat_ctx  newctx;
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   ierr = PetscNew(&newctx);CHKERRQ(ierr);
13   /* increase the reference count for BDDC preconditioner */
14   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
15   newctx->pc              = pc;
16   *fetidpmat_ctx          = newctx;
17   PetscFunctionReturn(0);
18 }
19 
20 #undef __FUNCT__
21 #define __FUNCT__ "PCBDDCCreateFETIDPPCContext"
22 PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
23 {
24   FETIDPPC_ctx   newctx;
25   PetscErrorCode ierr;
26 
27   PetscFunctionBegin;
28   ierr = PetscNew(&newctx);CHKERRQ(ierr);
29   /* increase the reference count for BDDC preconditioner */
30   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
31   newctx->pc              = pc;
32   *fetidppc_ctx           = newctx;
33   PetscFunctionReturn(0);
34 }
35 
36 #undef __FUNCT__
37 #define __FUNCT__ "PCBDDCDestroyFETIDPMat"
38 PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
39 {
40   FETIDPMat_ctx  mat_ctx;
41   PetscErrorCode ierr;
42 
43   PetscFunctionBegin;
44   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
45   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
46   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
47   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
48   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
49   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
50   ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr);
51   ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr);
52   ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr);
53   ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr);
54   ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr);
55   ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr);
56   ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr);
57   ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr);
58   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
59   ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr);
60   ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr);
61   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
62   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "PCBDDCDestroyFETIDPPC"
68 PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
69 {
70   FETIDPPC_ctx   pc_ctx;
71   PetscErrorCode ierr;
72 
73   PetscFunctionBegin;
74   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
75   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
76   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
77   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
78   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
79   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
80   ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr);
81   ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr);
82   ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr);
83   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "PCBDDCSetupFETIDPMatContext"
89 PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
90 {
91   PetscErrorCode ierr;
92   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
93   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
94   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
95   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
96   MPI_Comm       comm;
97   Mat            ScalingMat;
98   Vec            fetidp_global;
99   IS             IS_l2g_lambda;
100   IS             subset,subset_mult,subset_n;
101   PetscBool      skip_node,fully_redundant;
102   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
103   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
104   PetscMPIInt    rank,size,buf_size,neigh;
105   PetscScalar    scalar_value;
106   PetscInt       *vertex_indices;
107   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
108   const PetscInt *aux_global_numbering;
109   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
110   PetscScalar    *array,*scaling_factors,*vals_B_delta;
111   PetscInt       *aux_local_numbering_2;
112   PetscLayout    llay;
113   /* For communication of scaling factors */
114   PetscInt       *ptrs_buffer,neigh_position;
115   PetscScalar    **all_factors,*send_buffer,*recv_buffer;
116   MPI_Request    *send_reqs,*recv_reqs;
117   /* tests */
118   PetscBool      test_fetidp;
119   PetscViewer    viewer;
120   /* saddlepoint */
121   Mat                    oA;
122   ISLocalToGlobalMapping l2gmap_p;
123   PetscLayout            play;
124   IS                     gP,pP;
125   PetscInt               nPl,nPg,nPgl;
126 
127   PetscFunctionBegin;
128   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
129   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
130   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
131 
132   /* saddlepoint */
133   nPl      = 0;
134   nPg      = 0;
135   nPgl     = 0;
136   gP       = NULL;
137   pP       = NULL;
138   l2gmap_p = NULL;
139   play     = NULL;
140   oA       = NULL;
141   ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_A",(PetscObject*)&oA);CHKERRQ(ierr);
142   if (oA) { /* oA is the original matrix */
143     /* pressures in parallel layout of oA */
144     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr);
145     if (!pP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pP not present");
146 
147     /* subdomain pressures in global numbering */
148     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr);
149     if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
150     ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr);
151     ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr);
152     ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr);
153     ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr);
154     ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr);
155 
156     /* interface pressure matrix */
157     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr);
158     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */
159       IS Pg;
160 
161       ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr);
162       ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr);
163       ierr = ISDestroy(&Pg);CHKERRQ(ierr);
164       ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr);
165       ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr);
166       ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr);
167       ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr);
168       ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr);
169       ierr = PetscLayoutSetUp(play);CHKERRQ(ierr);
170     } else {
171       ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr);
172       ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr);
173       ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr);
174       ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr);
175       ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr);
176       ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr);
177       ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr);
178     }
179     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr);
180     ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr);
181 
182     /* import matrices for interface pressures coupling */
183     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr);
184     if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
185     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr);
186 
187     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr);
188     if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
189     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr);
190 
191     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
192     if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
193     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr);
194 
195     ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
196     if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
197     ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr);
198   }
199 
200   /* Default type of lagrange multipliers is non-redundant */
201   fully_redundant = fetidpmat_ctx->fully_redundant;
202 
203   /* Evaluate local and global number of lagrange multipliers */
204   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
205   n_local_lambda = 0;
206   partial_sum = 0;
207   n_boundary_dofs = 0;
208   s = 0;
209   /* Get Vertices used to define the BDDC */
210   n_vertices = pcbddc->n_vertices;
211   vertex_indices = pcbddc->local_primal_ref_node;
212   dual_size = pcis->n_B-n_vertices;
213   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
214   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
215   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
216 
217   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
218   for (i=0;i<pcis->n;i++){
219     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
220     if ( j > 0 ) {
221       n_boundary_dofs++;
222     }
223     skip_node = PETSC_FALSE;
224     if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */
225       skip_node = PETSC_TRUE;
226       s++;
227     }
228     if (j < 1) {
229       skip_node = PETSC_TRUE;
230     }
231     if ( !skip_node ) {
232       if (fully_redundant) {
233         /* fully redundant set of lagrange multipliers */
234         n_lambda_for_dof = (j*(j+1))/2;
235       } else {
236         n_lambda_for_dof = j;
237       }
238       n_local_lambda += j;
239       /* needed to evaluate global number of lagrange multipliers */
240       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
241       /* store some data needed */
242       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
243       aux_local_numbering_1[partial_sum] = i;
244       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
245       partial_sum++;
246     }
247   }
248   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
249 
250   /* compute global ordering of lagrange multipliers and associate l2g map */
251   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
252   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
253   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
254   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
255   ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr);
256   ierr = ISDestroy(&subset);CHKERRQ(ierr);
257 
258 #if defined(PETSC_USE_DEBUG)
259   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
260   ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
261   ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
262   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
263   i = (PetscInt)PetscRealPart(scalar_value);
264   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);
265 #endif
266 
267   /* init data for scaling factors exchange */
268   partial_sum = 0;
269   ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
270   ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr);
271   ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr);
272   ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr);
273   if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
274   for (i=1;i<pcis->n_neigh;i++) {
275     partial_sum += pcis->n_shared[i];
276     ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
277   }
278   ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
279   ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
280   ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
281   for (i=0;i<pcis->n-1;i++) {
282     j = mat_graph->count[i];
283     all_factors[i+1]=all_factors[i]+j;
284   }
285   /* scatter B scaling to N vec */
286   ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
287   ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
288   /* communications */
289   ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
290   for (i=1;i<pcis->n_neigh;i++) {
291     for (j=0;j<pcis->n_shared[i];j++) {
292       send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
293     }
294     ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
295     ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
296     ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
297     ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
298   }
299   ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
300   if (pcis->n_neigh > 0) {
301     ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
302   }
303   /* put values in correct places */
304   for (i=1;i<pcis->n_neigh;i++) {
305     for (j=0;j<pcis->n_shared[i];j++) {
306       k = pcis->shared[i][j];
307       neigh_position = 0;
308       while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
309       all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
310     }
311   }
312   if (pcis->n_neigh > 0) {
313     ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
314   }
315   ierr = PetscFree(send_reqs);CHKERRQ(ierr);
316   ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
317   ierr = PetscFree(send_buffer);CHKERRQ(ierr);
318   ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
319   ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
320 
321   /* Compute B and B_delta (local actions) */
322   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
323   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
324   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
325   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
326   ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
327   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
328   partial_sum=0;
329   cum = 0;
330   for (i=0;i<dual_size;i++) {
331     n_global_lambda = aux_global_numbering[cum];
332     j = mat_graph->count[aux_local_numbering_1[i]];
333     aux_sums[0]=0;
334     for (s=1;s<j;s++) {
335       aux_sums[s]=aux_sums[s-1]+j-s+1;
336     }
337     array = all_factors[aux_local_numbering_1[i]];
338     n_neg_values = 0;
339     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
340     n_pos_values = j - n_neg_values;
341     if (fully_redundant) {
342       for (s=0;s<n_neg_values;s++) {
343         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
344         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
345         vals_B_delta   [partial_sum+s]=-1.0;
346         scaling_factors[partial_sum+s]=array[s];
347       }
348       for (s=0;s<n_pos_values;s++) {
349         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
350         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
351         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
352         scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
353       }
354       partial_sum += j;
355     } else {
356       /* l2g_indices and default cols and vals of B_delta */
357       for (s=0;s<j;s++) {
358         l2g_indices    [partial_sum+s]=n_global_lambda+s;
359         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
360         vals_B_delta   [partial_sum+s]=0.0;
361       }
362       /* B_delta */
363       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
364         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
365       }
366       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
367         vals_B_delta   [partial_sum+n_neg_values]=1.0;
368       }
369       /* scaling as in Klawonn-Widlund 1999*/
370       for (s=0;s<n_neg_values;s++) {
371         scalar_value = 0.0;
372         for (k=0;k<s+1;k++) {
373           scalar_value += array[k];
374         }
375         scaling_factors[partial_sum+s] = -scalar_value;
376       }
377       for (s=0;s<n_pos_values;s++) {
378         scalar_value = 0.0;
379         for (k=s+n_neg_values;k<j;k++) {
380           scalar_value += array[k];
381         }
382         scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
383       }
384       partial_sum += j;
385     }
386     cum += aux_local_numbering_2[i];
387   }
388   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
389   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
390   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
391   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
392   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
393   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
394   ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
395   ierr = PetscFree(all_factors);CHKERRQ(ierr);
396 
397   /* Create local part of B_delta */
398   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
399   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
400   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
401   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
402   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
403   for (i=0;i<n_local_lambda;i++) {
404     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
405   }
406   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
407   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
408   ierr = MatAssemblyEnd  (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409 
410   if (fully_redundant) {
411     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
412     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
413     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
414     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
415     for (i=0;i<n_local_lambda;i++) {
416       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
417     }
418     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
419     ierr = MatAssemblyEnd  (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
420     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
421     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
422   } else {
423     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
424     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
425     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
426     ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
427     for (i=0;i<n_local_lambda;i++) {
428       ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
429     }
430     ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
431     ierr = MatAssemblyEnd  (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
432   }
433   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
434   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
435 
436   /* Layout of multipliers */
437   ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr);
438   ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr);
439   ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
440   ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr);
441   ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr);
442 
443   /* fetidpmat sizes */
444   fetidpmat_ctx->n += nPgl;
445   fetidpmat_ctx->N  = fetidpmat_ctx->n_lambda+nPg;
446 
447   /* Local work vector of multipliers */
448   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
449   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
450   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
451 
452   /* Global vector for FETI-DP linear system */
453   ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr);
454   ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr);
455   ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr);
456   ierr = VecSetUp(fetidp_global);CHKERRQ(ierr);
457 
458   /* Decide layout for fetidp dofs: if it is a saddle point problem
459      pressure is ordered first in the local part of the global vector
460      of the FETI-DP linear system */
461   if (nPg) {
462     IS             IS_l2g_p,ais;
463     PetscLayout    alay;
464     const PetscInt *idxs,*pranges,*aranges,*lranges;
465     PetscInt       *l2g_indices_p,rst;
466 
467     ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr);
468     ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr);
469     ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr);
470     ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr);
471     ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr);
472     ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
473     /* shift local to global indices for pressure */
474     for (i=0;i<nPl;i++) {
475       PetscInt owner;
476 
477       ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr);
478       l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
479     }
480     ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr);
481     ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr);
482 
483     /* local to global scatter for interface pressure */
484     ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr);
485     ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr);
486 
487     /* shift local to global indices for multipliers */
488     for (i=0;i<n_local_lambda;i++) {
489       PetscInt owner,ps;
490 
491       ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr);
492       ps = pranges[owner+1]-pranges[owner];
493       l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
494     }
495 
496     /* scatter from alldofs to interface pressures global fetidp vector */
497     ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr);
498     ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr);
499     ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr);
500     ierr = ISDestroy(&ais);CHKERRQ(ierr);
501   }
502   ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr);
503   ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr);
504   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
505   /* scatter from local to global multipliers */
506   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
507   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
508   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr);
509 
510   /* Create some vectors needed by fetidp */
511   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
512   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
513 
514   test_fetidp = PETSC_FALSE;
515   ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr);
516 
517   if (test_fetidp && !pcbddc->use_deluxe_scaling) {
518     Vec       test_vec,test_vec_p = NULL;
519     PetscReal real_value;
520 
521     ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
522     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
523     ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP SIZES--------------\n");CHKERRQ(ierr);
524     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%04d]: Sizes %D %D %D %D %D\n",rank,fetidpmat_ctx->n,fetidpmat_ctx->N,nPgl,nPg,nPl);CHKERRQ(ierr);
525     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
526     ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr);
527     ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr);
528     ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr);
529     if (fully_redundant) {
530       ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr);
531     } else {
532       ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr);
533     }
534     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
535 
536     /******************************************************************/
537     /* TEST A/B: Test numbering of global fetidp dofs             */
538     /******************************************************************/
539 
540     ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr);
541     ierr = VecSet(fetidp_global,1.0);CHKERRQ(ierr);
542     ierr = VecSet(test_vec,1.);CHKERRQ(ierr);
543     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
544     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
545     if (fetidpmat_ctx->l2g_p) {
546       ierr = VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);CHKERRQ(ierr);
547       ierr = VecSet(test_vec_p,1.);CHKERRQ(ierr);
548       ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
549       ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
550     }
551     scalar_value = -1.0;
552     ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
553     ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr);
554     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
555     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
556     if (fetidpmat_ctx->l2g_p) {
557       ierr = VecAXPY(test_vec_p,scalar_value,fetidpmat_ctx->vP);CHKERRQ(ierr);
558       ierr = VecNorm(test_vec_p,NORM_INFINITY,&real_value);CHKERRQ(ierr);
559       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc (p): % 1.14e\n",rank,real_value);CHKERRQ(ierr);
560     }
561     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
562     if (fully_redundant) {
563       ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
564       ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr);
565       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
566       ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
567       ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr);
568       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
569       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
570     }
571     if (fetidpmat_ctx->l2g_p) {
572       ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr);
573       ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
574       ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
575       ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
576 
577       ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
578       ierr = VecSet(fetidpmat_ctx->vP,-1.0);CHKERRQ(ierr);
579       ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
580       ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
581       ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
582       ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
583       ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
584       ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
585       ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr);
586       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob (p): % 1.14e\n",rank,PetscRealPart(scalar_value));CHKERRQ(ierr);
587     }
588     /******************************************************************/
589     /* TEST C: It should holds B_delta*w=0, w\in\widehat{W}           */
590     /* This is the meaning of the B matrix                            */
591     /******************************************************************/
592 
593     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
594     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
595     ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
596     ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
597     ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598     ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
599     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
600     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
601     /* Action of B_delta */
602     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
603     ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
604     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
605     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
606     ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
607     ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr);
608     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
609 
610     /******************************************************************/
611     /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W}     */
612     /* E_D = R_D^TR                                                   */
613     /* P_D = B_{D,delta}^T B_{delta}                                  */
614     /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
615     /******************************************************************/
616 
617     /* compute a random vector in \widetilde{W} */
618     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
619     scalar_value = 0.0;  /* set zero at vertices */
620     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
621     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
622     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
623     /* store w for final comparison */
624     ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr);
625     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627 
628     /* Jump operator P_D : results stored in pcis->vec1_B */
629 
630     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
631     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632     /* Action of B_delta */
633     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
634     ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
635     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
636     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
637     /* Action of B_Ddelta^T */
638     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
639     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
640     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
641 
642     /* Average operator E_D : results stored in pcis->vec2_B */
643     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
644     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
645     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr);
646     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
647     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
648 
649     /* test E_D=I-P_D */
650     scalar_value = 1.0;
651     ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr);
652     scalar_value = -1.0;
653     ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr);
654     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr);
655     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
656     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
657     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
658 
659     /******************************************************************/
660     /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W}          */
661     /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
662     /******************************************************************/
663 
664     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
665     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
666     scalar_value = 0.0;  /* set zero at vertices */
667     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
668     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
669 
670     /* Jump operator P_D : results stored in pcis->vec1_B */
671 
672     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
673     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674     /* Action of B_delta */
675     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
676     ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr);
677     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
678     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
679     /* Action of B_Ddelta^T */
680     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
681     ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
682     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
683     /* scaling */
684     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
685     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
686     ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr);
687     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
688 
689     if (!fully_redundant) {
690       /******************************************************************/
691       /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
692       /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
693       /******************************************************************/
694       ierr = VecDuplicate(fetidp_global,&test_vec);CHKERRQ(ierr);
695       ierr = VecSetRandom(fetidp_global,NULL);CHKERRQ(ierr);
696       if (fetidpmat_ctx->l2g_p) {
697         ierr = VecSet(fetidpmat_ctx->vP,0.);CHKERRQ(ierr);
698         ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
699         ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
700       }
701       /* Action of B_Ddelta^T */
702       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
703       ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
704       ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
705       /* Action of B_delta */
706       ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
707       ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
708       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
709       ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
710       scalar_value = -1.0;
711       ierr = VecAXPY(fetidp_global,scalar_value,test_vec);CHKERRQ(ierr);
712       ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
713       ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr);
714       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
715       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
716       ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
717     }
718     ierr = VecDestroy(&test_vec_p);CHKERRQ(ierr);
719   }
720   /* final cleanup */
721   ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr);
722   PetscFunctionReturn(0);
723 }
724 
725 #undef __FUNCT__
726 #define __FUNCT__ "PCBDDCSetupFETIDPPCContext"
727 PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
728 {
729   FETIDPMat_ctx  mat_ctx;
730   PC_IS          *pcis;
731   PetscBool      lumped = PETSC_FALSE;
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
736   /* get references from objects created when setting up feti mat context */
737   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
738   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
739   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
740   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
741   ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
742   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
743   /* create preconditioner */
744   pcis = (PC_IS*)fetidppc_ctx->pc->data;
745   /* Dirichlet preconditioner */
746   ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_lumped",&lumped,NULL);CHKERRQ(ierr);
747   if (!lumped) {
748     ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
749     ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
750   } else {
751     ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr);
752     fetidppc_ctx->S_j = pcis->A_BB;
753   }
754   /* saddle-point */
755   if (mat_ctx->xPg) {
756     Mat PAmat = NULL,PPmat = NULL;
757 
758     ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr);
759     fetidppc_ctx->xPg = mat_ctx->xPg;
760     ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr);
761     fetidppc_ctx->yPg = mat_ctx->yPg;
762     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);CHKERRQ(ierr);
763     ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PAmat",(PetscObject*)&PAmat);CHKERRQ(ierr);
764     ierr = KSPCreate(PetscObjectComm((PetscObject)fetidppc_ctx->pc),&fetidppc_ctx->kP);CHKERRQ(ierr);
765     ierr = KSPSetOperators(fetidppc_ctx->kP,PAmat,PPmat);CHKERRQ(ierr);
766     ierr = KSPSetOptionsPrefix(fetidppc_ctx->kP,"fetidp_pmass_");CHKERRQ(ierr);
767     ierr = KSPSetFromOptions(fetidppc_ctx->kP);CHKERRQ(ierr);
768   }
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "FETIDPMatMult"
774 PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
775 {
776   FETIDPMat_ctx  mat_ctx;
777   PC_BDDC        *pcbddc;
778   PC_IS          *pcis;
779   PetscErrorCode ierr;
780 
781   PetscFunctionBegin;
782   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
783   pcis = (PC_IS*)mat_ctx->pc->data;
784   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
785   /* Application of B_delta^T */
786   ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr);
787   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
788   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
789   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
790 
791   /* Add contribution from saddle point */
792   if (mat_ctx->l2g_p) {
793     ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
794     ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
795     if (pcbddc->switch_static) {
796       ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr);
797     }
798     ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr);
799   } else {
800     if (pcbddc->switch_static) {
801       ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
802     }
803   }
804   /* Application of \widetilde{S}^-1 */
805   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
806   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
807   ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr);
808   ierr = VecSet(y,0.0);CHKERRQ(ierr);
809   /* Application of B_delta */
810   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
811   /* Contribution from boundary pressures */
812   if (mat_ctx->C) {
813     const PetscScalar *lx;
814     PetscScalar       *ly;
815 
816     /* pressures ordered first in x and y */
817     ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr);
818     ierr = VecGetArray(y,&ly);CHKERRQ(ierr);
819     ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr);
820     ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr);
821     ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr);
822     ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr);
823     ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr);
824     ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr);
825     ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr);
826   }
827   /* Add contribution from saddle point */
828   if (mat_ctx->l2g_p) {
829     ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr);
830     if (pcbddc->switch_static) {
831       ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr);
832     }
833     ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
834     ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
835   }
836   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
837   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "FETIDPMatMultTranspose"
843 PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
844 {
845   FETIDPMat_ctx  mat_ctx;
846   PC_IS          *pcis;
847   PetscErrorCode ierr;
848 
849   PetscFunctionBegin;
850   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
851   pcis = (PC_IS*)mat_ctx->pc->data;
852   /* Application of B_delta^T */
853   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
854   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
855   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
856   /* Application of \widetilde{S}^-1 */
857   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
858   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr);
859   /* Application of B_delta */
860   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
861   ierr = VecSet(y,0.0);CHKERRQ(ierr);
862   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
863   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
864   PetscFunctionReturn(0);
865 }
866 
867 #undef __FUNCT__
868 #define __FUNCT__ "FETIDPPCApply"
869 PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
870 {
871   FETIDPPC_ctx   pc_ctx;
872   PC_IS          *pcis;
873   PetscErrorCode ierr;
874 
875   PetscFunctionBegin;
876   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
877   pcis = (PC_IS*)pc_ctx->pc->data;
878   /* Application of B_Ddelta^T */
879   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
880   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
881   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
882   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
883   /* Application of local Schur complement */
884   ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
885   /* Application of B_Ddelta */
886   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
887   ierr = VecSet(y,0.0);CHKERRQ(ierr);
888   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
889   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
890   PetscFunctionReturn(0);
891 }
892 
893 #undef __FUNCT__
894 #define __FUNCT__ "FETIDPPCApplyTranspose"
895 PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
896 {
897   FETIDPPC_ctx   pc_ctx;
898   PC_IS          *pcis;
899   PetscErrorCode ierr;
900 
901   PetscFunctionBegin;
902   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
903   pcis = (PC_IS*)pc_ctx->pc->data;
904   /* Application of B_Ddelta^T */
905   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
906   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
907   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
908   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
909   /* Application of local Schur complement */
910   ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
911   /* Application of B_Ddelta */
912   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
913   ierr = VecSet(y,0.0);CHKERRQ(ierr);
914   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
915   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
916   PetscFunctionReturn(0);
917 }
918