xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
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   newctx->lambda_local    = 0;
14   newctx->temp_solution_B = 0;
15   newctx->temp_solution_D = 0;
16   newctx->B_delta         = 0;
17   newctx->B_Ddelta        = 0; /* theoretically belongs to the FETIDP preconditioner */
18   newctx->l2g_lambda      = 0;
19   /* increase the reference count for BDDC preconditioner */
20   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
21   newctx->pc              = pc;
22   *fetidpmat_ctx          = newctx;
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "PCBDDCCreateFETIDPPCContext"
28 PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
29 {
30   FETIDPPC_ctx   newctx;
31   PetscErrorCode ierr;
32 
33   PetscFunctionBegin;
34   ierr = PetscNew(&newctx);CHKERRQ(ierr);
35   newctx->lambda_local    = 0;
36   newctx->B_Ddelta        = 0;
37   newctx->l2g_lambda      = 0;
38   /* increase the reference count for BDDC preconditioner */
39   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
40   newctx->pc              = pc;
41   *fetidppc_ctx           = newctx;
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "PCBDDCDestroyFETIDPMat"
47 PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
48 {
49   FETIDPMat_ctx  mat_ctx;
50   PetscErrorCode ierr;
51 
52   PetscFunctionBegin;
53   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
54   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
55   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
56   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
57   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
58   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
59   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
60   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
61   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "PCBDDCDestroyFETIDPPC"
67 PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
68 {
69   FETIDPPC_ctx   pc_ctx;
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
74   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
75   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
76   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
77   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
78   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
79   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "PCBDDCSetupFETIDPMatContext"
85 PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
86 {
87   PetscErrorCode ierr;
88   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
89   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
90   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
91   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
92   MPI_Comm       comm;
93   Mat            ScalingMat;
94   Vec            lambda_global;
95   IS             IS_l2g_lambda;
96   IS             subset,subset_mult,subset_n;
97   PetscBool      skip_node,fully_redundant;
98   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
99   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
100   PetscMPIInt    rank,size,buf_size,neigh;
101   PetscScalar    scalar_value;
102   PetscInt       *vertex_indices;
103   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
104   const PetscInt *aux_global_numbering;
105   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
106   PetscScalar    *array,*scaling_factors,*vals_B_delta;
107   PetscInt       *aux_local_numbering_2;
108   /* For communication of scaling factors */
109   PetscInt       *ptrs_buffer,neigh_position;
110   PetscScalar    **all_factors,*send_buffer,*recv_buffer;
111   MPI_Request    *send_reqs,*recv_reqs;
112   /* tests */
113   Vec            test_vec;
114   PetscBool      test_fetidp;
115   PetscViewer    viewer;
116 
117   PetscFunctionBegin;
118   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
119   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
120   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
121 
122   /* Default type of lagrange multipliers is non-redundant */
123   fully_redundant = PETSC_FALSE;
124   ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_fullyredundant",&fully_redundant,NULL);CHKERRQ(ierr);
125 
126   /* Evaluate local and global number of lagrange multipliers */
127   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
128   n_local_lambda = 0;
129   partial_sum = 0;
130   n_boundary_dofs = 0;
131   s = 0;
132   /* Get Vertices used to define the BDDC */
133   n_vertices = pcbddc->n_vertices;
134   vertex_indices = pcbddc->local_primal_ref_node;
135   dual_size = pcis->n_B-n_vertices;
136   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
137   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
138   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
139 
140   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
141   for (i=0;i<pcis->n;i++){
142     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
143     if ( j > 0 ) {
144       n_boundary_dofs++;
145     }
146     skip_node = PETSC_FALSE;
147     if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */
148       skip_node = PETSC_TRUE;
149       s++;
150     }
151     if (j < 1) {
152       skip_node = PETSC_TRUE;
153     }
154     if ( !skip_node ) {
155       if (fully_redundant) {
156         /* fully redundant set of lagrange multipliers */
157         n_lambda_for_dof = (j*(j+1))/2;
158       } else {
159         n_lambda_for_dof = j;
160       }
161       n_local_lambda += j;
162       /* needed to evaluate global number of lagrange multipliers */
163       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
164       /* store some data needed */
165       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
166       aux_local_numbering_1[partial_sum] = i;
167       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
168       partial_sum++;
169     }
170   }
171   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
172 
173   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
174   ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
175   ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
176   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
177   fetidpmat_ctx->n_lambda = (PetscInt)PetscRealPart(scalar_value);
178 
179   /* compute global ordering of lagrange multipliers and associate l2g map */
180   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
181   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
182   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
183   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
184   ierr = PCBDDCSubsetNumbering(subset,subset_mult,&i,&subset_n);CHKERRQ(ierr);
185   ierr = ISDestroy(&subset);CHKERRQ(ierr);
186   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);
187 
188   /* init data for scaling factors exchange */
189   partial_sum = 0;
190   j = 0;
191   ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
192   ierr = PetscMalloc1(pcis->n_neigh-1,&send_reqs);CHKERRQ(ierr);
193   ierr = PetscMalloc1(pcis->n_neigh-1,&recv_reqs);CHKERRQ(ierr);
194   ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr);
195   ptrs_buffer[0]=0;
196   for (i=1;i<pcis->n_neigh;i++) {
197     partial_sum += pcis->n_shared[i];
198     ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
199   }
200   ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
201   ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
202   ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
203   for (i=0;i<pcis->n-1;i++) {
204     j = mat_graph->count[i];
205     all_factors[i+1]=all_factors[i]+j;
206   }
207   /* scatter B scaling to N vec */
208   ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
209   ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
210   /* communications */
211   ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
212   for (i=1;i<pcis->n_neigh;i++) {
213     for (j=0;j<pcis->n_shared[i];j++) {
214       send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
215     }
216     ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
217     ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
218     ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
219     ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
220   }
221   ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
222   ierr = MPI_Waitall((pcis->n_neigh-1),recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
223   /* put values in correct places */
224   for (i=1;i<pcis->n_neigh;i++) {
225     for (j=0;j<pcis->n_shared[i];j++) {
226       k = pcis->shared[i][j];
227       neigh_position = 0;
228       while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
229       all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
230     }
231   }
232   ierr = MPI_Waitall((pcis->n_neigh-1),send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
233   ierr = PetscFree(send_reqs);CHKERRQ(ierr);
234   ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
235   ierr = PetscFree(send_buffer);CHKERRQ(ierr);
236   ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
237   ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
238 
239   /* Compute B and B_delta (local actions) */
240   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
241   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
242   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
243   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
244   ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
245   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
246   n_global_lambda=0;
247   partial_sum=0;
248   cum = 0;
249   for (i=0;i<dual_size;i++) {
250     n_global_lambda = aux_global_numbering[cum];
251     j = mat_graph->count[aux_local_numbering_1[i]];
252     aux_sums[0]=0;
253     for (s=1;s<j;s++) {
254       aux_sums[s]=aux_sums[s-1]+j-s+1;
255     }
256     array = all_factors[aux_local_numbering_1[i]];
257     n_neg_values = 0;
258     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
259     n_pos_values = j - n_neg_values;
260     if (fully_redundant) {
261       for (s=0;s<n_neg_values;s++) {
262         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
263         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
264         vals_B_delta   [partial_sum+s]=-1.0;
265         scaling_factors[partial_sum+s]=array[s];
266       }
267       for (s=0;s<n_pos_values;s++) {
268         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
269         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
270         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
271         scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
272       }
273       partial_sum += j;
274     } else {
275       /* l2g_indices and default cols and vals of B_delta */
276       for (s=0;s<j;s++) {
277         l2g_indices    [partial_sum+s]=n_global_lambda+s;
278         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
279         vals_B_delta   [partial_sum+s]=0.0;
280       }
281       /* B_delta */
282       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
283         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
284       }
285       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
286         vals_B_delta   [partial_sum+n_neg_values]=1.0;
287       }
288       /* scaling as in Klawonn-Widlund 1999*/
289       for (s=0;s<n_neg_values;s++) {
290         scalar_value = 0.0;
291         for (k=0;k<s+1;k++) {
292           scalar_value += array[k];
293         }
294         scaling_factors[partial_sum+s] = -scalar_value;
295       }
296       for (s=0;s<n_pos_values;s++) {
297         scalar_value = 0.0;
298         for (k=s+n_neg_values;k<j;k++) {
299           scalar_value += array[k];
300         }
301         scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
302       }
303       partial_sum += j;
304     }
305     cum += aux_local_numbering_2[i];
306   }
307   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
308   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
309   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
310   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
311   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
312   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
313   ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
314   ierr = PetscFree(all_factors);CHKERRQ(ierr);
315 
316   /* Local to global mapping of fetidpmat */
317   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
318   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
319   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
320   ierr = VecCreate(comm,&lambda_global);CHKERRQ(ierr);
321   ierr = VecSetSizes(lambda_global,PETSC_DECIDE,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
322   ierr = VecSetType(lambda_global,VECMPI);CHKERRQ(ierr);
323   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
324   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,(IS)0,lambda_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
325   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
326 
327   /* Create local part of B_delta */
328   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
329   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
330   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
331   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
332   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
333   for (i=0;i<n_local_lambda;i++) {
334     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
335   }
336   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
337   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
338   ierr = MatAssemblyEnd  (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
339 
340   if (fully_redundant) {
341     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
342     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
343     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
344     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
345     for (i=0;i<n_local_lambda;i++) {
346       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
347     }
348     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
349     ierr = MatAssemblyEnd  (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
350     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
351     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
352   } else {
353     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
354     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
355     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
356     ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
357     for (i=0;i<n_local_lambda;i++) {
358       ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
359     }
360     ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
361     ierr = MatAssemblyEnd  (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
362   }
363   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
364   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
365 
366   /* Create some vectors needed by fetidp */
367   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
368   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
369 
370   test_fetidp = PETSC_FALSE;
371   ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr);
372 
373   if (test_fetidp && !pcbddc->use_deluxe_scaling) {
374 
375     PetscReal real_value;
376 
377     ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
378     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
379     ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr);
380     ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr);
381     ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr);
382     if (fully_redundant) {
383       ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr);
384     } else {
385       ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr);
386     }
387     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
388 
389     /******************************************************************/
390     /* TEST A/B: Test numbering of global lambda dofs             */
391     /******************************************************************/
392 
393     ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr);
394     ierr = VecSet(lambda_global,1.0);CHKERRQ(ierr);
395     ierr = VecSet(test_vec,1.0);CHKERRQ(ierr);
396     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
397     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
398     scalar_value = -1.0;
399     ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
400     ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr);
401     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
402     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
403     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
404     if (fully_redundant) {
405       ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
406       ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr);
407       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
408       ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
409       ierr = VecSum(lambda_global,&scalar_value);CHKERRQ(ierr);
410       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
411       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
412     }
413 
414     /******************************************************************/
415     /* TEST C: It should holds B_delta*w=0, w\in\widehat{W}           */
416     /* This is the meaning of the B matrix                            */
417     /******************************************************************/
418 
419     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
420     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
421     ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
422     ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
423     ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
424     ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
425     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
426     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
427     /* Action of B_delta */
428     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
429     ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
430     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
431     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
432     ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
433     ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr);
434     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
435 
436     /******************************************************************/
437     /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W}     */
438     /* E_D = R_D^TR                                                   */
439     /* P_D = B_{D,delta}^T B_{delta}                                  */
440     /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
441     /******************************************************************/
442 
443     /* compute a random vector in \widetilde{W} */
444     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
445     scalar_value = 0.0;  /* set zero at vertices */
446     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
447     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
448     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
449     /* store w for final comparison */
450     ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr);
451     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
452     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
453 
454     /* Jump operator P_D : results stored in pcis->vec1_B */
455 
456     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
457     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
458     /* Action of B_delta */
459     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
460     ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
461     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
462     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
463     /* Action of B_Ddelta^T */
464     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
465     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
466     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
467 
468     /* Average operator E_D : results stored in pcis->vec2_B */
469     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
470     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
471     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr);
472     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
473     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
474 
475     /* test E_D=I-P_D */
476     scalar_value = 1.0;
477     ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr);
478     scalar_value = -1.0;
479     ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr);
480     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr);
481     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
482     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
483     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
484 
485     /******************************************************************/
486     /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W}          */
487     /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
488     /******************************************************************/
489 
490     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
491     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
492     scalar_value = 0.0;  /* set zero at vertices */
493     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
494     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
495 
496     /* Jump operator P_D : results stored in pcis->vec1_B */
497 
498     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
499     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
500     /* Action of B_delta */
501     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
502     ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
503     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
504     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
505     /* Action of B_Ddelta^T */
506     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
507     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
508     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
509     /* scaling */
510     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
511     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
512     ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr);
513     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
514 
515     if (!fully_redundant) {
516       /******************************************************************/
517       /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
518       /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
519       /******************************************************************/
520       ierr = VecDuplicate(lambda_global,&test_vec);CHKERRQ(ierr);
521       ierr = VecSetRandom(lambda_global,NULL);CHKERRQ(ierr);
522       /* Action of B_Ddelta^T */
523       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
524       ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
525       ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
526       /* Action of B_delta */
527       ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
528       ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
529       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
530       ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
531       scalar_value = -1.0;
532       ierr = VecAXPY(lambda_global,scalar_value,test_vec);CHKERRQ(ierr);
533       ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
534       ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr);
535       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
536       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
537       ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
538     }
539   }
540   /* final cleanup */
541   ierr = VecDestroy(&lambda_global);CHKERRQ(ierr);
542 
543   PetscFunctionReturn(0);
544 }
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "PCBDDCSetupFETIDPPCContext"
548 PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
549 {
550   FETIDPMat_ctx  mat_ctx;
551   PC_IS          *pcis;
552   PetscErrorCode ierr;
553 
554   PetscFunctionBegin;
555   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
556   /* get references from objects created when setting up feti mat context */
557   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
558   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
559   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
560   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
561   ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
562   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
563   /* create local Schur complement matrix */
564   pcis = (PC_IS*)fetidppc_ctx->pc->data;
565   ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
566   ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
567   PetscFunctionReturn(0);
568 }
569 
570 #undef __FUNCT__
571 #define __FUNCT__ "FETIDPMatMult"
572 PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
573 {
574   FETIDPMat_ctx  mat_ctx;
575   PC_IS          *pcis;
576   PetscErrorCode ierr;
577 
578   PetscFunctionBegin;
579   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
580   pcis = (PC_IS*)mat_ctx->pc->data;
581   /* Application of B_delta^T */
582   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
583   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
584   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
585   /* Application of \widetilde{S}^-1 */
586   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
587   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
588   /* Application of B_delta */
589   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
590   ierr = VecSet(y,0.0);CHKERRQ(ierr);
591   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
592   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595 
596 #undef __FUNCT__
597 #define __FUNCT__ "FETIDPMatMultTranspose"
598 PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
599 {
600   FETIDPMat_ctx  mat_ctx;
601   PC_IS          *pcis;
602   PetscErrorCode ierr;
603 
604   PetscFunctionBegin;
605   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
606   pcis = (PC_IS*)mat_ctx->pc->data;
607   /* Application of B_delta^T */
608   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
609   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
610   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
611   /* Application of \widetilde{S}^-1 */
612   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
613   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr);
614   /* Application of B_delta */
615   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
616   ierr = VecSet(y,0.0);CHKERRQ(ierr);
617   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
618   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
619   PetscFunctionReturn(0);
620 }
621 
622 #undef __FUNCT__
623 #define __FUNCT__ "FETIDPPCApply"
624 PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
625 {
626   FETIDPPC_ctx   pc_ctx;
627   PC_IS          *pcis;
628   PetscErrorCode ierr;
629 
630   PetscFunctionBegin;
631   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
632   pcis = (PC_IS*)pc_ctx->pc->data;
633   /* Application of B_Ddelta^T */
634   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
635   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
636   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
637   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
638   /* Application of local Schur complement */
639   ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
640   /* Application of B_Ddelta */
641   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
642   ierr = VecSet(y,0.0);CHKERRQ(ierr);
643   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
644   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
645   PetscFunctionReturn(0);
646 }
647 
648 #undef __FUNCT__
649 #define __FUNCT__ "FETIDPPCApplyTranspose"
650 PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
651 {
652   FETIDPPC_ctx   pc_ctx;
653   PC_IS          *pcis;
654   PetscErrorCode ierr;
655 
656   PetscFunctionBegin;
657   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
658   pcis = (PC_IS*)pc_ctx->pc->data;
659   /* Application of B_Ddelta^T */
660   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
661   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
662   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
663   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
664   /* Application of local Schur complement */
665   ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
666   /* Application of B_Ddelta */
667   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
668   ierr = VecSet(y,0.0);CHKERRQ(ierr);
669   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
670   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
671   PetscFunctionReturn(0);
672 }
673