Actual source code: bddcfetidp.c
petsc-3.10.5 2019-03-28
1: #include <../src/ksp/pc/impls/bddc/bddc.h>
2: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3: #include <petscblaslapack.h>
5: static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
6: {
7: BDdelta_DN ctx;
11: MatShellGetContext(A,(void**)&ctx);
12: MatMultTranspose(ctx->BD,x,ctx->work);
13: KSPSolveTranspose(ctx->kBD,ctx->work,y);
14: return(0);
15: }
17: static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
18: {
19: BDdelta_DN ctx;
23: MatShellGetContext(A,(void**)&ctx);
24: KSPSolve(ctx->kBD,x,ctx->work);
25: MatMult(ctx->BD,ctx->work,y);
26: return(0);
27: }
29: static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
30: {
31: BDdelta_DN ctx;
35: MatShellGetContext(A,(void**)&ctx);
36: MatDestroy(&ctx->BD);
37: KSPDestroy(&ctx->kBD);
38: VecDestroy(&ctx->work);
39: PetscFree(ctx);
40: return(0);
41: }
44: PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
45: {
46: FETIDPMat_ctx newctx;
50: PetscNew(&newctx);
51: /* increase the reference count for BDDC preconditioner */
52: PetscObjectReference((PetscObject)pc);
53: newctx->pc = pc;
54: *fetidpmat_ctx = newctx;
55: return(0);
56: }
58: PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
59: {
60: FETIDPPC_ctx newctx;
64: PetscNew(&newctx);
65: /* increase the reference count for BDDC preconditioner */
66: PetscObjectReference((PetscObject)pc);
67: newctx->pc = pc;
68: *fetidppc_ctx = newctx;
69: return(0);
70: }
72: PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
73: {
74: FETIDPMat_ctx mat_ctx;
78: MatShellGetContext(A,(void**)&mat_ctx);
79: VecDestroy(&mat_ctx->lambda_local);
80: VecDestroy(&mat_ctx->temp_solution_D);
81: VecDestroy(&mat_ctx->temp_solution_B);
82: MatDestroy(&mat_ctx->B_delta);
83: MatDestroy(&mat_ctx->B_Ddelta);
84: MatDestroy(&mat_ctx->B_BB);
85: MatDestroy(&mat_ctx->B_BI);
86: MatDestroy(&mat_ctx->Bt_BB);
87: MatDestroy(&mat_ctx->Bt_BI);
88: MatDestroy(&mat_ctx->C);
89: VecDestroy(&mat_ctx->rhs_flip);
90: VecDestroy(&mat_ctx->vP);
91: VecDestroy(&mat_ctx->xPg);
92: VecDestroy(&mat_ctx->yPg);
93: VecScatterDestroy(&mat_ctx->l2g_lambda);
94: VecScatterDestroy(&mat_ctx->l2g_lambda_only);
95: VecScatterDestroy(&mat_ctx->l2g_p);
96: VecScatterDestroy(&mat_ctx->g2g_p);
97: PCDestroy(&mat_ctx->pc); /* decrease PCBDDC reference count */
98: ISDestroy(&mat_ctx->pressure);
99: ISDestroy(&mat_ctx->lagrange);
100: PetscFree(mat_ctx);
101: return(0);
102: }
104: PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
105: {
106: FETIDPPC_ctx pc_ctx;
110: PCShellGetContext(pc,(void**)&pc_ctx);
111: VecDestroy(&pc_ctx->lambda_local);
112: MatDestroy(&pc_ctx->B_Ddelta);
113: VecScatterDestroy(&pc_ctx->l2g_lambda);
114: MatDestroy(&pc_ctx->S_j);
115: PCDestroy(&pc_ctx->pc); /* decrease PCBDDC reference count */
116: VecDestroy(&pc_ctx->xPg);
117: VecDestroy(&pc_ctx->yPg);
118: PetscFree(pc_ctx);
119: return(0);
120: }
122: PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx )
123: {
125: PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data;
126: PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data;
127: PCBDDCGraph mat_graph=pcbddc->mat_graph;
128: #if defined(PETSC_USE_DEBUG)
129: Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data;
130: #endif
131: MPI_Comm comm;
132: Mat ScalingMat,BD1,BD2;
133: Vec fetidp_global;
134: IS IS_l2g_lambda;
135: IS subset,subset_mult,subset_n,isvert;
136: PetscBool skip_node,fully_redundant;
137: PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum;
138: PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values;
139: PetscMPIInt rank,size,buf_size,neigh;
140: PetscScalar scalar_value;
141: const PetscInt *vertex_indices;
142: PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1;
143: const PetscInt *aux_global_numbering;
144: PetscInt *aux_sums,*cols_B_delta,*l2g_indices;
145: PetscScalar *array,*scaling_factors,*vals_B_delta;
146: PetscScalar **all_factors;
147: PetscInt *aux_local_numbering_2;
148: PetscLayout llay;
150: /* saddlepoint */
151: ISLocalToGlobalMapping l2gmap_p;
152: PetscLayout play;
153: IS gP,pP;
154: PetscInt nPl,nPg,nPgl;
157: PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);
158: MPI_Comm_rank(comm,&rank);
159: MPI_Comm_size(comm,&size);
161: /* saddlepoint */
162: nPl = 0;
163: nPg = 0;
164: nPgl = 0;
165: gP = NULL;
166: pP = NULL;
167: l2gmap_p = NULL;
168: play = NULL;
169: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);
170: if (pP) { /* saddle point */
171: /* subdomain pressures in global numbering */
172: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);
173: if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
174: ISGetLocalSize(gP,&nPl);
175: VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);
176: VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);
177: VecSetType(fetidpmat_ctx->vP,VECSTANDARD);
178: VecSetUp(fetidpmat_ctx->vP);
180: /* pressure matrix */
181: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);
182: if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
183: IS Pg;
185: ISRenumber(gP,NULL,&nPg,&Pg);
186: ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);
187: ISDestroy(&Pg);
188: PetscLayoutCreate(comm,&play);
189: PetscLayoutSetBlockSize(play,1);
190: PetscLayoutSetSize(play,nPg);
191: ISGetLocalSize(pP,&nPgl);
192: PetscLayoutSetLocalSize(play,nPgl);
193: PetscLayoutSetUp(play);
194: } else {
195: PetscObjectReference((PetscObject)fetidpmat_ctx->C);
196: MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);
197: PetscObjectReference((PetscObject)l2gmap_p);
198: MatGetSize(fetidpmat_ctx->C,&nPg,NULL);
199: MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);
200: MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);
201: PetscLayoutReference(llay,&play);
202: }
203: VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);
204: VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);
206: /* import matrices for pressures coupling */
207: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);
208: if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
209: PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);
211: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);
212: if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
213: PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);
215: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);
216: if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
217: PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);
219: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);
220: if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
221: PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);
223: PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);
224: if (fetidpmat_ctx->rhs_flip) {
225: PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);
226: }
227: }
229: /* Default type of lagrange multipliers is non-redundant */
230: fully_redundant = fetidpmat_ctx->fully_redundant;
232: /* Evaluate local and global number of lagrange multipliers */
233: VecSet(pcis->vec1_N,0.0);
234: n_local_lambda = 0;
235: partial_sum = 0;
236: n_boundary_dofs = 0;
237: s = 0;
239: /* Get Vertices used to define the BDDC */
240: PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);
241: ISGetLocalSize(isvert,&n_vertices);
242: ISGetIndices(isvert,&vertex_indices);
244: dual_size = pcis->n_B-n_vertices;
245: PetscMalloc1(dual_size,&dual_dofs_boundary_indices);
246: PetscMalloc1(dual_size,&aux_local_numbering_1);
247: PetscMalloc1(dual_size,&aux_local_numbering_2);
249: VecGetArray(pcis->vec1_N,&array);
250: for (i=0;i<pcis->n;i++){
251: j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
252: if (j > 0) n_boundary_dofs++;
253: skip_node = PETSC_FALSE;
254: if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
255: skip_node = PETSC_TRUE;
256: s++;
257: }
258: if (j < 1) skip_node = PETSC_TRUE;
259: if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
260: if (!skip_node) {
261: if (fully_redundant) {
262: /* fully redundant set of lagrange multipliers */
263: n_lambda_for_dof = (j*(j+1))/2;
264: } else {
265: n_lambda_for_dof = j;
266: }
267: n_local_lambda += j;
268: /* needed to evaluate global number of lagrange multipliers */
269: array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
270: /* store some data needed */
271: dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
272: aux_local_numbering_1[partial_sum] = i;
273: aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
274: partial_sum++;
275: }
276: }
277: VecRestoreArray(pcis->vec1_N,&array);
278: ISRestoreIndices(isvert,&vertex_indices);
279: PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);
280: dual_size = partial_sum;
282: /* compute global ordering of lagrange multipliers and associate l2g map */
283: ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);
284: ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);
285: ISDestroy(&subset_n);
286: ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);
287: ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);
288: ISDestroy(&subset);
290: #if defined(PETSC_USE_DEBUG)
291: VecSet(pcis->vec1_global,0.0);
292: VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
293: VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
294: VecSum(pcis->vec1_global,&scalar_value);
295: i = (PetscInt)PetscRealPart(scalar_value);
296: if (i != fetidpmat_ctx->n_lambda) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global number of multipliers mismatch! (%D != %D)",fetidpmat_ctx->n_lambda,i);
297: #endif
299: /* init data for scaling factors exchange */
300: if (!pcbddc->use_deluxe_scaling) {
301: PetscInt *ptrs_buffer,neigh_position;
302: PetscScalar *send_buffer,*recv_buffer;
303: MPI_Request *send_reqs,*recv_reqs;
305: partial_sum = 0;
306: PetscMalloc1(pcis->n_neigh,&ptrs_buffer);
307: PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);
308: PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);
309: PetscMalloc1(pcis->n+1,&all_factors);
310: if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
311: for (i=1;i<pcis->n_neigh;i++) {
312: partial_sum += pcis->n_shared[i];
313: ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
314: }
315: PetscMalloc1(partial_sum,&send_buffer);
316: PetscMalloc1(partial_sum,&recv_buffer);
317: PetscMalloc1(partial_sum,&all_factors[0]);
318: for (i=0;i<pcis->n-1;i++) {
319: j = mat_graph->count[i];
320: all_factors[i+1]=all_factors[i]+j;
321: }
323: /* scatter B scaling to N vec */
324: VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
325: VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);
326: /* communications */
327: VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);
328: for (i=1;i<pcis->n_neigh;i++) {
329: for (j=0;j<pcis->n_shared[i];j++) {
330: send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
331: }
332: PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);
333: PetscMPIIntCast(pcis->neigh[i],&neigh);
334: MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);
335: MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);
336: }
337: VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);
338: if (pcis->n_neigh > 0) {
339: MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);
340: }
341: /* put values in correct places */
342: for (i=1;i<pcis->n_neigh;i++) {
343: for (j=0;j<pcis->n_shared[i];j++) {
344: k = pcis->shared[i][j];
345: neigh_position = 0;
346: while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
347: all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
348: }
349: }
350: if (pcis->n_neigh > 0) {
351: MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);
352: }
353: PetscFree(send_reqs);
354: PetscFree(recv_reqs);
355: PetscFree(send_buffer);
356: PetscFree(recv_buffer);
357: PetscFree(ptrs_buffer);
358: }
360: /* Compute B and B_delta (local actions) */
361: PetscMalloc1(pcis->n_neigh,&aux_sums);
362: PetscMalloc1(n_local_lambda,&l2g_indices);
363: PetscMalloc1(n_local_lambda,&vals_B_delta);
364: PetscMalloc1(n_local_lambda,&cols_B_delta);
365: if (!pcbddc->use_deluxe_scaling) {
366: PetscMalloc1(n_local_lambda,&scaling_factors);
367: } else {
368: scaling_factors = NULL;
369: all_factors = NULL;
370: }
371: ISGetIndices(subset_n,&aux_global_numbering);
372: partial_sum=0;
373: cum = 0;
374: for (i=0;i<dual_size;i++) {
375: n_global_lambda = aux_global_numbering[cum];
376: j = mat_graph->count[aux_local_numbering_1[i]];
377: aux_sums[0]=0;
378: for (s=1;s<j;s++) {
379: aux_sums[s]=aux_sums[s-1]+j-s+1;
380: }
381: if (all_factors) array = all_factors[aux_local_numbering_1[i]];
382: n_neg_values = 0;
383: while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
384: n_pos_values = j - n_neg_values;
385: if (fully_redundant) {
386: for (s=0;s<n_neg_values;s++) {
387: l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
388: cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i];
389: vals_B_delta [partial_sum+s]=-1.0;
390: if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s];
391: }
392: for (s=0;s<n_pos_values;s++) {
393: l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
394: cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
395: vals_B_delta [partial_sum+s+n_neg_values]=1.0;
396: if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
397: }
398: partial_sum += j;
399: } else {
400: /* l2g_indices and default cols and vals of B_delta */
401: for (s=0;s<j;s++) {
402: l2g_indices [partial_sum+s]=n_global_lambda+s;
403: cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i];
404: vals_B_delta [partial_sum+s]=0.0;
405: }
406: /* B_delta */
407: if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */
408: vals_B_delta [partial_sum+n_neg_values-1]=-1.0;
409: }
410: if ( n_neg_values < j ) { /* there's a rank next to me to the right */
411: vals_B_delta [partial_sum+n_neg_values]=1.0;
412: }
413: /* scaling as in Klawonn-Widlund 1999 */
414: if (!pcbddc->use_deluxe_scaling) {
415: for (s=0;s<n_neg_values;s++) {
416: scalar_value = 0.0;
417: for (k=0;k<s+1;k++) scalar_value += array[k];
418: scaling_factors[partial_sum+s] = -scalar_value;
419: }
420: for (s=0;s<n_pos_values;s++) {
421: scalar_value = 0.0;
422: for (k=s+n_neg_values;k<j;k++) scalar_value += array[k];
423: scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
424: }
425: }
426: partial_sum += j;
427: }
428: cum += aux_local_numbering_2[i];
429: }
430: ISRestoreIndices(subset_n,&aux_global_numbering);
431: ISDestroy(&subset_mult);
432: ISDestroy(&subset_n);
433: PetscFree(aux_sums);
434: PetscFree(aux_local_numbering_1);
435: PetscFree(dual_dofs_boundary_indices);
436: if (all_factors) {
437: PetscFree(all_factors[0]);
438: PetscFree(all_factors);
439: }
441: /* Create local part of B_delta */
442: MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);
443: MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);
444: MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);
445: MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);
446: MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
447: for (i=0;i<n_local_lambda;i++) {
448: MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);
449: }
450: PetscFree(vals_B_delta);
451: MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);
452: MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);
454: BD1 = NULL;
455: BD2 = NULL;
456: if (fully_redundant) {
457: if (pcbddc->use_deluxe_scaling) SETERRQ(comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented");
458: MatCreate(PETSC_COMM_SELF,&ScalingMat);
459: MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);
460: MatSetType(ScalingMat,MATSEQAIJ);
461: MatSeqAIJSetPreallocation(ScalingMat,1,NULL);
462: for (i=0;i<n_local_lambda;i++) {
463: MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);
464: }
465: MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);
466: MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY);
467: MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);
468: MatDestroy(&ScalingMat);
469: } else {
470: MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);
471: MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);
472: if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
473: MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);
474: MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);
475: for (i=0;i<n_local_lambda;i++) {
476: MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);
477: }
478: MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);
479: MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);
480: } else {
481: /* scaling as in Klawonn-Widlund 1999 */
482: PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
483: PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
484: Mat T;
485: PetscScalar *W,lwork,*Bwork;
486: const PetscInt *idxs = NULL;
487: PetscInt cum,mss,*nnz;
488: PetscBLASInt *pivots,B_lwork,B_N,B_ierr;
490: if (!pcbddc->deluxe_singlemat) SETERRQ(comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
491: mss = 0;
492: PetscCalloc1(pcis->n_B,&nnz);
493: if (sub_schurs->is_Ej_all) {
494: ISGetIndices(sub_schurs->is_Ej_all,&idxs);
495: for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
496: PetscInt subset_size;
498: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
499: for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size;
500: mss = PetscMax(mss,subset_size);
501: cum += subset_size;
502: }
503: }
504: MatCreate(PETSC_COMM_SELF,&T);
505: MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);
506: MatSetType(T,MATSEQAIJ);
507: MatSeqAIJSetPreallocation(T,0,nnz);
508: PetscFree(nnz);
510: /* workspace allocation */
511: B_lwork = 0;
512: if (mss) {
513: PetscScalar dummy = 1;
515: B_lwork = -1;
516: PetscBLASIntCast(mss,&B_N);
517: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
518: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr));
519: PetscFPTrapPop();
520: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
521: PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);
522: }
523: PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork);
525: for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
526: PetscScalar *M;
527: PetscInt subset_size;
529: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
530: PetscBLASIntCast(subset_size,&B_N);
531: MatDenseGetArray(deluxe_ctx->seq_mat[i],&M);
532: PetscMemcpy(W,M,subset_size*subset_size*sizeof(PetscScalar));
533: MatDenseRestoreArray(deluxe_ctx->seq_mat[i],&M);
534: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
535: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr));
536: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
537: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
538: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
539: PetscFPTrapPop();
540: /* silent static analyzer */
541: if (!idxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"IDXS not present");
542: MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES);
543: cum += subset_size;
544: }
545: MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);
546: MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);
547: MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1);
548: MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2);
549: MatDestroy(&T);
550: PetscFree3(W,pivots,Bwork);
551: if (sub_schurs->is_Ej_all) {
552: ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);
553: }
554: }
555: }
556: PetscFree(scaling_factors);
557: PetscFree(cols_B_delta);
559: /* Layout of multipliers */
560: PetscLayoutCreate(comm,&llay);
561: PetscLayoutSetBlockSize(llay,1);
562: PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);
563: PetscLayoutSetUp(llay);
564: PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);
566: /* Local work vector of multipliers */
567: VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);
568: VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);
569: VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);
571: if (BD2) {
572: ISLocalToGlobalMapping l2g;
573: Mat T,TA,*pT;
574: IS is;
575: PetscInt nl,N;
576: BDdelta_DN ctx;
578: PetscLayoutGetLocalSize(llay,&nl);
579: PetscLayoutGetSize(llay,&N);
580: MatCreate(comm,&T);
581: MatSetSizes(T,nl,nl,N,N);
582: MatSetType(T,MATIS);
583: ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g);
584: MatSetLocalToGlobalMapping(T,l2g,l2g);
585: ISLocalToGlobalMappingDestroy(&l2g);
586: MatISSetLocalMat(T,BD2);
587: MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);
588: MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);
589: MatDestroy(&BD2);
590: MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&TA);
591: MatDestroy(&T);
592: ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is);
593: MatCreateSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT);
594: MatDestroy(&TA);
595: ISDestroy(&is);
596: BD2 = pT[0];
597: PetscFree(pT);
599: /* B_Ddelta for non-redundant multipliers with deluxe scaling */
600: PetscNew(&ctx);
601: MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL);
602: MatShellSetContext(fetidpmat_ctx->B_Ddelta,(void *)ctx);
603: MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred);
604: MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred);
605: MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred);
606: MatSetUp(fetidpmat_ctx->B_Ddelta);
608: PetscObjectReference((PetscObject)BD1);
609: ctx->BD = BD1;
610: KSPCreate(PETSC_COMM_SELF,&ctx->kBD);
611: KSPSetOperators(ctx->kBD,BD2,BD2);
612: VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work);
613: fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
614: }
615: MatDestroy(&BD1);
616: MatDestroy(&BD2);
618: /* fetidpmat sizes */
619: fetidpmat_ctx->n += nPgl;
620: fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg;
622: /* Global vector for FETI-DP linear system */
623: VecCreate(comm,&fetidp_global);
624: VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);
625: VecSetType(fetidp_global,VECMPI);
626: VecSetUp(fetidp_global);
628: /* Decide layout for fetidp dofs: if it is a saddle point problem
629: pressure is ordered first in the local part of the global vector
630: of the FETI-DP linear system */
631: if (nPg) {
632: Vec v;
633: IS IS_l2g_p,ais;
634: PetscLayout alay;
635: const PetscInt *idxs,*pranges,*aranges,*lranges;
636: PetscInt *l2g_indices_p,rst;
637: PetscMPIInt rank;
639: PetscMalloc1(nPl,&l2g_indices_p);
640: VecGetLayout(fetidp_global,&alay);
641: PetscLayoutGetRanges(alay,&aranges);
642: PetscLayoutGetRanges(play,&pranges);
643: PetscLayoutGetRanges(llay,&lranges);
645: MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global),&rank);
646: ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),pranges[rank+1]-pranges[rank],aranges[rank],1,&fetidpmat_ctx->pressure);
647: PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure,"F_P");
648: ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),lranges[rank+1]-lranges[rank],aranges[rank]+pranges[rank+1]-pranges[rank],1,&fetidpmat_ctx->lagrange);
649: PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange,"F_L");
650: ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);
651: /* shift local to global indices for pressure */
652: for (i=0;i<nPl;i++) {
653: PetscInt owner;
655: PetscLayoutFindOwner(play,idxs[i],&owner);
656: l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
657: }
658: ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);
659: ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);
661: /* local to global scatter for pressure */
662: VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);
663: ISDestroy(&IS_l2g_p);
665: /* scatter for lagrange multipliers only */
666: VecCreate(comm,&v);
667: VecSetType(v,VECSTANDARD);
668: VecSetLayout(v,llay);
669: VecSetUp(v);
670: ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&ais);
671: VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,v,ais,&fetidpmat_ctx->l2g_lambda_only);
672: ISDestroy(&ais);
673: VecDestroy(&v);
675: /* shift local to global indices for multipliers */
676: for (i=0;i<n_local_lambda;i++) {
677: PetscInt owner,ps;
679: PetscLayoutFindOwner(llay,l2g_indices[i],&owner);
680: ps = pranges[owner+1]-pranges[owner];
681: l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
682: }
684: /* scatter from alldofs to pressures global fetidp vector */
685: PetscLayoutGetRange(alay,&rst,NULL);
686: ISCreateStride(comm,nPgl,rst,1,&ais);
687: VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);
688: ISDestroy(&ais);
689: }
690: PetscLayoutDestroy(&llay);
691: PetscLayoutDestroy(&play);
692: ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);
694: /* scatter from local to global multipliers */
695: VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);
696: ISDestroy(&IS_l2g_lambda);
697: ISLocalToGlobalMappingDestroy(&l2gmap_p);
698: VecDestroy(&fetidp_global);
700: /* Create some work vectors needed by fetidp */
701: VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);
702: VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);
703: return(0);
704: }
707: PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
708: {
709: FETIDPMat_ctx mat_ctx;
710: PC_BDDC *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data;
711: PC_IS *pcis = (PC_IS*)fetidppc_ctx->pc->data;
712: PetscBool lumped = PETSC_FALSE;
716: MatShellGetContext(fetimat,(void**)&mat_ctx);
717: /* get references from objects created when setting up feti mat context */
718: PetscObjectReference((PetscObject)mat_ctx->lambda_local);
719: fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
720: PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);
721: fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
722: if (mat_ctx->deluxe_nonred) {
723: PC pc,mpc;
724: BDdelta_DN ctx;
725: MatSolverType solver;
726: const char *prefix;
728: MatShellGetContext(mat_ctx->B_Ddelta,&ctx);
729: KSPSetType(ctx->kBD,KSPPREONLY);
730: KSPGetPC(ctx->kBD,&mpc);
731: KSPGetPC(pcbddc->ksp_D,&pc);
732: PCSetType(mpc,PCLU);
733: PCFactorGetMatSolverType(pc,(MatSolverType*)&solver);
734: if (solver) {
735: PCFactorSetMatSolverType(mpc,solver);
736: }
737: MatGetOptionsPrefix(fetimat,&prefix);
738: KSPSetOptionsPrefix(ctx->kBD,prefix);
739: KSPAppendOptionsPrefix(ctx->kBD,"bddelta_");
740: KSPSetFromOptions(ctx->kBD);
741: }
743: if (mat_ctx->l2g_lambda_only) {
744: PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only);
745: fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
746: } else {
747: PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);
748: fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
749: }
750: /* Dirichlet preconditioner */
751: PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL);
752: if (!lumped) {
753: IS iV;
754: PetscBool discrete_harmonic = PETSC_FALSE;
756: PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iV",(PetscObject*)&iV);
757: if (iV) {
758: PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL);
759: }
760: if (discrete_harmonic) {
761: KSP sksp;
762: PC pc;
763: PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
764: Mat A_II,A_IB,A_BI;
765: IS iP = NULL;
766: PetscBool isshell,reuse = PETSC_FALSE;
767: KSPType ksptype;
768: const char *prefix;
770: /*
771: We constructs a Schur complement for
773: | A_II A_ID |
774: | A_DI A_DD |
776: instead of
778: | A_II B^t_II A_ID |
779: | B_II -C_II B_ID |
780: | A_DI B^t_ID A_DD |
782: */
783: if (sub_schurs && sub_schurs->reuse_solver) {
784: PetscObjectQuery((PetscObject)sub_schurs->A,"__KSPFETIDP_iP",(PetscObject*)&iP);
785: if (iP) reuse = PETSC_TRUE;
786: }
787: if (!reuse) {
788: IS aB;
789: PetscInt nb;
790: ISGetLocalSize(pcis->is_B_local,&nb);
791: ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB);
792: MatCreateSubMatrix(pcis->A_II,iV,iV,MAT_INITIAL_MATRIX,&A_II);
793: MatCreateSubMatrix(pcis->A_IB,iV,aB,MAT_INITIAL_MATRIX,&A_IB);
794: MatCreateSubMatrix(pcis->A_BI,aB,iV,MAT_INITIAL_MATRIX,&A_BI);
795: ISDestroy(&aB);
796: } else {
797: MatCreateSubMatrix(sub_schurs->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&A_IB);
798: MatCreateSubMatrix(sub_schurs->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&A_BI);
799: PetscObjectReference((PetscObject)pcis->A_II);
800: A_II = pcis->A_II;
801: }
802: MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j);
804: /* propagate settings of solver */
805: MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp);
806: KSPGetType(pcis->ksp_D,&ksptype);
807: KSPSetType(sksp,ksptype);
808: KSPGetPC(pcis->ksp_D,&pc);
809: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell);
810: if (!isshell) {
811: MatSolverType solver;
812: PCType pctype;
814: PCGetType(pc,&pctype);
815: PCFactorGetMatSolverType(pc,(MatSolverType*)&solver);
816: KSPGetPC(sksp,&pc);
817: PCSetType(pc,pctype);
818: if (solver) {
819: PCFactorSetMatSolverType(pc,solver);
820: }
821: } else {
822: KSPGetPC(sksp,&pc);
823: PCSetType(pc,PCLU);
824: }
825: MatDestroy(&A_II);
826: MatDestroy(&A_IB);
827: MatDestroy(&A_BI);
828: MatGetOptionsPrefix(fetimat,&prefix);
829: KSPSetOptionsPrefix(sksp,prefix);
830: KSPAppendOptionsPrefix(sksp,"harmonic_");
831: KSPSetFromOptions(sksp);
832: if (reuse) {
833: KSPSetPC(sksp,sub_schurs->reuse_solver->interior_solver);
834: PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver,(PetscObject)sksp,0);
835: }
836: } else { /* default Dirichlet preconditioner is pde-harmonic */
837: MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);
838: MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);
839: }
840: } else {
841: PetscObjectReference((PetscObject)pcis->A_BB);
842: fetidppc_ctx->S_j = pcis->A_BB;
843: }
844: /* saddle-point */
845: if (mat_ctx->xPg) {
846: PetscObjectReference((PetscObject)mat_ctx->xPg);
847: fetidppc_ctx->xPg = mat_ctx->xPg;
848: PetscObjectReference((PetscObject)mat_ctx->yPg);
849: fetidppc_ctx->yPg = mat_ctx->yPg;
850: }
851: return(0);
852: }
854: PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
855: {
856: FETIDPMat_ctx mat_ctx;
857: PC_BDDC *pcbddc;
858: PC_IS *pcis;
862: MatShellGetContext(fetimat,(void**)&mat_ctx);
863: pcis = (PC_IS*)mat_ctx->pc->data;
864: pcbddc = (PC_BDDC*)mat_ctx->pc->data;
865: /* Application of B_delta^T */
866: VecSet(pcis->vec1_B,0.);
867: VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
868: VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
869: MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);
871: /* Add contribution from saddle point */
872: if (mat_ctx->l2g_p) {
873: VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
874: VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
875: if (pcbddc->switch_static) {
876: if (trans) {
877: MatMultTranspose(mat_ctx->B_BI,mat_ctx->vP,pcis->vec1_D);
878: } else {
879: MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);
880: }
881: }
882: if (trans) {
883: MatMultTransposeAdd(mat_ctx->B_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);
884: } else {
885: MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);
886: }
887: } else {
888: if (pcbddc->switch_static) {
889: VecSet(pcis->vec1_D,0.0);
890: }
891: }
892: /* Application of \widetilde{S}^-1 */
893: PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
894: PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,trans);
895: PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));
896: VecSet(y,0.0);
897: /* Application of B_delta */
898: MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);
899: /* Contribution from boundary pressures */
900: if (mat_ctx->C) {
901: const PetscScalar *lx;
902: PetscScalar *ly;
904: /* pressure ordered first in the local part of x and y */
905: VecGetArrayRead(x,&lx);
906: VecGetArray(y,&ly);
907: VecPlaceArray(mat_ctx->xPg,lx);
908: VecPlaceArray(mat_ctx->yPg,ly);
909: if (trans) {
910: MatMultTranspose(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);
911: } else {
912: MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);
913: }
914: VecResetArray(mat_ctx->xPg);
915: VecResetArray(mat_ctx->yPg);
916: VecRestoreArrayRead(x,&lx);
917: VecRestoreArray(y,&ly);
918: }
919: /* Add contribution from saddle point */
920: if (mat_ctx->l2g_p) {
921: if (trans) {
922: MatMultTranspose(mat_ctx->Bt_BB,pcis->vec1_B,mat_ctx->vP);
923: } else {
924: MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);
925: }
926: if (pcbddc->switch_static) {
927: if (trans) {
928: MatMultTransposeAdd(mat_ctx->Bt_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);
929: } else {
930: MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);
931: }
932: }
933: VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);
934: VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);
935: }
936: VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
937: VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
938: return(0);
939: }
941: PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
942: {
946: FETIDPMatMult_Kernel(fetimat,x,y,PETSC_FALSE);
947: return(0);
948: }
950: PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
951: {
955: FETIDPMatMult_Kernel(fetimat,x,y,PETSC_TRUE);
956: return(0);
957: }
959: PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
960: {
961: FETIDPPC_ctx pc_ctx;
962: PC_IS *pcis;
966: PCShellGetContext(fetipc,(void**)&pc_ctx);
967: pcis = (PC_IS*)pc_ctx->pc->data;
968: /* Application of B_Ddelta^T */
969: VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
970: VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
971: VecSet(pcis->vec2_B,0.0);
972: MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);
973: /* Application of local Schur complement */
974: if (trans) {
975: MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);
976: } else {
977: MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);
978: }
979: /* Application of B_Ddelta */
980: MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);
981: VecSet(y,0.0);
982: VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
983: VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);
984: return(0);
985: }
987: PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
988: {
992: FETIDPPCApply_Kernel(pc,x,y,PETSC_FALSE);
993: return(0);
994: }
996: PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
997: {
1001: FETIDPPCApply_Kernel(pc,x,y,PETSC_TRUE);
1002: return(0);
1003: }
1005: PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
1006: {
1007: FETIDPPC_ctx pc_ctx;
1008: PetscBool iascii;
1009: PetscViewer sviewer;
1010: PetscErrorCode ierr;
1013: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1014: if (iascii) {
1015: PetscMPIInt rank;
1016: PetscBool isschur,isshell;
1018: PCShellGetContext(pc,(void**)&pc_ctx);
1019: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
1020: PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur);
1021: if (isschur) {
1022: PetscViewerASCIIPrintf(viewer," Dirichlet preconditioner (just from rank 0)\n");
1023: } else {
1024: PetscViewerASCIIPrintf(viewer," Lumped preconditioner (just from rank 0)\n");
1025: }
1026: PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);
1027: if (!rank) {
1028: PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);
1029: PetscViewerASCIIPushTab(sviewer);
1030: MatView(pc_ctx->S_j,sviewer);
1031: PetscViewerASCIIPopTab(sviewer);
1032: PetscViewerPopFormat(sviewer);
1033: }
1034: PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell);
1035: if (isshell) {
1036: BDdelta_DN ctx;
1037: PetscViewerASCIIPrintf(viewer," FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n");
1038: MatShellGetContext(pc_ctx->B_Ddelta,&ctx);
1039: if (!rank) {
1040: PetscInt tl;
1042: PetscViewerASCIIGetTab(sviewer,&tl);
1043: PetscObjectSetTabLevel((PetscObject)ctx->kBD,tl);
1044: KSPView(ctx->kBD,sviewer);
1045: PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);
1046: MatView(ctx->BD,sviewer);
1047: PetscViewerPopFormat(sviewer);
1048: }
1049: }
1050: PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);
1051: PetscViewerFlush(viewer);
1052: }
1053: return(0);
1054: }