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