xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 3bbff08aae55fef3c6bc0ec24cdbd560153bb73a)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3674ae819SStefano Zampini 
4674ae819SStefano Zampini #undef __FUNCT__
5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPMatContext"
6674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
7674ae819SStefano Zampini {
8674ae819SStefano Zampini   FETIDPMat_ctx  newctx;
9674ae819SStefano Zampini   PetscErrorCode ierr;
10674ae819SStefano Zampini 
11674ae819SStefano Zampini   PetscFunctionBegin;
12854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
13674ae819SStefano Zampini   newctx->lambda_local    = 0;
14674ae819SStefano Zampini   newctx->temp_solution_B = 0;
15674ae819SStefano Zampini   newctx->temp_solution_D = 0;
16674ae819SStefano Zampini   newctx->B_delta         = 0;
17674ae819SStefano Zampini   newctx->B_Ddelta        = 0; /* theoretically belongs to the FETIDP preconditioner */
18674ae819SStefano Zampini   newctx->l2g_lambda      = 0;
19674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
20674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
21674ae819SStefano Zampini   newctx->pc              = pc;
22674ae819SStefano Zampini   *fetidpmat_ctx          = newctx;
23674ae819SStefano Zampini   PetscFunctionReturn(0);
24674ae819SStefano Zampini }
25674ae819SStefano Zampini 
26674ae819SStefano Zampini #undef __FUNCT__
27674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext"
28674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
29674ae819SStefano Zampini {
30674ae819SStefano Zampini   FETIDPPC_ctx   newctx;
31674ae819SStefano Zampini   PetscErrorCode ierr;
32674ae819SStefano Zampini 
33674ae819SStefano Zampini   PetscFunctionBegin;
34854ce69bSBarry Smith   ierr = PetscNew(&newctx);CHKERRQ(ierr);
35674ae819SStefano Zampini   newctx->lambda_local    = 0;
36674ae819SStefano Zampini   newctx->B_Ddelta        = 0;
37674ae819SStefano Zampini   newctx->l2g_lambda      = 0;
38674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
39674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr);
40674ae819SStefano Zampini   newctx->pc              = pc;
41674ae819SStefano Zampini   *fetidppc_ctx           = newctx;
42674ae819SStefano Zampini   PetscFunctionReturn(0);
43674ae819SStefano Zampini }
44674ae819SStefano Zampini 
45674ae819SStefano Zampini #undef __FUNCT__
46674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat"
47674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
48674ae819SStefano Zampini {
49674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
50674ae819SStefano Zampini   PetscErrorCode ierr;
51674ae819SStefano Zampini 
52674ae819SStefano Zampini   PetscFunctionBegin;
53674ae819SStefano Zampini   ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr);
54674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr);
55674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr);
56674ae819SStefano Zampini   ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr);
57674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr);
58674ae819SStefano Zampini   ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr);
59674ae819SStefano Zampini   ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr);
60674ae819SStefano Zampini   ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
61674ae819SStefano Zampini   ierr = PetscFree(mat_ctx);CHKERRQ(ierr);
62674ae819SStefano Zampini   PetscFunctionReturn(0);
63674ae819SStefano Zampini }
64674ae819SStefano Zampini 
65674ae819SStefano Zampini #undef __FUNCT__
66674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC"
67674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
68674ae819SStefano Zampini {
69674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
70674ae819SStefano Zampini   PetscErrorCode ierr;
71674ae819SStefano Zampini 
72674ae819SStefano Zampini   PetscFunctionBegin;
73674ae819SStefano Zampini   ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr);
74674ae819SStefano Zampini   ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr);
75674ae819SStefano Zampini   ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr);
76674ae819SStefano Zampini   ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr);
7768668abeSStefano Zampini   ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr);
78674ae819SStefano Zampini   ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */
79674ae819SStefano Zampini   ierr = PetscFree(pc_ctx);CHKERRQ(ierr);
80674ae819SStefano Zampini   PetscFunctionReturn(0);
81674ae819SStefano Zampini }
82674ae819SStefano Zampini 
83674ae819SStefano Zampini #undef __FUNCT__
84674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext"
85674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
86674ae819SStefano Zampini {
87674ae819SStefano Zampini   PetscErrorCode ierr;
88674ae819SStefano Zampini   PC_IS          *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
89674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
90674ae819SStefano Zampini   PCBDDCGraph    mat_graph=pcbddc->mat_graph;
91674ae819SStefano Zampini   Mat_IS         *matis  = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
92674ae819SStefano Zampini   MPI_Comm       comm;
93674ae819SStefano Zampini   Mat            ScalingMat;
94674ae819SStefano Zampini   Vec            lambda_global;
95674ae819SStefano Zampini   IS             IS_l2g_lambda;
96dc456d91SStefano Zampini   IS             subset,subset_mult,subset_n;
97674ae819SStefano Zampini   PetscBool      skip_node,fully_redundant;
98674ae819SStefano Zampini   PetscInt       i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
99dc456d91SStefano Zampini   PetscInt       cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
10076ec1555SBarry Smith   PetscMPIInt    rank,size,buf_size,neigh;
101674ae819SStefano Zampini   PetscScalar    scalar_value;
102674ae819SStefano Zampini   PetscInt       *vertex_indices;
103dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices,*aux_local_numbering_1;
104dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
105674ae819SStefano Zampini   PetscInt       *aux_sums,*cols_B_delta,*l2g_indices;
106674ae819SStefano Zampini   PetscScalar    *array,*scaling_factors,*vals_B_delta;
107674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
108674ae819SStefano Zampini   /* For communication of scaling factors */
109674ae819SStefano Zampini   PetscInt       *ptrs_buffer,neigh_position;
110674ae819SStefano Zampini   PetscScalar    **all_factors,*send_buffer,*recv_buffer;
111674ae819SStefano Zampini   MPI_Request    *send_reqs,*recv_reqs;
112674ae819SStefano Zampini   /* tests */
113674ae819SStefano Zampini   Vec            test_vec;
114674ae819SStefano Zampini   PetscBool      test_fetidp;
115674ae819SStefano Zampini   PetscViewer    viewer;
116674ae819SStefano Zampini 
117674ae819SStefano Zampini   PetscFunctionBegin;
118674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr);
119674ae819SStefano Zampini   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12076ec1555SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
121674ae819SStefano Zampini 
122674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
123674ae819SStefano Zampini   fully_redundant = PETSC_FALSE;
124674ae819SStefano Zampini   ierr = PetscOptionsGetBool(NULL,"-fetidp_fullyredundant",&fully_redundant,NULL);CHKERRQ(ierr);
125674ae819SStefano Zampini 
126674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
127674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr);
128674ae819SStefano Zampini   n_local_lambda = 0;
129674ae819SStefano Zampini   partial_sum = 0;
130674ae819SStefano Zampini   n_boundary_dofs = 0;
131674ae819SStefano Zampini   s = 0;
132674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
133b371cd4fSStefano Zampini   n_vertices = pcbddc->n_vertices;
1340e6343abSStefano Zampini   vertex_indices = pcbddc->local_primal_ref_node;
135674ae819SStefano Zampini   dual_size = pcis->n_B-n_vertices;
136785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr);
137785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr);
138785e854fSJed Brown   ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr);
139674ae819SStefano Zampini 
140674ae819SStefano Zampini   ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
141674ae819SStefano Zampini   for (i=0;i<pcis->n;i++){
142674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
143674ae819SStefano Zampini     if ( j > 0 ) {
144674ae819SStefano Zampini       n_boundary_dofs++;
145674ae819SStefano Zampini     }
146674ae819SStefano Zampini     skip_node = PETSC_FALSE;
147674ae819SStefano Zampini     if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */
148674ae819SStefano Zampini       skip_node = PETSC_TRUE;
149674ae819SStefano Zampini       s++;
150674ae819SStefano Zampini     }
151674ae819SStefano Zampini     if (j < 1) {
152674ae819SStefano Zampini       skip_node = PETSC_TRUE;
153674ae819SStefano Zampini     }
154674ae819SStefano Zampini     if ( !skip_node ) {
155674ae819SStefano Zampini       if (fully_redundant) {
156674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
157674ae819SStefano Zampini         n_lambda_for_dof = (j*(j+1))/2;
158674ae819SStefano Zampini       } else {
159674ae819SStefano Zampini         n_lambda_for_dof = j;
160674ae819SStefano Zampini       }
161674ae819SStefano Zampini       n_local_lambda += j;
162674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
163674ae819SStefano Zampini       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
164674ae819SStefano Zampini       /* store some data needed */
165674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
166674ae819SStefano Zampini       aux_local_numbering_1[partial_sum] = i;
167674ae819SStefano Zampini       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
168674ae819SStefano Zampini       partial_sum++;
169674ae819SStefano Zampini     }
170674ae819SStefano Zampini   }
171674ae819SStefano Zampini   ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
172674ae819SStefano Zampini 
173674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
174674ae819SStefano Zampini   ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
175674ae819SStefano Zampini   ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
176674ae819SStefano Zampini   ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr);
1775b08dc53SStefano Zampini   fetidpmat_ctx->n_lambda = (PetscInt)PetscRealPart(scalar_value);
178674ae819SStefano Zampini 
179674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
180dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
181*3bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr);
182dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
183dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr);
184dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(subset,subset_mult,&i,&subset_n);CHKERRQ(ierr);
185dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
186674ae819SStefano Zampini   if (i != fetidpmat_ctx->n_lambda) {
187674ae819SStefano Zampini     SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in %s: global number of multipliers mismatch! (%d!=%d)\n",__FUNCT__,fetidpmat_ctx->n_lambda,i);
188674ae819SStefano Zampini   }
189674ae819SStefano Zampini 
190674ae819SStefano Zampini   /* init data for scaling factors exchange */
191674ae819SStefano Zampini   partial_sum = 0;
192674ae819SStefano Zampini   j = 0;
193785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr);
194854ce69bSBarry Smith   ierr = PetscMalloc1(pcis->n_neigh-1,&send_reqs);CHKERRQ(ierr);
195854ce69bSBarry Smith   ierr = PetscMalloc1(pcis->n_neigh-1,&recv_reqs);CHKERRQ(ierr);
196785e854fSJed Brown   ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr);
197674ae819SStefano Zampini   ptrs_buffer[0]=0;
198674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
199674ae819SStefano Zampini     partial_sum += pcis->n_shared[i];
200674ae819SStefano Zampini     ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
201674ae819SStefano Zampini   }
202785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr);
203785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr);
204785e854fSJed Brown   ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr);
205674ae819SStefano Zampini   for (i=0;i<pcis->n-1;i++) {
206674ae819SStefano Zampini     j = mat_graph->count[i];
207674ae819SStefano Zampini     all_factors[i+1]=all_factors[i]+j;
208674ae819SStefano Zampini   }
209674ae819SStefano Zampini   /* scatter B scaling to N vec */
210674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
211674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
212674ae819SStefano Zampini   /* communications */
2132b095fd8SStefano Zampini   ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
214674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
215674ae819SStefano Zampini     for (j=0;j<pcis->n_shared[i];j++) {
216674ae819SStefano Zampini       send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
217674ae819SStefano Zampini     }
218674ae819SStefano Zampini     ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr);
219674ae819SStefano Zampini     ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr);
220674ae819SStefano Zampini     ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr);
221674ae819SStefano Zampini     ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr);
222674ae819SStefano Zampini   }
2232b095fd8SStefano Zampini   ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr);
224674ae819SStefano Zampini   ierr = MPI_Waitall((pcis->n_neigh-1),recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
225674ae819SStefano Zampini   /* put values in correct places */
226674ae819SStefano Zampini   for (i=1;i<pcis->n_neigh;i++) {
227674ae819SStefano Zampini     for (j=0;j<pcis->n_shared[i];j++) {
228674ae819SStefano Zampini       k = pcis->shared[i][j];
229674ae819SStefano Zampini       neigh_position = 0;
230674ae819SStefano Zampini       while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
231674ae819SStefano Zampini       all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
232674ae819SStefano Zampini     }
233674ae819SStefano Zampini   }
234674ae819SStefano Zampini   ierr = MPI_Waitall((pcis->n_neigh-1),send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
235674ae819SStefano Zampini   ierr = PetscFree(send_reqs);CHKERRQ(ierr);
236674ae819SStefano Zampini   ierr = PetscFree(recv_reqs);CHKERRQ(ierr);
237674ae819SStefano Zampini   ierr = PetscFree(send_buffer);CHKERRQ(ierr);
238674ae819SStefano Zampini   ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
239674ae819SStefano Zampini   ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr);
240674ae819SStefano Zampini 
241674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
242785e854fSJed Brown   ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr);
243785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr);
244785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr);
245785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr);
246785e854fSJed Brown   ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr);
247dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
248674ae819SStefano Zampini   n_global_lambda=0;
249674ae819SStefano Zampini   partial_sum=0;
250dc456d91SStefano Zampini   cum = 0;
251674ae819SStefano Zampini   for (i=0;i<dual_size;i++) {
252dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
253674ae819SStefano Zampini     j = mat_graph->count[aux_local_numbering_1[i]];
254674ae819SStefano Zampini     aux_sums[0]=0;
255674ae819SStefano Zampini     for (s=1;s<j;s++) {
256674ae819SStefano Zampini       aux_sums[s]=aux_sums[s-1]+j-s+1;
257674ae819SStefano Zampini     }
258674ae819SStefano Zampini     array = all_factors[aux_local_numbering_1[i]];
259674ae819SStefano Zampini     n_neg_values = 0;
2602a7da448SStefano Zampini     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
261674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
262674ae819SStefano Zampini     if (fully_redundant) {
263674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
264674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
265674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
266674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=-1.0;
267674ae819SStefano Zampini         scaling_factors[partial_sum+s]=array[s];
268674ae819SStefano Zampini       }
269674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
270674ae819SStefano Zampini         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
271674ae819SStefano Zampini         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
272674ae819SStefano Zampini         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
273674ae819SStefano Zampini         scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
274674ae819SStefano Zampini       }
275674ae819SStefano Zampini       partial_sum += j;
276674ae819SStefano Zampini     } else {
277674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
278674ae819SStefano Zampini       for (s=0;s<j;s++) {
279674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=n_global_lambda+s;
280674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
281674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=0.0;
282674ae819SStefano Zampini       }
283674ae819SStefano Zampini       /* B_delta */
284674ae819SStefano Zampini       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
285674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
286674ae819SStefano Zampini       }
287674ae819SStefano Zampini       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
288674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values]=1.0;
289674ae819SStefano Zampini       }
290674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999*/
291674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
292674ae819SStefano Zampini         scalar_value = 0.0;
293674ae819SStefano Zampini         for (k=0;k<s+1;k++) {
294674ae819SStefano Zampini           scalar_value += array[k];
295674ae819SStefano Zampini         }
296674ae819SStefano Zampini         scaling_factors[partial_sum+s] = -scalar_value;
297674ae819SStefano Zampini       }
298674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
299674ae819SStefano Zampini         scalar_value = 0.0;
300674ae819SStefano Zampini         for (k=s+n_neg_values;k<j;k++) {
301674ae819SStefano Zampini           scalar_value += array[k];
302674ae819SStefano Zampini         }
303674ae819SStefano Zampini         scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
304674ae819SStefano Zampini       }
305674ae819SStefano Zampini       partial_sum += j;
306674ae819SStefano Zampini     }
307dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
308674ae819SStefano Zampini   }
309dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr);
310dc456d91SStefano Zampini   ierr = ISDestroy(&subset_mult);CHKERRQ(ierr);
311dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
312674ae819SStefano Zampini   ierr = PetscFree(aux_sums);CHKERRQ(ierr);
313674ae819SStefano Zampini   ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr);
314674ae819SStefano Zampini   ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr);
315674ae819SStefano Zampini   ierr = PetscFree(all_factors[0]);CHKERRQ(ierr);
316674ae819SStefano Zampini   ierr = PetscFree(all_factors);CHKERRQ(ierr);
317674ae819SStefano Zampini 
318674ae819SStefano Zampini   /* Local to global mapping of fetidpmat */
319674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
320674ae819SStefano Zampini   ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
321674ae819SStefano Zampini   ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr);
322674ae819SStefano Zampini   ierr = VecCreate(comm,&lambda_global);CHKERRQ(ierr);
323674ae819SStefano Zampini   ierr = VecSetSizes(lambda_global,PETSC_DECIDE,fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
324674ae819SStefano Zampini   ierr = VecSetType(lambda_global,VECMPI);CHKERRQ(ierr);
325674ae819SStefano Zampini   ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr);
326674ae819SStefano Zampini   ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,(IS)0,lambda_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr);
327674ae819SStefano Zampini   ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr);
328674ae819SStefano Zampini 
329674ae819SStefano Zampini   /* Create local part of B_delta */
330302440fdSBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr);
331674ae819SStefano Zampini   ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
332674ae819SStefano Zampini   ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr);
333674ae819SStefano Zampini   ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr);
334674ae819SStefano Zampini   ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
335674ae819SStefano Zampini   for (i=0;i<n_local_lambda;i++) {
336674ae819SStefano Zampini     ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr);
337674ae819SStefano Zampini   }
338674ae819SStefano Zampini   ierr = PetscFree(vals_B_delta);CHKERRQ(ierr);
339674ae819SStefano Zampini   ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
340674ae819SStefano Zampini   ierr = MatAssemblyEnd  (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
341674ae819SStefano Zampini 
342674ae819SStefano Zampini   if (fully_redundant) {
343302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr);
344674ae819SStefano Zampini     ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr);
345674ae819SStefano Zampini     ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr);
346674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr);
347674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
348674ae819SStefano Zampini       ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
349674ae819SStefano Zampini     }
350674ae819SStefano Zampini     ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
351674ae819SStefano Zampini     ierr = MatAssemblyEnd  (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
352674ae819SStefano Zampini     ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
353674ae819SStefano Zampini     ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr);
354674ae819SStefano Zampini   } else {
355302440fdSBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr);
356674ae819SStefano Zampini     ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr);
357674ae819SStefano Zampini     ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr);
358674ae819SStefano Zampini     ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr);
359674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
360674ae819SStefano Zampini       ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr);
361674ae819SStefano Zampini     }
362674ae819SStefano Zampini     ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
363674ae819SStefano Zampini     ierr = MatAssemblyEnd  (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
364674ae819SStefano Zampini   }
365674ae819SStefano Zampini   ierr = PetscFree(scaling_factors);CHKERRQ(ierr);
366674ae819SStefano Zampini   ierr = PetscFree(cols_B_delta);CHKERRQ(ierr);
367674ae819SStefano Zampini 
368674ae819SStefano Zampini   /* Create some vectors needed by fetidp */
369674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr);
370674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr);
371674ae819SStefano Zampini 
372674ae819SStefano Zampini   test_fetidp = PETSC_FALSE;
373674ae819SStefano Zampini   ierr = PetscOptionsGetBool(NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr);
374674ae819SStefano Zampini 
375674ae819SStefano Zampini   if (test_fetidp && !pcbddc->use_deluxe_scaling) {
376674ae819SStefano Zampini 
3775b08dc53SStefano Zampini     PetscReal real_value;
3785b08dc53SStefano Zampini 
379674ae819SStefano Zampini     ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
380674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
381674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr);
382674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr);
383674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr);
384674ae819SStefano Zampini     if (fully_redundant) {
385674ae819SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr);
386674ae819SStefano Zampini     } else {
387674ae819SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr);
388674ae819SStefano Zampini     }
389674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
390674ae819SStefano Zampini 
391674ae819SStefano Zampini     /******************************************************************/
392674ae819SStefano Zampini     /* TEST A/B: Test numbering of global lambda dofs             */
393674ae819SStefano Zampini     /******************************************************************/
394674ae819SStefano Zampini 
395674ae819SStefano Zampini     ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr);
396674ae819SStefano Zampini     ierr = VecSet(lambda_global,1.0);CHKERRQ(ierr);
397674ae819SStefano Zampini     ierr = VecSet(test_vec,1.0);CHKERRQ(ierr);
398674ae819SStefano Zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
399674ae819SStefano Zampini     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
400674ae819SStefano Zampini     scalar_value = -1.0;
401674ae819SStefano Zampini     ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
4025b08dc53SStefano Zampini     ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr);
403674ae819SStefano Zampini     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
4045b08dc53SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
405674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
406674ae819SStefano Zampini     if (fully_redundant) {
407674ae819SStefano Zampini       ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
408674ae819SStefano Zampini       ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr);
409674ae819SStefano Zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
410674ae819SStefano Zampini       ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
411674ae819SStefano Zampini       ierr = VecSum(lambda_global,&scalar_value);CHKERRQ(ierr);
4125b08dc53SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr);
413674ae819SStefano Zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
414674ae819SStefano Zampini     }
415674ae819SStefano Zampini 
416674ae819SStefano Zampini     /******************************************************************/
417674ae819SStefano Zampini     /* TEST C: It should holds B_delta*w=0, w\in\widehat{W}           */
418674ae819SStefano Zampini     /* This is the meaning of the B matrix                            */
419674ae819SStefano Zampini     /******************************************************************/
420674ae819SStefano Zampini 
421674ae819SStefano Zampini     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
422674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
423674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
424674ae819SStefano Zampini     ierr = VecScatterEnd  (matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
425674ae819SStefano Zampini     ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
426674ae819SStefano Zampini     ierr = VecScatterEnd  (matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
427674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
428674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
429674ae819SStefano Zampini     /* Action of B_delta */
430674ae819SStefano Zampini     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
431674ae819SStefano Zampini     ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
432674ae819SStefano Zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
433674ae819SStefano Zampini     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4345b08dc53SStefano Zampini     ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
4355b08dc53SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr);
436674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
437674ae819SStefano Zampini 
438674ae819SStefano Zampini     /******************************************************************/
439674ae819SStefano Zampini     /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W}     */
440674ae819SStefano Zampini     /* E_D = R_D^TR                                                   */
441674ae819SStefano Zampini     /* P_D = B_{D,delta}^T B_{delta}                                  */
442674ae819SStefano Zampini     /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
443674ae819SStefano Zampini     /******************************************************************/
444674ae819SStefano Zampini 
445674ae819SStefano Zampini     /* compute a random vector in \widetilde{W} */
446674ae819SStefano Zampini     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
447674ae819SStefano Zampini     scalar_value = 0.0;  /* set zero at vertices */
448674ae819SStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
449674ae819SStefano Zampini     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
450674ae819SStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
451674ae819SStefano Zampini     /* store w for final comparison */
452674ae819SStefano Zampini     ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr);
453674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
454674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
455674ae819SStefano Zampini 
456674ae819SStefano Zampini     /* Jump operator P_D : results stored in pcis->vec1_B */
457674ae819SStefano Zampini 
458674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
459674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
460674ae819SStefano Zampini     /* Action of B_delta */
461674ae819SStefano Zampini     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
462674ae819SStefano Zampini     ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
463674ae819SStefano Zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
464674ae819SStefano Zampini     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
465674ae819SStefano Zampini     /* Action of B_Ddelta^T */
466674ae819SStefano Zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
467674ae819SStefano Zampini     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
468674ae819SStefano Zampini     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
469674ae819SStefano Zampini 
470674ae819SStefano Zampini     /* Average operator E_D : results stored in pcis->vec2_B */
471674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
472674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
473674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr);
474674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
475674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
476674ae819SStefano Zampini 
477674ae819SStefano Zampini     /* test E_D=I-P_D */
478674ae819SStefano Zampini     scalar_value = 1.0;
479674ae819SStefano Zampini     ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr);
480674ae819SStefano Zampini     scalar_value = -1.0;
481674ae819SStefano Zampini     ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr);
4825b08dc53SStefano Zampini     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr);
483674ae819SStefano Zampini     ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
4845b08dc53SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr);
485674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
486674ae819SStefano Zampini 
487674ae819SStefano Zampini     /******************************************************************/
488674ae819SStefano Zampini     /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W}          */
489674ae819SStefano Zampini     /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
490674ae819SStefano Zampini     /******************************************************************/
491674ae819SStefano Zampini 
492674ae819SStefano Zampini     ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr);
493674ae819SStefano Zampini     ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr);
494674ae819SStefano Zampini     scalar_value = 0.0;  /* set zero at vertices */
495674ae819SStefano Zampini     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
496674ae819SStefano Zampini     ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr);
497674ae819SStefano Zampini 
498674ae819SStefano Zampini     /* Jump operator P_D : results stored in pcis->vec1_B */
499674ae819SStefano Zampini 
500674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
501674ae819SStefano Zampini     ierr = VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
502674ae819SStefano Zampini     /* Action of B_delta */
503674ae819SStefano Zampini     ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
504674ae819SStefano Zampini     ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr);
505674ae819SStefano Zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
506674ae819SStefano Zampini     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
507674ae819SStefano Zampini     /* Action of B_Ddelta^T */
508674ae819SStefano Zampini     ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
509674ae819SStefano Zampini     ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
510674ae819SStefano Zampini     ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
511674ae819SStefano Zampini     /* scaling */
512674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
5135b08dc53SStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
5145b08dc53SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr);
515674ae819SStefano Zampini     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
516674ae819SStefano Zampini 
517674ae819SStefano Zampini     if (!fully_redundant) {
518674ae819SStefano Zampini       /******************************************************************/
519674ae819SStefano Zampini       /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
520674ae819SStefano Zampini       /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
521674ae819SStefano Zampini       /******************************************************************/
522674ae819SStefano Zampini       ierr = VecDuplicate(lambda_global,&test_vec);CHKERRQ(ierr);
523674ae819SStefano Zampini       ierr = VecSetRandom(lambda_global,NULL);CHKERRQ(ierr);
524674ae819SStefano Zampini       /* Action of B_Ddelta^T */
525674ae819SStefano Zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
526674ae819SStefano Zampini       ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
527674ae819SStefano Zampini       ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
528674ae819SStefano Zampini       /* Action of B_delta */
529674ae819SStefano Zampini       ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr);
530674ae819SStefano Zampini       ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
531674ae819SStefano Zampini       ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
532674ae819SStefano Zampini       ierr = VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
533674ae819SStefano Zampini       scalar_value = -1.0;
534674ae819SStefano Zampini       ierr = VecAXPY(lambda_global,scalar_value,test_vec);CHKERRQ(ierr);
5355b08dc53SStefano Zampini       ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr);
5365b08dc53SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr);
537674ae819SStefano Zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
538674ae819SStefano Zampini       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
539674ae819SStefano Zampini       ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
540674ae819SStefano Zampini     }
541674ae819SStefano Zampini   }
542674ae819SStefano Zampini   /* final cleanup */
543674ae819SStefano Zampini   ierr = VecDestroy(&lambda_global);CHKERRQ(ierr);
544674ae819SStefano Zampini 
545674ae819SStefano Zampini   PetscFunctionReturn(0);
546674ae819SStefano Zampini }
547674ae819SStefano Zampini 
548674ae819SStefano Zampini #undef __FUNCT__
549674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext"
550674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
551674ae819SStefano Zampini {
552674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
553ed6c3d69SStefano Zampini   PC_IS          *pcis;
554674ae819SStefano Zampini   PetscErrorCode ierr;
555674ae819SStefano Zampini 
556674ae819SStefano Zampini   PetscFunctionBegin;
557674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
558674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
559674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr);
560674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
561674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr);
562674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
563674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr);
564674ae819SStefano Zampini   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
565ed6c3d69SStefano Zampini   /* create local Schur complement matrix */
566ed6c3d69SStefano Zampini   pcis = (PC_IS*)fetidppc_ctx->pc->data;
567ed6c3d69SStefano Zampini   ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr);
568ed6c3d69SStefano Zampini   ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr);
569674ae819SStefano Zampini   PetscFunctionReturn(0);
570674ae819SStefano Zampini }
571674ae819SStefano Zampini 
572674ae819SStefano Zampini #undef __FUNCT__
573674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult"
574674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
575674ae819SStefano Zampini {
576674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
577674ae819SStefano Zampini   PC_IS          *pcis;
578674ae819SStefano Zampini   PetscErrorCode ierr;
579674ae819SStefano Zampini 
580674ae819SStefano Zampini   PetscFunctionBegin;
581674ae819SStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
582674ae819SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
583674ae819SStefano Zampini   /* Application of B_delta^T */
584674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
585674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
586674ae819SStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
587674ae819SStefano Zampini   /* Application of \widetilde{S}^-1 */
588674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
589dc359a40SStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr);
590674ae819SStefano Zampini   /* Application of B_delta */
591674ae819SStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
592674ae819SStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
593674ae819SStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
594674ae819SStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
595674ae819SStefano Zampini   PetscFunctionReturn(0);
596674ae819SStefano Zampini }
597674ae819SStefano Zampini 
598674ae819SStefano Zampini #undef __FUNCT__
599edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose"
600edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
601edf7251bSStefano Zampini {
602edf7251bSStefano Zampini   FETIDPMat_ctx  mat_ctx;
603edf7251bSStefano Zampini   PC_IS          *pcis;
604edf7251bSStefano Zampini   PetscErrorCode ierr;
605edf7251bSStefano Zampini 
606edf7251bSStefano Zampini   PetscFunctionBegin;
607edf7251bSStefano Zampini   ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr);
608edf7251bSStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
609edf7251bSStefano Zampini   /* Application of B_delta^T */
610edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
611edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
612edf7251bSStefano Zampini   ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr);
613edf7251bSStefano Zampini   /* Application of \widetilde{S}^-1 */
614edf7251bSStefano Zampini   ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr);
615edf7251bSStefano Zampini   ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr);
616edf7251bSStefano Zampini   /* Application of B_delta */
617edf7251bSStefano Zampini   ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr);
618edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
619edf7251bSStefano Zampini   ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
620edf7251bSStefano Zampini   ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
621edf7251bSStefano Zampini   PetscFunctionReturn(0);
622edf7251bSStefano Zampini }
623edf7251bSStefano Zampini 
624edf7251bSStefano Zampini #undef __FUNCT__
625674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply"
626674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
627674ae819SStefano Zampini {
628674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
629674ae819SStefano Zampini   PC_IS          *pcis;
630674ae819SStefano Zampini   PetscErrorCode ierr;
631674ae819SStefano Zampini 
632674ae819SStefano Zampini   PetscFunctionBegin;
633302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
634674ae819SStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
635674ae819SStefano Zampini   /* Application of B_Ddelta^T */
636674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
637674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
638674ae819SStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
639674ae819SStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
640ed6c3d69SStefano Zampini   /* Application of local Schur complement */
641ed6c3d69SStefano Zampini   ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
642edf7251bSStefano Zampini   /* Application of B_Ddelta */
643edf7251bSStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
644edf7251bSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
645edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
646edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
647edf7251bSStefano Zampini   PetscFunctionReturn(0);
648edf7251bSStefano Zampini }
649edf7251bSStefano Zampini 
650edf7251bSStefano Zampini #undef __FUNCT__
651edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose"
652edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
653edf7251bSStefano Zampini {
654edf7251bSStefano Zampini   FETIDPPC_ctx   pc_ctx;
655edf7251bSStefano Zampini   PC_IS          *pcis;
656edf7251bSStefano Zampini   PetscErrorCode ierr;
657edf7251bSStefano Zampini 
658edf7251bSStefano Zampini   PetscFunctionBegin;
659302440fdSBarry Smith   ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr);
660edf7251bSStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
661edf7251bSStefano Zampini   /* Application of B_Ddelta^T */
662edf7251bSStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
663edf7251bSStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
664edf7251bSStefano Zampini   ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr);
665edf7251bSStefano Zampini   ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr);
666ed6c3d69SStefano Zampini   /* Application of local Schur complement */
667ed6c3d69SStefano Zampini   ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr);
668674ae819SStefano Zampini   /* Application of B_Ddelta */
669674ae819SStefano Zampini   ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr);
670674ae819SStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
671674ae819SStefano Zampini   ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
672674ae819SStefano Zampini   ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
673674ae819SStefano Zampini   PetscFunctionReturn(0);
674674ae819SStefano Zampini }
675