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