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