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