Actual source code: bddcfetidp.c
petsc-3.6.1 2015-08-06
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,"-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->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
175: VecScatterEnd(matis->ctx,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(matis->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) {
187: SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in %s: global number of multipliers mismatch! (%d!=%d)\n",__FUNCT__,fetidpmat_ctx->n_lambda,i);
188: }
190: /* init data for scaling factors exchange */
191: partial_sum = 0;
192: j = 0;
193: PetscMalloc1(pcis->n_neigh,&ptrs_buffer);
194: PetscMalloc1(pcis->n_neigh-1,&send_reqs);
195: PetscMalloc1(pcis->n_neigh-1,&recv_reqs);
196: PetscMalloc1(pcis->n,&all_factors);
197: ptrs_buffer[0]=0;
198: for (i=1;i<pcis->n_neigh;i++) {
199: partial_sum += pcis->n_shared[i];
200: ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
201: }
202: PetscMalloc1(partial_sum,&send_buffer);
203: PetscMalloc1(partial_sum,&recv_buffer);
204: PetscMalloc1(partial_sum,&all_factors[0]);
205: for (i=0;i<pcis->n-1;i++) {
206: j = mat_graph->count[i];
207: all_factors[i+1]=all_factors[i]+j;
208: }
209: /* scatter B scaling to N vec */
210: VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
211: VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
212: /* communications */
213: VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);
214: for (i=1;i<pcis->n_neigh;i++) {
215: for (j=0;j<pcis->n_shared[i];j++) {
216: send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
217: }
218: PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);
219: PetscMPIIntCast(pcis->neigh[i],&neigh);
220: MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);
221: MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);
222: }
223: VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);
224: MPI_Waitall((pcis->n_neigh-1),recv_reqs,MPI_STATUSES_IGNORE);
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: MPI_Waitall((pcis->n_neigh-1),send_reqs,MPI_STATUSES_IGNORE);
235: PetscFree(send_reqs);
236: PetscFree(recv_reqs);
237: PetscFree(send_buffer);
238: PetscFree(recv_buffer);
239: PetscFree(ptrs_buffer);
241: /* Compute B and B_delta (local actions) */
242: PetscMalloc1(pcis->n_neigh,&aux_sums);
243: PetscMalloc1(n_local_lambda,&l2g_indices);
244: PetscMalloc1(n_local_lambda,&vals_B_delta);
245: PetscMalloc1(n_local_lambda,&cols_B_delta);
246: PetscMalloc1(n_local_lambda,&scaling_factors);
247: ISGetIndices(subset_n,&aux_global_numbering);
248: n_global_lambda=0;
249: partial_sum=0;
250: cum = 0;
251: for (i=0;i<dual_size;i++) {
252: n_global_lambda = aux_global_numbering[cum];
253: j = mat_graph->count[aux_local_numbering_1[i]];
254: aux_sums[0]=0;
255: for (s=1;s<j;s++) {
256: aux_sums[s]=aux_sums[s-1]+j-s+1;
257: }
258: array = all_factors[aux_local_numbering_1[i]];
259: n_neg_values = 0;
260: while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
261: n_pos_values = j - n_neg_values;
262: if (fully_redundant) {
263: for (s=0;s<n_neg_values;s++) {
264: l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
265: cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i];
266: vals_B_delta [partial_sum+s]=-1.0;
267: scaling_factors[partial_sum+s]=array[s];
268: }
269: for (s=0;s<n_pos_values;s++) {
270: l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
271: cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
272: vals_B_delta [partial_sum+s+n_neg_values]=1.0;
273: scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
274: }
275: partial_sum += j;
276: } else {
277: /* l2g_indices and default cols and vals of B_delta */
278: for (s=0;s<j;s++) {
279: l2g_indices [partial_sum+s]=n_global_lambda+s;
280: cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i];
281: vals_B_delta [partial_sum+s]=0.0;
282: }
283: /* B_delta */
284: if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
285: vals_B_delta [partial_sum+n_neg_values-1]=-1.0;
286: }
287: if ( n_neg_values < j ) { /* there's a rank next to me to the right */
288: vals_B_delta [partial_sum+n_neg_values]=1.0;
289: }
290: /* scaling as in Klawonn-Widlund 1999*/
291: for (s=0;s<n_neg_values;s++) {
292: scalar_value = 0.0;
293: for (k=0;k<s+1;k++) {
294: scalar_value += array[k];
295: }
296: scaling_factors[partial_sum+s] = -scalar_value;
297: }
298: for (s=0;s<n_pos_values;s++) {
299: scalar_value = 0.0;
300: for (k=s+n_neg_values;k<j;k++) {
301: scalar_value += array[k];
302: }
303: scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
304: }
305: partial_sum += j;
306: }
307: cum += aux_local_numbering_2[i];
308: }
309: ISRestoreIndices(subset_n,&aux_global_numbering);
310: ISDestroy(&subset_mult);
311: ISDestroy(&subset_n);
312: PetscFree(aux_sums);
313: PetscFree(aux_local_numbering_1);
314: PetscFree(dual_dofs_boundary_indices);
315: PetscFree(all_factors[0]);
316: PetscFree(all_factors);
318: /* Local to global mapping of fetidpmat */
319: VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);
320: VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);
321: VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);
322: VecCreate(comm,&lambda_global);
323: VecSetSizes(lambda_global,PETSC_DECIDE,fetidpmat_ctx->n_lambda);
324: VecSetType(lambda_global,VECMPI);
325: ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);
326: VecScatterCreate(fetidpmat_ctx->lambda_local,(IS)0,lambda_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);
327: ISDestroy(&IS_l2g_lambda);
329: /* Create local part of B_delta */
330: MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);
331: MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);
332: MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);
333: MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);
334: MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
335: for (i=0;i<n_local_lambda;i++) {
336: MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);
337: }
338: PetscFree(vals_B_delta);
339: MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);
340: MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);
342: if (fully_redundant) {
343: MatCreate(PETSC_COMM_SELF,&ScalingMat);
344: MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);
345: MatSetType(ScalingMat,MATSEQAIJ);
346: MatSeqAIJSetPreallocation(ScalingMat,1,NULL);
347: for (i=0;i<n_local_lambda;i++) {
348: MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);
349: }
350: MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);
351: MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);
352: MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);
353: MatDestroy(&ScalingMat);
354: } else {
355: MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);
356: MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);
357: MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);
358: MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);
359: for (i=0;i<n_local_lambda;i++) {
360: MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);
361: }
362: MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);
363: MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);
364: }
365: PetscFree(scaling_factors);
366: PetscFree(cols_B_delta);
368: /* Create some vectors needed by fetidp */
369: VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);
370: VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);
372: test_fetidp = PETSC_FALSE;
373: PetscOptionsGetBool(NULL,"-fetidp_check",&test_fetidp,NULL);
375: if (test_fetidp && !pcbddc->use_deluxe_scaling) {
377: PetscReal real_value;
379: PetscViewerASCIIGetStdout(comm,&viewer);
380: PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
381: PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");
382: PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");
383: PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");
384: if (fully_redundant) {
385: PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");
386: } else {
387: PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");
388: }
389: PetscViewerFlush(viewer);
391: /******************************************************************/
392: /* TEST A/B: Test numbering of global lambda dofs */
393: /******************************************************************/
395: VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);
396: VecSet(lambda_global,1.0);
397: VecSet(test_vec,1.0);
398: VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
399: VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
400: scalar_value = -1.0;
401: VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);
402: VecNorm(test_vec,NORM_INFINITY,&real_value);
403: VecDestroy(&test_vec);
404: PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);
405: PetscViewerFlush(viewer);
406: if (fully_redundant) {
407: VecSet(lambda_global,0.0);
408: VecSet(fetidpmat_ctx->lambda_local,0.5);
409: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
410: VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
411: VecSum(lambda_global,&scalar_value);
412: PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);
413: PetscViewerFlush(viewer);
414: }
416: /******************************************************************/
417: /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */
418: /* This is the meaning of the B matrix */
419: /******************************************************************/
421: VecSetRandom(pcis->vec1_N,NULL);
422: VecSet(pcis->vec1_global,0.0);
423: VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
424: VecScatterEnd (matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
425: VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
426: VecScatterEnd (matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
427: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
428: VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
429: /* Action of B_delta */
430: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
431: VecSet(lambda_global,0.0);
432: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
433: VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
434: VecNorm(lambda_global,NORM_INFINITY,&real_value);
435: PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);
436: PetscViewerFlush(viewer);
438: /******************************************************************/
439: /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */
440: /* E_D = R_D^TR */
441: /* P_D = B_{D,delta}^T B_{delta} */
442: /* eq.44 Mandel Tezaur and Dohrmann 2005 */
443: /******************************************************************/
445: /* compute a random vector in \widetilde{W} */
446: VecSetRandom(pcis->vec1_N,NULL);
447: scalar_value = 0.0; /* set zero at vertices */
448: VecGetArray(pcis->vec1_N,&array);
449: for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
450: VecRestoreArray(pcis->vec1_N,&array);
451: /* store w for final comparison */
452: VecDuplicate(pcis->vec1_B,&test_vec);
453: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
454: VecScatterEnd (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
456: /* Jump operator P_D : results stored in pcis->vec1_B */
458: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
459: VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
460: /* Action of B_delta */
461: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
462: VecSet(lambda_global,0.0);
463: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
464: VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
465: /* Action of B_Ddelta^T */
466: VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
467: VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
468: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
470: /* Average operator E_D : results stored in pcis->vec2_B */
471: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
472: VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
473: PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);
474: VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
475: VecScatterEnd (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
477: /* test E_D=I-P_D */
478: scalar_value = 1.0;
479: VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);
480: scalar_value = -1.0;
481: VecAXPY(pcis->vec1_B,scalar_value,test_vec);
482: VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);
483: VecDestroy(&test_vec);
484: PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);
485: PetscViewerFlush(viewer);
487: /******************************************************************/
488: /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W} */
489: /* eq.48 Mandel Tezaur and Dohrmann 2005 */
490: /******************************************************************/
492: VecSetRandom(pcis->vec1_N,NULL);
493: VecGetArray(pcis->vec1_N,&array);
494: scalar_value = 0.0; /* set zero at vertices */
495: for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; }
496: VecRestoreArray(pcis->vec1_N,&array);
498: /* Jump operator P_D : results stored in pcis->vec1_B */
500: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
501: VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
502: /* Action of B_delta */
503: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
504: VecSet(lambda_global,0.0);
505: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
506: VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);
507: /* Action of B_Ddelta^T */
508: VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
509: VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
510: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
511: /* scaling */
512: PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);
513: VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);
514: PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);
515: PetscViewerFlush(viewer);
517: if (!fully_redundant) {
518: /******************************************************************/
519: /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */
520: /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */
521: /******************************************************************/
522: VecDuplicate(lambda_global,&test_vec);
523: VecSetRandom(lambda_global,NULL);
524: /* Action of B_Ddelta^T */
525: VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
526: VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
527: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
528: /* Action of B_delta */
529: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
530: VecSet(test_vec,0.0);
531: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
532: VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
533: scalar_value = -1.0;
534: VecAXPY(lambda_global,scalar_value,test_vec);
535: VecNorm(lambda_global,NORM_INFINITY,&real_value);
536: PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);
537: PetscViewerFlush(viewer);
538: PetscViewerFlush(viewer);
539: VecDestroy(&test_vec);
540: }
541: }
542: /* final cleanup */
543: VecDestroy(&lambda_global);
545: return(0);
546: }
550: PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
551: {
552: FETIDPMat_ctx mat_ctx;
553: PC_IS *pcis;
557: MatShellGetContext(fetimat,(void**)&mat_ctx);
558: /* get references from objects created when setting up feti mat context */
559: PetscObjectReference((PetscObject)mat_ctx->lambda_local);
560: fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
561: PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);
562: fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
563: PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);
564: fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
565: /* create local Schur complement matrix */
566: pcis = (PC_IS*)fetidppc_ctx->pc->data;
567: MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);
568: MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);
569: return(0);
570: }
574: PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
575: {
576: FETIDPMat_ctx mat_ctx;
577: PC_IS *pcis;
581: MatShellGetContext(fetimat,(void**)&mat_ctx);
582: pcis = (PC_IS*)mat_ctx->pc->data;
583: /* Application of B_delta^T */
584: VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
585: VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
586: MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
587: /* Application of \widetilde{S}^-1 */
588: VecSet(pcis->vec1_D,0.0);
589: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);
590: /* Application of B_delta */
591: MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
592: VecSet(y,0.0);
593: VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
594: VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
595: return(0);
596: }
600: PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
601: {
602: FETIDPMat_ctx mat_ctx;
603: PC_IS *pcis;
607: MatShellGetContext(fetimat,(void**)&mat_ctx);
608: pcis = (PC_IS*)mat_ctx->pc->data;
609: /* Application of B_delta^T */
610: VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
611: VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
612: MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
613: /* Application of \widetilde{S}^-1 */
614: VecSet(pcis->vec1_D,0.0);
615: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);
616: /* Application of B_delta */
617: MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
618: VecSet(y,0.0);
619: VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
620: VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
621: return(0);
622: }
626: PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y)
627: {
628: FETIDPPC_ctx pc_ctx;
629: PC_IS *pcis;
633: PCShellGetContext(fetipc,(void**)&pc_ctx);
634: pcis = (PC_IS*)pc_ctx->pc->data;
635: /* Application of B_Ddelta^T */
636: VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
637: VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
638: VecSet(pcis->vec2_B,0.0);
639: MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);
640: /* Application of local Schur complement */
641: MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);
642: /* Application of B_Ddelta */
643: MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);
644: VecSet(y,0.0);
645: VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
646: VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
647: return(0);
648: }
652: PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y)
653: {
654: FETIDPPC_ctx pc_ctx;
655: PC_IS *pcis;
659: PCShellGetContext(fetipc,(void**)&pc_ctx);
660: pcis = (PC_IS*)pc_ctx->pc->data;
661: /* Application of B_Ddelta^T */
662: VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
663: VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
664: VecSet(pcis->vec2_B,0.0);
665: MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);
666: /* Application of local Schur complement */
667: MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);
668: /* Application of B_Ddelta */
669: MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);
670: VecSet(y,0.0);
671: VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
672: VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
673: return(0);
674: }