Actual source code: bddcfetidp.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: #include <../src/ksp/pc/impls/bddc/bddc.h>
  2: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>

  6: PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
  7: {
  8:   FETIDPMat_ctx  newctx;

 12:   PetscNew(&newctx);
 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:   PetscObjectReference((PetscObject)pc);
 21:   newctx->pc              = pc;
 22:   *fetidpmat_ctx          = newctx;
 23:   return(0);
 24: }

 28: PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
 29: {
 30:   FETIDPPC_ctx   newctx;

 34:   PetscNew(&newctx);
 35:   newctx->lambda_local    = 0;
 36:   newctx->B_Ddelta        = 0;
 37:   newctx->l2g_lambda      = 0;
 38:   /* increase the reference count for BDDC preconditioner */
 39:   PetscObjectReference((PetscObject)pc);
 40:   newctx->pc              = pc;
 41:   *fetidppc_ctx           = newctx;
 42:   return(0);
 43: }

 47: PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
 48: {
 49:   FETIDPMat_ctx  mat_ctx;

 53:   MatShellGetContext(A,(void**)&mat_ctx);
 54:   VecDestroy(&mat_ctx->lambda_local);
 55:   VecDestroy(&mat_ctx->temp_solution_D);
 56:   VecDestroy(&mat_ctx->temp_solution_B);
 57:   MatDestroy(&mat_ctx->B_delta);
 58:   MatDestroy(&mat_ctx->B_Ddelta);
 59:   VecScatterDestroy(&mat_ctx->l2g_lambda);
 60:   PCDestroy(&mat_ctx->pc); /* decrease PCBDDC reference count */
 61:   PetscFree(mat_ctx);
 62:   return(0);
 63: }

 67: PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
 68: {
 69:   FETIDPPC_ctx   pc_ctx;

 73:   PCShellGetContext(pc,(void**)&pc_ctx);
 74:   VecDestroy(&pc_ctx->lambda_local);
 75:   MatDestroy(&pc_ctx->B_Ddelta);
 76:   VecScatterDestroy(&pc_ctx->l2g_lambda);
 77:   MatDestroy(&pc_ctx->S_j);
 78:   PCDestroy(&pc_ctx->pc); /* decrease PCBDDC reference count */
 79:   PetscFree(pc_ctx);
 80:   return(0);
 81: }

 85: PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
 86: {
 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;

118:   PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);
119:   MPI_Comm_rank(comm,&rank);
120:   MPI_Comm_size(comm,&size);

122:   /* Default type of lagrange multipliers is non-redundant */
123:   fully_redundant = PETSC_FALSE;
124:   PetscOptionsGetBool(NULL,NULL,"-fetidp_fullyredundant",&fully_redundant,NULL);

126:   /* Evaluate local and global number of lagrange multipliers */
127:   VecSet(pcis->vec1_N,0.0);
128:   n_local_lambda = 0;
129:   partial_sum = 0;
130:   n_boundary_dofs = 0;
131:   s = 0;
132:   /* Get Vertices used to define the BDDC */
133:   n_vertices = pcbddc->n_vertices;
134:   vertex_indices = pcbddc->local_primal_ref_node;
135:   dual_size = pcis->n_B-n_vertices;
136:   PetscMalloc1(dual_size,&dual_dofs_boundary_indices);
137:   PetscMalloc1(dual_size,&aux_local_numbering_1);
138:   PetscMalloc1(dual_size,&aux_local_numbering_2);

140:   VecGetArray(pcis->vec1_N,&array);
141:   for (i=0;i<pcis->n;i++){
142:     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
143:     if ( j > 0 ) {
144:       n_boundary_dofs++;
145:     }
146:     skip_node = PETSC_FALSE;
147:     if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */
148:       skip_node = PETSC_TRUE;
149:       s++;
150:     }
151:     if (j < 1) {
152:       skip_node = PETSC_TRUE;
153:     }
154:     if ( !skip_node ) {
155:       if (fully_redundant) {
156:         /* fully redundant set of lagrange multipliers */
157:         n_lambda_for_dof = (j*(j+1))/2;
158:       } else {
159:         n_lambda_for_dof = j;
160:       }
161:       n_local_lambda += j;
162:       /* needed to evaluate global number of lagrange multipliers */
163:       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
164:       /* store some data needed */
165:       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
166:       aux_local_numbering_1[partial_sum] = i;
167:       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
168:       partial_sum++;
169:     }
170:   }
171:   VecRestoreArray(pcis->vec1_N,&array);

173:   VecSet(pcis->vec1_global,0.0);
174:   VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
175:   VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
176:   VecSum(pcis->vec1_global,&scalar_value);
177:   fetidpmat_ctx->n_lambda = (PetscInt)PetscRealPart(scalar_value);

179:   /* compute global ordering of lagrange multipliers and associate l2g map */
180:   ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);
181:   ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);
182:   ISDestroy(&subset_n);
183:   ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);
184:   PCBDDCSubsetNumbering(subset,subset_mult,&i,&subset_n);
185:   ISDestroy(&subset);
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);

188:   /* init data for scaling factors exchange */
189:   partial_sum = 0;
190:   PetscMalloc1(pcis->n_neigh,&ptrs_buffer);
191:   PetscMalloc1(pcis->n_neigh-1,&send_reqs);
192:   PetscMalloc1(pcis->n_neigh-1,&recv_reqs);
193:   PetscMalloc1(pcis->n,&all_factors);
194:   ptrs_buffer[0]=0;
195:   for (i=1;i<pcis->n_neigh;i++) {
196:     partial_sum += pcis->n_shared[i];
197:     ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
198:   }
199:   PetscMalloc1(partial_sum,&send_buffer);
200:   PetscMalloc1(partial_sum,&recv_buffer);
201:   PetscMalloc1(partial_sum,&all_factors[0]);
202:   for (i=0;i<pcis->n-1;i++) {
203:     j = mat_graph->count[i];
204:     all_factors[i+1]=all_factors[i]+j;
205:   }
206:   /* scatter B scaling to N vec */
207:   VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
208:   VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
209:   /* communications */
210:   VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);
211:   for (i=1;i<pcis->n_neigh;i++) {
212:     for (j=0;j<pcis->n_shared[i];j++) {
213:       send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
214:     }
215:     PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);
216:     PetscMPIIntCast(pcis->neigh[i],&neigh);
217:     MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);
218:     MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);
219:   }
220:   VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);
221:   MPI_Waitall((pcis->n_neigh-1),recv_reqs,MPI_STATUSES_IGNORE);
222:   /* put values in correct places */
223:   for (i=1;i<pcis->n_neigh;i++) {
224:     for (j=0;j<pcis->n_shared[i];j++) {
225:       k = pcis->shared[i][j];
226:       neigh_position = 0;
227:       while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
228:       all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
229:     }
230:   }
231:   MPI_Waitall((pcis->n_neigh-1),send_reqs,MPI_STATUSES_IGNORE);
232:   PetscFree(send_reqs);
233:   PetscFree(recv_reqs);
234:   PetscFree(send_buffer);
235:   PetscFree(recv_buffer);
236:   PetscFree(ptrs_buffer);

238:   /* Compute B and B_delta (local actions) */
239:   PetscMalloc1(pcis->n_neigh,&aux_sums);
240:   PetscMalloc1(n_local_lambda,&l2g_indices);
241:   PetscMalloc1(n_local_lambda,&vals_B_delta);
242:   PetscMalloc1(n_local_lambda,&cols_B_delta);
243:   PetscMalloc1(n_local_lambda,&scaling_factors);
244:   ISGetIndices(subset_n,&aux_global_numbering);
245:   partial_sum=0;
246:   cum = 0;
247:   for (i=0;i<dual_size;i++) {
248:     n_global_lambda = aux_global_numbering[cum];
249:     j = mat_graph->count[aux_local_numbering_1[i]];
250:     aux_sums[0]=0;
251:     for (s=1;s<j;s++) {
252:       aux_sums[s]=aux_sums[s-1]+j-s+1;
253:     }
254:     array = all_factors[aux_local_numbering_1[i]];
255:     n_neg_values = 0;
256:     while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
257:     n_pos_values = j - n_neg_values;
258:     if (fully_redundant) {
259:       for (s=0;s<n_neg_values;s++) {
260:         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
261:         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
262:         vals_B_delta   [partial_sum+s]=-1.0;
263:         scaling_factors[partial_sum+s]=array[s];
264:       }
265:       for (s=0;s<n_pos_values;s++) {
266:         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
267:         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
268:         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
269:         scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
270:       }
271:       partial_sum += j;
272:     } else {
273:       /* l2g_indices and default cols and vals of B_delta */
274:       for (s=0;s<j;s++) {
275:         l2g_indices    [partial_sum+s]=n_global_lambda+s;
276:         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
277:         vals_B_delta   [partial_sum+s]=0.0;
278:       }
279:       /* B_delta */
280:       if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
281:         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
282:       }
283:       if ( n_neg_values < j ) { /* there's a rank next to me to the right */
284:         vals_B_delta   [partial_sum+n_neg_values]=1.0;
285:       }
286:       /* scaling as in Klawonn-Widlund 1999*/
287:       for (s=0;s<n_neg_values;s++) {
288:         scalar_value = 0.0;
289:         for (k=0;k<s+1;k++) {
290:           scalar_value += array[k];
291:         }
292:         scaling_factors[partial_sum+s] = -scalar_value;
293:       }
294:       for (s=0;s<n_pos_values;s++) {
295:         scalar_value = 0.0;
296:         for (k=s+n_neg_values;k<j;k++) {
297:           scalar_value += array[k];
298:         }
299:         scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
300:       }
301:       partial_sum += j;
302:     }
303:     cum += aux_local_numbering_2[i];
304:   }
305:   ISRestoreIndices(subset_n,&aux_global_numbering);
306:   ISDestroy(&subset_mult);
307:   ISDestroy(&subset_n);
308:   PetscFree(aux_sums);
309:   PetscFree(aux_local_numbering_1);
310:   PetscFree(dual_dofs_boundary_indices);
311:   PetscFree(all_factors[0]);
312:   PetscFree(all_factors);

314:   /* Local to global mapping of fetidpmat */
315:   VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);
316:   VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);
317:   VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);
318:   VecCreate(comm,&lambda_global);
319:   VecSetSizes(lambda_global,PETSC_DECIDE,fetidpmat_ctx->n_lambda);
320:   VecSetType(lambda_global,VECMPI);
321:   ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);
322:   VecScatterCreate(fetidpmat_ctx->lambda_local,(IS)0,lambda_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);
323:   ISDestroy(&IS_l2g_lambda);

325:   /* Create local part of B_delta */
326:   MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);
327:   MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);
328:   MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);
329:   MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);
330:   MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
331:   for (i=0;i<n_local_lambda;i++) {
332:     MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);
333:   }
334:   PetscFree(vals_B_delta);
335:   MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);
336:   MatAssemblyEnd  (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);

338:   if (fully_redundant) {
339:     MatCreate(PETSC_COMM_SELF,&ScalingMat);
340:     MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);
341:     MatSetType(ScalingMat,MATSEQAIJ);
342:     MatSeqAIJSetPreallocation(ScalingMat,1,NULL);
343:     for (i=0;i<n_local_lambda;i++) {
344:       MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);
345:     }
346:     MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);
347:     MatAssemblyEnd  (ScalingMat,MAT_FINAL_ASSEMBLY);
348:     MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);
349:     MatDestroy(&ScalingMat);
350:   } else {
351:     MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);
352:     MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);
353:     MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);
354:     MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);
355:     for (i=0;i<n_local_lambda;i++) {
356:       MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);
357:     }
358:     MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);
359:     MatAssemblyEnd  (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);
360:   }
361:   PetscFree(scaling_factors);
362:   PetscFree(cols_B_delta);

364:   /* Create some vectors needed by fetidp */
365:   VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);
366:   VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);

368:   test_fetidp = PETSC_FALSE;
369:   PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);

371:   if (test_fetidp && !pcbddc->use_deluxe_scaling) {

373:     PetscReal real_value;

375:     PetscViewerASCIIGetStdout(comm,&viewer);
376:     PetscViewerASCIIPushSynchronized(viewer);
377:     PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");
378:     PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");
379:     PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");
380:     if (fully_redundant) {
381:       PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");
382:     } else {
383:       PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");
384:     }
385:     PetscViewerFlush(viewer);

387:     /******************************************************************/
388:     /* TEST A/B: Test numbering of global lambda dofs             */
389:     /******************************************************************/

391:     VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);
392:     VecSet(lambda_global,1.0);
393:     VecSet(test_vec,1.0);
394:     VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
395:     VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
396:     scalar_value = -1.0;
397:     VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);
398:     VecNorm(test_vec,NORM_INFINITY,&real_value);
399:     VecDestroy(&test_vec);
400:     PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);
401:     PetscViewerFlush(viewer);
402:     if (fully_redundant) {
403:       VecSet(lambda_global,0.0);
404:       VecSet(fetidpmat_ctx->lambda_local,0.5);
405:       VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
406:       VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
407:       VecSum(lambda_global,&scalar_value);
408:       PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);
409:       PetscViewerFlush(viewer);
410:     }

412:     /******************************************************************/
413:     /* TEST C: It should holds B_delta*w=0, w\in\widehat{W}           */
414:     /* This is the meaning of the B matrix                            */
415:     /******************************************************************/

417:     VecSetRandom(pcis->vec1_N,NULL);
418:     VecSet(pcis->vec1_global,0.0);
419:     VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
420:     VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
421:     VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
422:     VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
423:     VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
424:     VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
425:     /* Action of B_delta */
426:     MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
427:     VecSet(lambda_global,0.0);
428:     VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
429:     VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
430:     VecNorm(lambda_global,NORM_INFINITY,&real_value);
431:     PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);
432:     PetscViewerFlush(viewer);

434:     /******************************************************************/
435:     /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W}     */
436:     /* E_D = R_D^TR                                                   */
437:     /* P_D = B_{D,delta}^T B_{delta}                                  */
438:     /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
439:     /******************************************************************/

441:     /* compute a random vector in \widetilde{W} */
442:     VecSetRandom(pcis->vec1_N,NULL);
443:     scalar_value = 0.0;  /* set zero at vertices */
444:     VecGetArray(pcis->vec1_N,&array);
445:     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
446:     VecRestoreArray(pcis->vec1_N,&array);
447:     /* store w for final comparison */
448:     VecDuplicate(pcis->vec1_B,&test_vec);
449:     VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
450:     VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);

452:     /* Jump operator P_D : results stored in pcis->vec1_B */

454:     VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
455:     VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
456:     /* Action of B_delta */
457:     MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
458:     VecSet(lambda_global,0.0);
459:     VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
460:     VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
461:     /* Action of B_Ddelta^T */
462:     VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
463:     VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
464:     MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);

466:     /* Average operator E_D : results stored in pcis->vec2_B */
467:     VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
468:     VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
469:     PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);
470:     VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
471:     VecScatterEnd  (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);

473:     /* test E_D=I-P_D */
474:     scalar_value = 1.0;
475:     VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);
476:     scalar_value = -1.0;
477:     VecAXPY(pcis->vec1_B,scalar_value,test_vec);
478:     VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);
479:     VecDestroy(&test_vec);
480:     PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);
481:     PetscViewerFlush(viewer);

483:     /******************************************************************/
484:     /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W}          */
485:     /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
486:     /******************************************************************/

488:     VecSetRandom(pcis->vec1_N,NULL);
489:     VecGetArray(pcis->vec1_N,&array);
490:     scalar_value = 0.0;  /* set zero at vertices */
491:     for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
492:     VecRestoreArray(pcis->vec1_N,&array);

494:     /* Jump operator P_D : results stored in pcis->vec1_B */

496:     VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
497:     VecScatterEnd  (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
498:     /* Action of B_delta */
499:     MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
500:     VecSet(lambda_global,0.0);
501:     VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
502:     VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
503:     /* Action of B_Ddelta^T */
504:     VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
505:     VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
506:     MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
507:     /* scaling */
508:     PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);
509:     VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);
510:     PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);
511:     PetscViewerFlush(viewer);

513:     if (!fully_redundant) {
514:       /******************************************************************/
515:       /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
516:       /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
517:       /******************************************************************/
518:       VecDuplicate(lambda_global,&test_vec);
519:       VecSetRandom(lambda_global,NULL);
520:       /* Action of B_Ddelta^T */
521:       VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
522:       VecScatterEnd  (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
523:       MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
524:       /* Action of B_delta */
525:       MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
526:       VecSet(test_vec,0.0);
527:       VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
528:       VecScatterEnd  (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
529:       scalar_value = -1.0;
530:       VecAXPY(lambda_global,scalar_value,test_vec);
531:       VecNorm(lambda_global,NORM_INFINITY,&real_value);
532:       PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);
533:       PetscViewerFlush(viewer);
534:       PetscViewerFlush(viewer);
535:       VecDestroy(&test_vec);
536:     }
537:   }
538:   /* final cleanup */
539:   VecDestroy(&lambda_global);

541:   return(0);
542: }

546: PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
547: {
548:   FETIDPMat_ctx  mat_ctx;
549:   PC_IS          *pcis;

553:   MatShellGetContext(fetimat,(void**)&mat_ctx);
554:   /* get references from objects created when setting up feti mat context */
555:   PetscObjectReference((PetscObject)mat_ctx->lambda_local);
556:   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
557:   PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);
558:   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
559:   PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);
560:   fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
561:   /* create local Schur complement matrix */
562:   pcis = (PC_IS*)fetidppc_ctx->pc->data;
563:   MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);
564:   MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);
565:   return(0);
566: }

570: PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
571: {
572:   FETIDPMat_ctx  mat_ctx;
573:   PC_IS          *pcis;

577:   MatShellGetContext(fetimat,(void**)&mat_ctx);
578:   pcis = (PC_IS*)mat_ctx->pc->data;
579:   /* Application of B_delta^T */
580:   VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
581:   VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
582:   MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
583:   /* Application of \widetilde{S}^-1 */
584:   VecSet(pcis->vec1_D,0.0);
585:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
586:   /* Application of B_delta */
587:   MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
588:   VecSet(y,0.0);
589:   VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
590:   VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
591:   return(0);
592: }

596: PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
597: {
598:   FETIDPMat_ctx  mat_ctx;
599:   PC_IS          *pcis;

603:   MatShellGetContext(fetimat,(void**)&mat_ctx);
604:   pcis = (PC_IS*)mat_ctx->pc->data;
605:   /* Application of B_delta^T */
606:   VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
607:   VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
608:   MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
609:   /* Application of \widetilde{S}^-1 */
610:   VecSet(pcis->vec1_D,0.0);
611:   PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);
612:   /* Application of B_delta */
613:   MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
614:   VecSet(y,0.0);
615:   VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
616:   VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
617:   return(0);
618: }

622: PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
623: {
624:   FETIDPPC_ctx   pc_ctx;
625:   PC_IS          *pcis;

629:   PCShellGetContext(fetipc,(void**)&pc_ctx);
630:   pcis = (PC_IS*)pc_ctx->pc->data;
631:   /* Application of B_Ddelta^T */
632:   VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
633:   VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
634:   VecSet(pcis->vec2_B,0.0);
635:   MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);
636:   /* Application of local Schur complement */
637:   MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);
638:   /* Application of B_Ddelta */
639:   MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);
640:   VecSet(y,0.0);
641:   VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
642:   VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
643:   return(0);
644: }

648: PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
649: {
650:   FETIDPPC_ctx   pc_ctx;
651:   PC_IS          *pcis;

655:   PCShellGetContext(fetipc,(void**)&pc_ctx);
656:   pcis = (PC_IS*)pc_ctx->pc->data;
657:   /* Application of B_Ddelta^T */
658:   VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
659:   VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
660:   VecSet(pcis->vec2_B,0.0);
661:   MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);
662:   /* Application of local Schur complement */
663:   MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);
664:   /* Application of B_Ddelta */
665:   MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);
666:   VecSet(y,0.0);
667:   VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
668:   VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
669:   return(0);
670: }