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