Actual source code: bddcschurs.c
petsc-3.8.3 2017-12-09
1: #include <../src/ksp/pc/impls/bddc/bddc.h>
2: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3: #include <../src/mat/impls/dense/seq/dense.h>
4: #include <petscblaslapack.h>
6: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
7: static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
8: static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
9: static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
11: /* if v2 is not present, correction is done in-place */
12: PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
13: {
14: PetscScalar *array;
15: PetscScalar *array2;
19: if (!ctx->benign_n) return(0);
20: if (sol && full) {
21: PetscInt n_I,size_schur;
23: /* get sizes */
24: MatGetSize(ctx->benign_csAIB,&size_schur,NULL);
25: VecGetSize(v,&n_I);
26: n_I = n_I - size_schur;
27: /* get schur sol from array */
28: VecGetArray(v,&array);
29: VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);
30: VecRestoreArray(v,&array);
31: /* apply interior sol correction */
32: MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);
33: VecResetArray(ctx->benign_dummy_schur_vec);
34: MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);
35: }
36: if (v2) {
37: PetscInt nl;
39: VecGetArrayRead(v,(const PetscScalar**)&array);
40: VecGetLocalSize(v2,&nl);
41: VecGetArray(v2,&array2);
42: PetscMemcpy(array2,array,nl*sizeof(PetscScalar));
43: } else {
44: VecGetArray(v,&array);
45: array2 = array;
46: }
47: if (!sol) { /* change rhs */
48: PetscInt n;
49: for (n=0;n<ctx->benign_n;n++) {
50: PetscScalar sum = 0.;
51: const PetscInt *cols;
52: PetscInt nz,i;
54: ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);
55: ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);
56: for (i=0;i<nz-1;i++) sum += array[cols[i]];
57: #if defined(PETSC_USE_COMPLEX)
58: sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
59: #else
60: sum = -sum/nz;
61: #endif
62: for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
63: ctx->benign_save_vals[n] = array2[cols[nz-1]];
64: array2[cols[nz-1]] = sum;
65: ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);
66: }
67: } else {
68: PetscInt n;
69: for (n=0;n<ctx->benign_n;n++) {
70: PetscScalar sum = 0.;
71: const PetscInt *cols;
72: PetscInt nz,i;
73: ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);
74: ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);
75: for (i=0;i<nz-1;i++) sum += array[cols[i]];
76: #if defined(PETSC_USE_COMPLEX)
77: sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
78: #else
79: sum = -sum/nz;
80: #endif
81: for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
82: array2[cols[nz-1]] = ctx->benign_save_vals[n];
83: ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);
84: }
85: }
86: if (v2) {
87: VecRestoreArrayRead(v,(const PetscScalar**)&array);
88: VecRestoreArray(v2,&array2);
89: } else {
90: VecRestoreArray(v,&array);
91: }
92: if (!sol && full) {
93: Vec usedv;
94: PetscInt n_I,size_schur;
96: /* get sizes */
97: MatGetSize(ctx->benign_csAIB,&size_schur,NULL);
98: VecGetSize(v,&n_I);
99: n_I = n_I - size_schur;
100: /* compute schur rhs correction */
101: if (v2) {
102: usedv = v2;
103: } else {
104: usedv = v;
105: }
106: /* apply schur rhs correction */
107: MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);
108: VecGetArrayRead(usedv,(const PetscScalar**)&array);
109: VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);
110: VecRestoreArrayRead(usedv,(const PetscScalar**)&array);
111: MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);
112: VecResetArray(ctx->benign_dummy_schur_vec);
113: }
114: return(0);
115: }
117: static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
118: {
119: PCBDDCReuseSolvers ctx;
120: PetscBool copy = PETSC_FALSE;
121: PetscErrorCode ierr;
124: PCShellGetContext(pc,(void **)&ctx);
125: if (full) {
126: #if defined(PETSC_HAVE_MUMPS)
127: MatMumpsSetIcntl(ctx->F,26,-1);
128: #endif
129: #if defined(PETSC_HAVE_MKL_PARDISO)
130: MatMkl_PardisoSetCntl(ctx->F,70,0);
131: #endif
132: copy = ctx->has_vertices;
133: } else { /* interior solver */
134: #if defined(PETSC_HAVE_MUMPS)
135: MatMumpsSetIcntl(ctx->F,26,0);
136: #endif
137: #if defined(PETSC_HAVE_MKL_PARDISO)
138: MatMkl_PardisoSetCntl(ctx->F,70,1);
139: #endif
140: copy = PETSC_TRUE;
141: }
142: /* copy rhs into factored matrix workspace */
143: if (copy) {
144: PetscInt n;
145: PetscScalar *array,*array_solver;
147: VecGetLocalSize(rhs,&n);
148: VecGetArrayRead(rhs,(const PetscScalar**)&array);
149: VecGetArray(ctx->rhs,&array_solver);
150: PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));
151: VecRestoreArray(ctx->rhs,&array_solver);
152: VecRestoreArrayRead(rhs,(const PetscScalar**)&array);
154: PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);
155: if (transpose) {
156: MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);
157: } else {
158: MatSolve(ctx->F,ctx->rhs,ctx->sol);
159: }
160: PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);
162: /* get back data to caller worskpace */
163: VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);
164: VecGetArray(sol,&array);
165: PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));
166: VecRestoreArray(sol,&array);
167: VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);
168: } else {
169: if (ctx->benign_n) {
170: PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);
171: if (transpose) {
172: MatSolveTranspose(ctx->F,ctx->rhs,sol);
173: } else {
174: MatSolve(ctx->F,ctx->rhs,sol);
175: }
176: PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);
177: } else {
178: if (transpose) {
179: MatSolveTranspose(ctx->F,rhs,sol);
180: } else {
181: MatSolve(ctx->F,rhs,sol);
182: }
183: }
184: }
185: /* restore defaults */
186: #if defined(PETSC_HAVE_MUMPS)
187: MatMumpsSetIcntl(ctx->F,26,-1);
188: #endif
189: #if defined(PETSC_HAVE_MKL_PARDISO)
190: MatMkl_PardisoSetCntl(ctx->F,70,0);
191: #endif
192: return(0);
193: }
195: static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
196: {
197: PetscErrorCode ierr;
200: PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);
201: return(0);
202: }
204: static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
205: {
206: PetscErrorCode ierr;
209: PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);
210: return(0);
211: }
213: static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
214: {
215: PetscErrorCode ierr;
218: PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);
219: return(0);
220: }
222: static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
223: {
224: PetscErrorCode ierr;
227: PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);
228: return(0);
229: }
231: static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
232: {
233: PetscInt i;
237: MatDestroy(&reuse->F);
238: VecDestroy(&reuse->sol);
239: VecDestroy(&reuse->rhs);
240: PCDestroy(&reuse->interior_solver);
241: PCDestroy(&reuse->correction_solver);
242: ISDestroy(&reuse->is_R);
243: ISDestroy(&reuse->is_B);
244: VecScatterDestroy(&reuse->correction_scatter_B);
245: VecDestroy(&reuse->sol_B);
246: VecDestroy(&reuse->rhs_B);
247: for (i=0;i<reuse->benign_n;i++) {
248: ISDestroy(&reuse->benign_zerodiag_subs[i]);
249: }
250: PetscFree(reuse->benign_zerodiag_subs);
251: PetscFree(reuse->benign_save_vals);
252: MatDestroy(&reuse->benign_csAIB);
253: MatDestroy(&reuse->benign_AIIm1ones);
254: VecDestroy(&reuse->benign_corr_work);
255: VecDestroy(&reuse->benign_dummy_schur_vec);
256: return(0);
257: }
259: static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
260: {
261: Mat B, C, D, Bd, Cd, AinvBd;
262: KSP ksp;
263: PC pc;
264: PetscBool isLU, isILU, isCHOL, Bdense, Cdense;
265: PetscReal fill = 2.0;
266: PetscInt n_I;
267: PetscMPIInt size;
271: MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);
272: if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
273: if (reuse == MAT_REUSE_MATRIX) {
274: PetscBool Sdense;
276: PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);
277: if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
278: }
279: MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);
280: MatSchurComplementGetKSP(M, &ksp);
281: KSPGetPC(ksp, &pc);
282: PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);
283: PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);
284: PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);
285: PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);
286: PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);
287: MatGetSize(B,&n_I,NULL);
288: if (n_I) {
289: if (!Bdense) {
290: MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);
291: } else {
292: Bd = B;
293: }
295: if (isLU || isILU || isCHOL) {
296: Mat fact;
297: KSPSetUp(ksp);
298: PCFactorGetMatrix(pc, &fact);
299: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
300: MatMatSolve(fact, Bd, AinvBd);
301: } else {
302: PetscBool ex = PETSC_TRUE;
304: if (ex) {
305: Mat Ainvd;
307: PCComputeExplicitOperator(pc, &Ainvd);
308: MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);
309: MatDestroy(&Ainvd);
310: } else {
311: Vec sol,rhs;
312: PetscScalar *arrayrhs,*arraysol;
313: PetscInt i,nrhs,n;
315: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
316: MatGetSize(Bd,&n,&nrhs);
317: MatDenseGetArray(Bd,&arrayrhs);
318: MatDenseGetArray(AinvBd,&arraysol);
319: KSPGetSolution(ksp,&sol);
320: KSPGetRhs(ksp,&rhs);
321: for (i=0;i<nrhs;i++) {
322: VecPlaceArray(rhs,arrayrhs+i*n);
323: VecPlaceArray(sol,arraysol+i*n);
324: KSPSolve(ksp,rhs,sol);
325: VecResetArray(rhs);
326: VecResetArray(sol);
327: }
328: MatDenseRestoreArray(Bd,&arrayrhs);
329: MatDenseRestoreArray(AinvBd,&arrayrhs);
330: }
331: }
332: if (!Bdense & !issym) {
333: MatDestroy(&Bd);
334: }
336: if (!issym) {
337: if (!Cdense) {
338: MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);
339: } else {
340: Cd = C;
341: }
342: MatMatMult(Cd, AinvBd, reuse, fill, S);
343: if (!Cdense) {
344: MatDestroy(&Cd);
345: }
346: } else {
347: MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);
348: if (!Bdense) {
349: MatDestroy(&Bd);
350: }
351: }
352: MatDestroy(&AinvBd);
353: }
355: if (D) {
356: Mat Dd;
357: PetscBool Ddense;
359: PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);
360: if (!Ddense) {
361: MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);
362: } else {
363: Dd = D;
364: }
365: if (n_I) {
366: MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);
367: } else {
368: if (reuse == MAT_INITIAL_MATRIX) {
369: MatDuplicate(Dd,MAT_COPY_VALUES,S);
370: } else {
371: MatCopy(Dd,*S,SAME_NONZERO_PATTERN);
372: }
373: }
374: if (!Ddense) {
375: MatDestroy(&Dd);
376: }
377: } else {
378: MatScale(*S,-1.0);
379: }
380: return(0);
381: }
383: PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal)
384: {
385: Mat F,A_II,A_IB,A_BI,A_BB,AE_II;
386: Mat S_all;
387: Mat global_schur_subsets,work_mat,*submats;
388: ISLocalToGlobalMapping l2gmap_subsets;
389: IS is_I,is_I_layer;
390: IS all_subsets,all_subsets_mult,all_subsets_n;
391: PetscInt *nnz,*all_local_idx_N;
392: PetscInt *auxnum1,*auxnum2;
393: PetscInt i,subset_size,max_subset_size;
394: PetscInt n_B,extra,local_size,global_size;
395: PetscBLASInt B_N,B_ierr,B_lwork,*pivots;
396: PetscScalar *Bwork;
397: PetscSubcomm subcomm;
398: PetscMPIInt color,rank;
399: MPI_Comm comm_n;
400: PetscBool deluxe = PETSC_TRUE, use_getr = PETSC_FALSE;
401: PetscErrorCode ierr;
404: MatDestroy(&sub_schurs->A);
405: MatDestroy(&sub_schurs->S);
407: sub_schurs->is_hermitian = PETSC_TRUE;
408: sub_schurs->is_posdef = PETSC_TRUE;
409: if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
410: PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
411: PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);
412: PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);
413: PetscOptionsEnd();
415: /* convert matrix if needed */
416: if (Ain) {
417: PetscBool isseqaij;
418: PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);
419: if (isseqaij) {
420: PetscObjectReference((PetscObject)Ain);
421: sub_schurs->A = Ain;
422: } else {
423: MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);
424: }
425: }
427: PetscObjectReference((PetscObject)Sin);
428: sub_schurs->S = Sin;
429: if (sub_schurs->schur_explicit) {
430: sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
431: }
433: /* preliminary checks */
434: if (!sub_schurs->schur_explicit && compute_Stilda) SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO");
436: /* restrict work on active processes */
437: color = 0;
438: if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
439: MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
440: PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);
441: PetscSubcommSetNumber(subcomm,2);
442: PetscSubcommSetTypeGeneral(subcomm,color,rank);
443: PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);
444: PetscSubcommDestroy(&subcomm);
445: if (!sub_schurs->n_subs) {
446: PetscCommDestroy(&comm_n);
447: return(0);
448: }
449: /* PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL); */
451: /* get Schur complement matrices */
452: if (!sub_schurs->schur_explicit) {
453: Mat tA_IB,tA_BI,tA_BB;
454: PetscBool isseqsbaij;
455: MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);
456: PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);
457: if (isseqsbaij) {
458: MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);
459: MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);
460: MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);
461: } else {
462: PetscObjectReference((PetscObject)tA_BB);
463: A_BB = tA_BB;
464: PetscObjectReference((PetscObject)tA_IB);
465: A_IB = tA_IB;
466: PetscObjectReference((PetscObject)tA_BI);
467: A_BI = tA_BI;
468: }
469: } else {
470: A_II = NULL;
471: A_IB = NULL;
472: A_BI = NULL;
473: A_BB = NULL;
474: }
475: S_all = NULL;
477: /* determine interior problems */
478: ISGetLocalSize(sub_schurs->is_I,&i);
479: if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
480: PetscBT touched;
481: const PetscInt* idx_B;
482: PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
484: if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
485: /* get sizes */
486: ISGetLocalSize(sub_schurs->is_I,&n_I);
487: ISGetLocalSize(sub_schurs->is_B,&n_B);
489: PetscMalloc1(n_I+n_B,&local_numbering);
490: PetscBTCreate(n_I+n_B,&touched);
491: PetscBTMemzero(n_I+n_B,touched);
493: /* all boundary dofs must be skipped when adding layers */
494: ISGetIndices(sub_schurs->is_B,&idx_B);
495: for (j=0;j<n_B;j++) {
496: PetscBTSet(touched,idx_B[j]);
497: }
498: PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));
499: ISRestoreIndices(sub_schurs->is_B,&idx_B);
501: /* add prescribed number of layers of dofs */
502: n_local_dofs = n_B;
503: n_prev_added = n_B;
504: for (layer=0;layer<nlayers;layer++) {
505: PetscInt n_added;
506: if (n_local_dofs == n_I+n_B) break;
507: if (n_local_dofs > n_I+n_B) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error querying layer %D. Out of bound access (%D > %D)",layer,n_local_dofs,n_I+n_B);
508: PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);
509: n_prev_added = n_added;
510: n_local_dofs += n_added;
511: if (!n_added) break;
512: }
513: PetscBTDestroy(&touched);
515: /* IS for I layer dofs in original numbering */
516: ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);
517: PetscFree(local_numbering);
518: ISSort(is_I_layer);
519: /* IS for I layer dofs in I numbering */
520: if (!sub_schurs->schur_explicit) {
521: ISLocalToGlobalMapping ItoNmap;
522: ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);
523: ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);
524: ISLocalToGlobalMappingDestroy(&ItoNmap);
526: /* II block */
527: MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);
528: }
529: } else {
530: PetscInt n_I;
532: /* IS for I dofs in original numbering */
533: PetscObjectReference((PetscObject)sub_schurs->is_I);
534: is_I_layer = sub_schurs->is_I;
536: /* IS for I dofs in I numbering (strided 1) */
537: if (!sub_schurs->schur_explicit) {
538: ISGetSize(sub_schurs->is_I,&n_I);
539: ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);
541: /* II block is the same */
542: PetscObjectReference((PetscObject)A_II);
543: AE_II = A_II;
544: }
545: }
547: /* Get info on subset sizes and sum of all subsets sizes */
548: max_subset_size = 0;
549: local_size = 0;
550: for (i=0;i<sub_schurs->n_subs;i++) {
551: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
552: max_subset_size = PetscMax(subset_size,max_subset_size);
553: local_size += subset_size;
554: }
556: /* Work arrays for local indices */
557: extra = 0;
558: ISGetLocalSize(sub_schurs->is_B,&n_B);
559: if (sub_schurs->schur_explicit && is_I_layer) {
560: ISGetLocalSize(is_I_layer,&extra);
561: }
562: PetscMalloc1(n_B+extra,&all_local_idx_N);
563: if (extra) {
564: const PetscInt *idxs;
565: ISGetIndices(is_I_layer,&idxs);
566: PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));
567: ISRestoreIndices(is_I_layer,&idxs);
568: }
569: PetscMalloc1(local_size,&nnz);
570: PetscMalloc1(sub_schurs->n_subs,&auxnum1);
571: PetscMalloc1(sub_schurs->n_subs,&auxnum2);
573: /* Get local indices in local numbering */
574: local_size = 0;
575: for (i=0;i<sub_schurs->n_subs;i++) {
576: PetscInt j;
577: const PetscInt *idxs;
579: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
580: ISGetIndices(sub_schurs->is_subs[i],&idxs);
581: /* start (smallest in global ordering) and multiplicity */
582: auxnum1[i] = idxs[0];
583: auxnum2[i] = subset_size;
584: /* subset indices in local numbering */
585: PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));
586: ISRestoreIndices(sub_schurs->is_subs[i],&idxs);
587: for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
588: local_size += subset_size;
589: }
591: /* allocate extra workspace needed only for GETRI */
592: if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
593: PetscScalar lwork,dummyscalar = 0.;
594: PetscBLASInt dummyint = 0;
596: use_getr = PETSC_TRUE;
597: B_lwork = -1;
598: PetscBLASIntCast(local_size,&B_N);
599: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
600: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
601: PetscFPTrapPop();
602: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
603: PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);
604: PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);
605: } else {
606: Bwork = NULL;
607: pivots = NULL;
608: }
610: /* prepare parallel matrices for summing up properly schurs on subsets */
611: ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);
612: ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);
613: ISDestroy(&all_subsets_n);
614: ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);
615: ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);
616: ISDestroy(&all_subsets);
617: ISDestroy(&all_subsets_mult);
618: ISGetLocalSize(all_subsets_n,&i);
619: if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
620: ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);
621: MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);
622: ISLocalToGlobalMappingDestroy(&l2gmap_subsets);
623: MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);
624: MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);
625: MatSetType(global_schur_subsets,MATMPIAIJ);
627: /* subset indices in local boundary numbering */
628: if (!sub_schurs->is_Ej_all) {
629: PetscInt *all_local_idx_B;
631: PetscMalloc1(local_size,&all_local_idx_B);
632: ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);
633: if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D\n",subset_size,local_size);
634: ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);
635: }
637: if (change) {
638: ISLocalToGlobalMapping BtoS;
639: IS change_primal_B;
640: IS change_primal_all;
642: if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
643: if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
644: PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);
645: for (i=0;i<sub_schurs->n_subs;i++) {
646: ISLocalToGlobalMapping NtoS;
647: ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);
648: ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);
649: ISLocalToGlobalMappingDestroy(&NtoS);
650: }
651: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);
652: ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);
653: ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);
654: ISLocalToGlobalMappingDestroy(&BtoS);
655: ISDestroy(&change_primal_B);
656: PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);
657: for (i=0;i<sub_schurs->n_subs;i++) {
658: Mat change_sub;
660: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
661: KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);
662: KSPSetType(sub_schurs->change[i],KSPPREONLY);
663: if (!sub_schurs->change_with_qr) {
664: MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);
665: } else {
666: Mat change_subt;
667: MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);
668: MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);
669: MatDestroy(&change_subt);
670: }
671: KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);
672: MatDestroy(&change_sub);
673: KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);
674: KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");
675: }
676: ISDestroy(&change_primal_all);
677: }
679: /* Local matrix of all local Schur on subsets (transposed) */
680: if (!sub_schurs->S_Ej_all) {
681: MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);
682: MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);
683: MatSetType(sub_schurs->S_Ej_all,MATAIJ);
684: MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);
685: }
687: /* Compute Schur complements explicitly */
688: F = NULL;
689: if (!sub_schurs->schur_explicit) { /* this code branch is used when MatFactor with Schur complemnent support is not present; it is not very efficient, unless the economic version of the scaling is required */
690: Mat S_Ej_expl;
691: PetscScalar *work;
692: PetscInt j,*dummy_idx;
693: PetscBool Sdense;
695: PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
696: local_size = 0;
697: for (i=0;i<sub_schurs->n_subs;i++) {
698: IS is_subset_B;
699: Mat AE_EE,AE_IE,AE_EI,S_Ej;
701: /* subsets in original and boundary numbering */
702: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);
703: /* EE block */
704: MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);
705: /* IE block */
706: MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);
707: /* EI block */
708: if (sub_schurs->is_hermitian) {
709: MatCreateTranspose(AE_IE,&AE_EI);
710: } else {
711: MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);
712: }
713: ISDestroy(&is_subset_B);
714: MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);
715: MatDestroy(&AE_EE);
716: MatDestroy(&AE_IE);
717: MatDestroy(&AE_EI);
718: if (AE_II == A_II) { /* we can reuse the same ksp */
719: KSP ksp;
720: MatSchurComplementGetKSP(sub_schurs->S,&ksp);
721: MatSchurComplementSetKSP(S_Ej,ksp);
722: } else { /* build new ksp object which inherits ksp and pc types from the original one */
723: KSP origksp,schurksp;
724: PC origpc,schurpc;
725: KSPType ksp_type;
726: PetscInt n_internal;
727: PetscBool ispcnone;
729: MatSchurComplementGetKSP(sub_schurs->S,&origksp);
730: MatSchurComplementGetKSP(S_Ej,&schurksp);
731: KSPGetType(origksp,&ksp_type);
732: KSPSetType(schurksp,ksp_type);
733: KSPGetPC(schurksp,&schurpc);
734: KSPGetPC(origksp,&origpc);
735: PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);
736: if (!ispcnone) {
737: PCType pc_type;
738: PCGetType(origpc,&pc_type);
739: PCSetType(schurpc,pc_type);
740: } else {
741: PCSetType(schurpc,PCLU);
742: }
743: ISGetSize(is_I,&n_internal);
744: if (n_internal) { /* UMFPACK gives error with 0 sized problems */
745: MatSolverPackage solver=NULL;
746: PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);
747: if (solver) {
748: PCFactorSetMatSolverPackage(schurpc,solver);
749: }
750: }
751: KSPSetUp(schurksp);
752: }
753: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
754: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);
755: PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);
756: PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);
757: if (Sdense) {
758: for (j=0;j<subset_size;j++) {
759: dummy_idx[j]=local_size+j;
760: }
761: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
762: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
763: MatDestroy(&S_Ej);
764: MatDestroy(&S_Ej_expl);
765: local_size += subset_size;
766: }
767: PetscFree2(dummy_idx,work);
768: /* free */
769: ISDestroy(&is_I);
770: MatDestroy(&AE_II);
771: PetscFree(all_local_idx_N);
772: } else {
773: Mat A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
774: Vec Dall = NULL;
775: IS is_A_all,*is_p_r = NULL;
776: PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
777: PetscInt n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
778: PetscBool economic,solver_S,S_lower_triangular = PETSC_FALSE;
779: PetscBool schur_has_vertices,factor_workaround;
780: char solver[256];
782: /* get sizes */
783: n_I = 0;
784: if (is_I_layer) {
785: ISGetLocalSize(is_I_layer,&n_I);
786: }
787: economic = PETSC_FALSE;
788: ISGetLocalSize(sub_schurs->is_I,&cum);
789: if (cum != n_I) economic = PETSC_TRUE;
790: MatGetLocalSize(sub_schurs->A,&n,NULL);
791: size_active_schur = local_size;
793: /* import scaling vector (wrong formulation if we have 3D edges) */
794: if (scaling && compute_Stilda) {
795: const PetscScalar *array;
796: PetscScalar *array2;
797: const PetscInt *idxs;
798: PetscInt i;
800: ISGetIndices(sub_schurs->is_Ej_all,&idxs);
801: VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);
802: VecGetArrayRead(scaling,&array);
803: VecGetArray(Dall,&array2);
804: for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
805: VecRestoreArray(Dall,&array2);
806: VecRestoreArrayRead(scaling,&array);
807: ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);
808: deluxe = PETSC_FALSE;
809: }
811: /* size active schurs does not count any dirichlet or vertex dof on the interface */
812: factor_workaround = PETSC_FALSE;
813: schur_has_vertices = PETSC_FALSE;
814: cum = n_I+size_active_schur;
815: if (sub_schurs->is_dir) {
816: const PetscInt* idxs;
817: PetscInt n_dir;
819: ISGetLocalSize(sub_schurs->is_dir,&n_dir);
820: ISGetIndices(sub_schurs->is_dir,&idxs);
821: PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));
822: ISRestoreIndices(sub_schurs->is_dir,&idxs);
823: cum += n_dir;
824: factor_workaround = PETSC_TRUE;
825: }
826: /* include the primal vertices in the Schur complement */
827: if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
828: PetscInt n_v;
830: ISGetLocalSize(sub_schurs->is_vertices,&n_v);
831: if (n_v) {
832: const PetscInt* idxs;
834: ISGetIndices(sub_schurs->is_vertices,&idxs);
835: PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));
836: ISRestoreIndices(sub_schurs->is_vertices,&idxs);
837: cum += n_v;
838: factor_workaround = PETSC_TRUE;
839: schur_has_vertices = PETSC_TRUE;
840: }
841: }
842: size_schur = cum - n_I;
843: ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);
844: /* get working mat (TODO: factorize without actually permuting) */
845: if (cum == n) {
846: ISSetPermutation(is_A_all);
847: MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);
848: } else {
849: MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);
850: }
851: MatSetOptionsPrefix(A,sub_schurs->prefix);
852: MatAppendOptionsPrefix(A,"sub_schurs_");
854: /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
855: n^2 instead of n^1.5 or something. This is a workaround */
856: if (benign_n) {
857: Vec v;
858: ISLocalToGlobalMapping N_to_reor;
859: IS is_p0,is_p0_p;
860: PetscScalar *cs_AIB,*AIIm1_data;
861: PetscInt sizeA;
863: ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);
864: ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);
865: ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);
866: ISDestroy(&is_p0);
867: MatCreateVecs(A,&v,NULL);
868: VecGetSize(v,&sizeA);
869: MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);
870: MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);
871: MatDenseGetArray(cs_AIB_mat,&cs_AIB);
872: MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
873: PetscMalloc1(benign_n,&is_p_r);
874: /* compute colsum of A_IB restricted to pressures */
875: for (i=0;i<benign_n;i++) {
876: Vec benign_AIIm1_ones;
877: PetscScalar *array;
878: const PetscInt *idxs;
879: PetscInt j,nz;
881: ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);
882: ISGetLocalSize(is_p_r[i],&nz);
883: ISGetIndices(is_p_r[i],&idxs);
884: for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
885: ISRestoreIndices(is_p_r[i],&idxs);
886: VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);
887: MatMult(A,benign_AIIm1_ones,v);
888: VecDestroy(&benign_AIIm1_ones);
889: VecGetArray(v,&array);
890: for (j=0;j<size_schur;j++) {
891: #if defined(PETSC_USE_COMPLEX)
892: cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
893: #else
894: cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
895: #endif
896: }
897: VecRestoreArray(v,&array);
898: }
899: MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
900: MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
901: VecDestroy(&v);
902: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);
903: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
904: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
905: MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);
906: ISDestroy(&is_p0_p);
907: ISLocalToGlobalMappingDestroy(&N_to_reor);
908: }
909: MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);
910: MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);
911: MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);
912: #if defined(PETSC_HAVE_MUMPS)
913: PetscStrcpy(solver,MATSOLVERMUMPS);
914: #else
915: PetscStrcpy(solver,MATSOLVERMKL_PARDISO);
916: #endif
917: PetscOptionsGetString(NULL,((PetscObject)A)->prefix,"-mat_solver_package",solver,256,NULL);
919: /* when using the benign subspace trick, the local Schur complements are SPD */
920: if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
922: if (n_I) { /* TODO add ordering from options */
923: IS is_schur;
925: if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
926: MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);
927: } else {
928: MatGetFactor(A,solver,MAT_FACTOR_LU,&F);
929: }
930: MatSetErrorIfFailure(A,PETSC_TRUE);
932: /* subsets ordered last */
933: ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);
934: MatFactorSetSchurIS(F,is_schur);
935: ISDestroy(&is_schur);
937: /* factorization step */
938: if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
939: MatCholeskyFactorSymbolic(F,A,NULL,NULL);
940: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
941: MatMumpsSetIcntl(F,19,2);
942: #endif
943: MatCholeskyFactorNumeric(F,A,NULL);
944: S_lower_triangular = PETSC_TRUE;
945: } else {
946: MatLUFactorSymbolic(F,A,NULL,NULL,NULL);
947: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
948: MatMumpsSetIcntl(F,19,3);
949: #endif
950: MatLUFactorNumeric(F,A,NULL);
951: }
952: MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");
954: /* get explicit Schur Complement computed during numeric factorization */
955: MatFactorGetSchurComplement(F,&S_all,NULL);
956: MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);
957: MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);
959: /* we can reuse the solvers if we are not using the economic version */
960: reuse_solvers = (PetscBool)(reuse_solvers && !economic);
961: factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
962: solver_S = PETSC_TRUE;
964: /* update the Schur complement with the change of basis on the pressures */
965: if (benign_n) {
966: PetscScalar *S_data,*cs_AIB,*AIIm1_data;
967: Vec v;
968: PetscInt sizeA;
970: MatDenseGetArray(S_all,&S_data);
971: MatCreateVecs(A,&v,NULL);
972: VecGetSize(v,&sizeA);
973: #if defined(PETSC_HAVE_MUMPS)
974: MatMumpsSetIcntl(F,26,0);
975: #endif
976: #if defined(PETSC_HAVE_MKL_PARDISO)
977: MatMkl_PardisoSetCntl(F,70,1);
978: #endif
979: MatDenseGetArray(cs_AIB_mat,&cs_AIB);
980: MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
981: for (i=0;i<benign_n;i++) {
982: Vec benign_AIIm1_ones;
983: PetscScalar *array,sum = 0.,one = 1.;
984: const PetscInt *idxs;
985: PetscInt j,nz;
986: PetscBLASInt B_k,B_n;
988: VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);
989: VecCopy(benign_AIIm1_ones,v);
990: MatSolve(F,v,benign_AIIm1_ones);
991: ISGetLocalSize(is_p_r[i],&nz);
992: ISGetIndices(is_p_r[i],&idxs);
993: /* p0 dof (eliminated) is excluded from the sum */
994: for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
995: ISRestoreIndices(is_p_r[i],&idxs);
996: MatMult(A,benign_AIIm1_ones,v);
997: VecGetArray(v,&array);
998: /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
999: /* cs_AIB already scaled by 1./nz */
1000: B_k = 1;
1001: PetscBLASIntCast(size_schur,&B_n);
1002: PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1003: sum = 1.;
1004: PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1005: VecRestoreArray(v,&array);
1006: /* set p0 entry of AIIm1_ones to zero */
1007: ISGetIndices(is_p_r[i],&idxs);
1008: for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1009: ISRestoreIndices(is_p_r[i],&idxs);
1010: VecDestroy(&benign_AIIm1_ones);
1011: }
1012: /* restore defaults */
1013: #if defined(PETSC_HAVE_MUMPS)
1014: MatMumpsSetIcntl(F,26,-1);
1015: #endif
1016: #if defined(PETSC_HAVE_MKL_PARDISO)
1017: MatMkl_PardisoSetCntl(F,70,0);
1018: #endif
1019: MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
1020: MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
1021: VecDestroy(&v);
1022: MatDenseRestoreArray(S_all,&S_data);
1023: }
1024: if (!reuse_solvers) {
1025: for (i=0;i<benign_n;i++) {
1026: ISDestroy(&is_p_r[i]);
1027: }
1028: PetscFree(is_p_r);
1029: MatDestroy(&cs_AIB_mat);
1030: MatDestroy(&benign_AIIm1_ones_mat);
1031: }
1032: } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1033: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);
1034: reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1035: factor_workaround = PETSC_FALSE;
1036: solver_S = PETSC_FALSE;
1037: }
1039: if (reuse_solvers) {
1040: Mat A_II,Afake;
1041: Vec vec1_B;
1042: PCBDDCReuseSolvers msolv_ctx;
1043: PetscInt n_R;
1045: if (sub_schurs->reuse_solver) {
1046: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1047: } else {
1048: PetscNew(&sub_schurs->reuse_solver);
1049: }
1050: msolv_ctx = sub_schurs->reuse_solver;
1051: MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);
1052: PetscObjectReference((PetscObject)F);
1053: msolv_ctx->F = F;
1054: MatCreateVecs(F,&msolv_ctx->sol,NULL);
1055: /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1056: {
1057: PetscScalar *array;
1058: PetscInt n;
1060: VecGetLocalSize(msolv_ctx->sol,&n);
1061: VecGetArray(msolv_ctx->sol,&array);
1062: VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);
1063: VecRestoreArray(msolv_ctx->sol,&array);
1064: }
1065: msolv_ctx->has_vertices = schur_has_vertices;
1067: /* interior solver */
1068: PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);
1069: PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);
1070: PCSetType(msolv_ctx->interior_solver,PCSHELL);
1071: PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);
1072: PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);
1073: PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);
1075: /* correction solver */
1076: PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);
1077: PCSetType(msolv_ctx->correction_solver,PCSHELL);
1078: PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);
1079: PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);
1080: PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);
1082: /* scatter and vecs for Schur complement solver */
1083: MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);
1084: MatCreateVecs(sub_schurs->S,&vec1_B,NULL);
1085: if (!schur_has_vertices) {
1086: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);
1087: VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);
1088: PetscObjectReference((PetscObject)is_A_all);
1089: msolv_ctx->is_R = is_A_all;
1090: } else {
1091: IS is_B_all;
1092: const PetscInt* idxs;
1093: PetscInt dual,n_v,n;
1095: ISGetLocalSize(sub_schurs->is_vertices,&n_v);
1096: dual = size_schur - n_v;
1097: ISGetLocalSize(is_A_all,&n);
1098: ISGetIndices(is_A_all,&idxs);
1099: ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);
1100: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);
1101: ISDestroy(&is_B_all);
1102: ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);
1103: VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);
1104: ISDestroy(&is_B_all);
1105: ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);
1106: ISRestoreIndices(is_A_all,&idxs);
1107: }
1108: ISGetLocalSize(msolv_ctx->is_R,&n_R);
1109: MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);
1110: MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);
1111: MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);
1112: PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);
1113: MatDestroy(&Afake);
1114: VecDestroy(&vec1_B);
1116: /* communicate benign info to solver context */
1117: if (benign_n) {
1118: PetscScalar *array;
1120: msolv_ctx->benign_n = benign_n;
1121: msolv_ctx->benign_zerodiag_subs = is_p_r;
1122: PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);
1123: msolv_ctx->benign_csAIB = cs_AIB_mat;
1124: MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);
1125: VecGetArray(msolv_ctx->benign_corr_work,&array);
1126: VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);
1127: VecRestoreArray(msolv_ctx->benign_corr_work,&array);
1128: msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1129: }
1130: } else {
1131: if (sub_schurs->reuse_solver) {
1132: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1133: }
1134: PetscFree(sub_schurs->reuse_solver);
1135: }
1136: MatDestroy(&A);
1137: ISDestroy(&is_A_all);
1139: /* Work arrays */
1140: PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
1142: /* matrices for deluxe scaling and adaptive selection */
1143: if (compute_Stilda) {
1144: if (!sub_schurs->sum_S_Ej_tilda_all) {
1145: MatDuplicate(sub_schurs->sum_S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);
1146: }
1147: if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1148: MatDuplicate(sub_schurs->sum_S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);
1149: }
1150: }
1152: /* S_Ej_all */
1153: cum = cum2 = 0;
1154: MatDenseGetArray(S_all,&S_data);
1155: for (i=0;i<sub_schurs->n_subs;i++) {
1156: PetscInt j;
1158: /* get S_E */
1159: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1160: if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1161: PetscInt k;
1162: for (k=0;k<subset_size;k++) {
1163: for (j=k;j<subset_size;j++) {
1164: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1165: work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1166: }
1167: }
1168: } else { /* just copy to workspace */
1169: PetscInt k;
1170: for (k=0;k<subset_size;k++) {
1171: for (j=0;j<subset_size;j++) {
1172: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1173: }
1174: }
1175: }
1176: /* insert S_E values */
1177: for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1178: if (sub_schurs->change) {
1179: Mat change_sub,SEj,T;
1181: /* change basis */
1182: KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1183: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1184: if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1185: Mat T2;
1186: MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1187: MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1188: MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1189: MatDestroy(&T2);
1190: } else {
1191: MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1192: }
1193: MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1194: MatDestroy(&T);
1195: MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);
1196: MatDestroy(&SEj);
1197: }
1198: if (deluxe) {
1199: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1200: /* if adaptivity is requested, invert S_E blocks */
1201: if (compute_Stilda) {
1202: PetscBLASIntCast(subset_size,&B_N);
1203: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1204: if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1205: PetscInt k;
1206: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1207: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1208: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1209: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1210: for (k=0;k<subset_size;k++) {
1211: for (j=k;j<subset_size;j++) {
1212: work[j*subset_size+k] = work[k*subset_size+j];
1213: }
1214: }
1215: } else {
1216: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1217: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1218: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1219: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1220: }
1221: PetscFPTrapPop();
1222: MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1223: }
1224: } else if (compute_Stilda) { /* not using deluxe */
1225: Mat SEj;
1226: Vec D;
1227: PetscScalar *array;
1229: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1230: VecGetArray(Dall,&array);
1231: VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);
1232: VecRestoreArray(Dall,&array);
1233: VecShift(D,-1.);
1234: MatDiagonalScale(SEj,D,D);
1235: MatDestroy(&SEj);
1236: VecDestroy(&D);
1237: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1238: }
1239: cum += subset_size;
1240: cum2 += subset_size*(size_schur + 1);
1241: }
1242: MatDenseRestoreArray(S_all,&S_data);
1244: if (solver_S) {
1245: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1246: }
1248: schur_factor = NULL;
1249: if (compute_Stilda && size_active_schur) {
1251: if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1252: PetscInt j;
1253: for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1254: MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);
1255: } else {
1256: Mat S_all_inv=NULL;
1257: if (solver_S) {
1258: /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1259: The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */
1260: if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1261: PetscScalar *data;
1262: PetscInt nd = 0;
1264: if (!sub_schurs->is_posdef) {
1265: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1266: }
1267: MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1268: MatDenseGetArray(S_all_inv,&data);
1269: if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1270: ISGetLocalSize(sub_schurs->is_dir,&nd);
1271: }
1273: /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1274: if (schur_has_vertices) {
1275: Mat M;
1276: PetscScalar *tdata;
1277: PetscInt nv = 0, news;
1279: ISGetLocalSize(sub_schurs->is_vertices,&nv);
1280: news = size_active_schur + nv;
1281: PetscCalloc1(news*news,&tdata);
1282: for (i=0;i<size_active_schur;i++) {
1283: PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));
1284: PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));
1285: }
1286: for (i=0;i<nv;i++) {
1287: PetscInt k = i+size_active_schur;
1288: PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));
1289: }
1291: MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);
1292: MatSetOption(M,MAT_SPD,PETSC_TRUE);
1293: MatCholeskyFactor(M,NULL,NULL);
1294: /* save the factors */
1295: cum = 0;
1296: PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);
1297: for (i=0;i<size_active_schur;i++) {
1298: PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));
1299: cum += size_active_schur - i;
1300: }
1301: for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1302: MatSeqDenseInvertFactors_Private(M);
1303: /* move back just the active dofs to the Schur complement */
1304: for (i=0;i<size_active_schur;i++) {
1305: PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));
1306: }
1307: PetscFree(tdata);
1308: MatDestroy(&M);
1309: } else { /* we can factorize and invert just the activedofs */
1310: Mat M;
1312: MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);
1313: MatSeqDenseSetLDA(M,size_schur);
1314: MatSetOption(M,MAT_SPD,PETSC_TRUE);
1315: MatCholeskyFactor(M,NULL,NULL);
1316: MatSeqDenseInvertFactors_Private(M);
1317: MatDestroy(&M);
1318: for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1319: }
1320: MatDenseRestoreArray(S_all_inv,&data);
1321: } else { /* use MatFactor calls to invert S */
1322: MatFactorInvertSchurComplement(F);
1323: MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1324: }
1325: } else { /* we need to invert explicitly since we are not using MatFactor for S */
1326: PetscObjectReference((PetscObject)S_all);
1327: S_all_inv = S_all;
1328: MatDenseGetArray(S_all_inv,&S_data);
1329: PetscBLASIntCast(size_schur,&B_N);
1330: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1331: if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1332: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1333: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1334: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1335: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1336: } else {
1337: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1338: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1339: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1340: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1341: }
1342: PetscFPTrapPop();
1343: MatDenseRestoreArray(S_all_inv,&S_data);
1344: }
1345: /* S_Ej_tilda_all */
1346: cum = cum2 = 0;
1347: MatDenseGetArray(S_all_inv,&S_data);
1348: for (i=0;i<sub_schurs->n_subs;i++) {
1349: PetscInt j;
1351: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1352: /* get (St^-1)_E */
1353: /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1354: will be properly accessed later during adaptive selection */
1355: if (S_lower_triangular) {
1356: PetscInt k;
1357: if (sub_schurs->change) {
1358: for (k=0;k<subset_size;k++) {
1359: for (j=k;j<subset_size;j++) {
1360: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1361: work[j*subset_size+k] = work[k*subset_size+j];
1362: }
1363: }
1364: } else {
1365: for (k=0;k<subset_size;k++) {
1366: for (j=k;j<subset_size;j++) {
1367: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1368: }
1369: }
1370: }
1371: } else {
1372: PetscInt k;
1373: for (k=0;k<subset_size;k++) {
1374: for (j=0;j<subset_size;j++) {
1375: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1376: }
1377: }
1378: }
1379: if (sub_schurs->change) {
1380: Mat change_sub,SEj,T;
1382: /* change basis */
1383: KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1384: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1385: if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1386: Mat T2;
1387: MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1388: MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1389: MatDestroy(&T2);
1390: MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1391: } else {
1392: MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1393: }
1394: MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1395: MatDestroy(&T);
1396: /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1397: MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);
1398: MatDestroy(&SEj);
1399: }
1400: for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1401: MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1402: cum += subset_size;
1403: cum2 += subset_size*(size_schur + 1);
1404: }
1405: MatDenseRestoreArray(S_all_inv,&S_data);
1406: if (solver_S) {
1407: if (schur_has_vertices) {
1408: MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);
1409: } else {
1410: MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);
1411: }
1412: }
1413: MatDestroy(&S_all_inv);
1414: }
1416: /* move back factors if needed */
1417: if (schur_has_vertices) {
1418: Mat S_tmp;
1419: PetscInt nd = 0;
1421: if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1422: MatFactorGetSchurComplement(F,&S_tmp,NULL);
1423: if (sub_schurs->is_posdef) {
1424: PetscScalar *data;
1426: MatDenseGetArray(S_tmp,&data);
1427: PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));
1429: if (S_lower_triangular) {
1430: cum = 0;
1431: for (i=0;i<size_active_schur;i++) {
1432: PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));
1433: cum += size_active_schur-i;
1434: }
1435: } else {
1436: PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));
1437: }
1438: if (sub_schurs->is_dir) {
1439: ISGetLocalSize(sub_schurs->is_dir,&nd);
1440: for (i=0;i<nd;i++) {
1441: data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1442: }
1443: }
1444: /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1445: set the diagonal entry of the Schur factor to a very large value */
1446: for (i=size_active_schur+nd;i<size_schur;i++) {
1447: data[i*(size_schur+1)] = infty;
1448: }
1449: MatDenseRestoreArray(S_tmp,&data);
1450: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1451: MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);
1452: }
1453: } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1454: PetscScalar *data;
1455: PetscInt nd = 0;
1457: if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1458: ISGetLocalSize(sub_schurs->is_dir,&nd);
1459: }
1460: MatFactorGetSchurComplement(F,&S_all,NULL);
1461: MatDenseGetArray(S_all,&data);
1462: for (i=0;i<size_active_schur;i++) {
1463: PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));
1464: }
1465: for (i=size_active_schur+nd;i<size_schur;i++) {
1466: PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));
1467: data[i*(size_schur+1)] = infty;
1468: }
1469: MatDenseRestoreArray(S_all,&data);
1470: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1471: }
1472: PetscFree2(dummy_idx,work);
1473: PetscFree(schur_factor);
1474: VecDestroy(&Dall);
1475: }
1476: PetscFree(nnz);
1477: MatDestroy(&F);
1478: ISDestroy(&is_I_layer);
1479: MatDestroy(&S_all);
1480: MatDestroy(&A_BB);
1481: MatDestroy(&A_IB);
1482: MatDestroy(&A_BI);
1483: MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1484: MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1485: if (compute_Stilda) {
1486: MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1487: MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1488: if (deluxe) {
1489: MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1490: MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1491: }
1492: }
1494: /* Global matrix of all assembled Schur on subsets */
1495: MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);
1496: MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);
1497: MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);
1499: /* Get local part of (\sum_j S_Ej) */
1500: MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);
1501: if (!sub_schurs->sum_S_Ej_all) {
1502: MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);
1503: } else {
1504: MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);
1505: }
1507: /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1508: if (compute_Stilda) {
1509: MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);
1510: MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);
1511: MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);
1512: MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);
1513: if (deluxe) {
1514: MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);
1515: MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);
1516: MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);
1517: MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);
1518: } else {
1519: PetscScalar *array;
1520: PetscInt cum;
1522: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1523: cum = 0;
1524: for (i=0;i<sub_schurs->n_subs;i++) {
1525: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1526: PetscBLASIntCast(subset_size,&B_N);
1527: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1528: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1529: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1530: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1531: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1532: PetscFPTrapPop();
1533: cum += subset_size*subset_size;
1534: }
1535: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1536: PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);
1537: MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1538: sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1539: }
1540: }
1541: MatDestroySubMatrices(1,&submats);
1543: /* free workspace */
1544: PetscFree(submats);
1545: PetscFree2(Bwork,pivots);
1546: MatDestroy(&global_schur_subsets);
1547: MatDestroy(&work_mat);
1548: ISDestroy(&all_subsets_n);
1549: PetscCommDestroy(&comm_n);
1550: return(0);
1551: }
1553: PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1554: {
1555: IS *faces,*edges,*all_cc,vertices;
1556: PetscInt i,n_faces,n_edges,n_all_cc;
1557: PetscBool is_sorted;
1558: PetscErrorCode ierr;
1561: ISSorted(is_I,&is_sorted);
1562: if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1563: ISSorted(is_B,&is_sorted);
1564: if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1566: /* reset any previous data */
1567: PCBDDCSubSchursReset(sub_schurs);
1569: /* get index sets for faces and edges (already sorted by global ordering) */
1570: PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1571: n_all_cc = n_faces+n_edges;
1572: PetscBTCreate(n_all_cc,&sub_schurs->is_edge);
1573: PetscMalloc1(n_all_cc,&all_cc);
1574: for (i=0;i<n_faces;i++) {
1575: if (copycc) {
1576: ISDuplicate(faces[i],&all_cc[i]);
1577: } else {
1578: PetscObjectReference((PetscObject)faces[i]);
1579: all_cc[i] = faces[i];
1580: }
1581: }
1582: for (i=0;i<n_edges;i++) {
1583: if (copycc) {
1584: ISDuplicate(edges[i],&all_cc[n_faces+i]);
1585: } else {
1586: PetscObjectReference((PetscObject)edges[i]);
1587: all_cc[n_faces+i] = edges[i];
1588: }
1589: PetscBTSet(sub_schurs->is_edge,n_faces+i);
1590: }
1591: PetscObjectReference((PetscObject)vertices);
1592: sub_schurs->is_vertices = vertices;
1593: PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1594: sub_schurs->is_dir = NULL;
1595: PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);
1597: /* Determine if MatFactor can be used */
1598: sub_schurs->schur_explicit = PETSC_FALSE;
1599: #if defined(PETSC_HAVE_MUMPS)
1600: sub_schurs->schur_explicit = PETSC_TRUE;
1601: #endif
1602: #if defined(PETSC_HAVE_MKL_PARDISO)
1603: sub_schurs->schur_explicit = PETSC_TRUE;
1604: #endif
1606: PetscObjectReference((PetscObject)is_I);
1607: sub_schurs->is_I = is_I;
1608: PetscObjectReference((PetscObject)is_B);
1609: sub_schurs->is_B = is_B;
1610: PetscObjectReference((PetscObject)graph->l2gmap);
1611: sub_schurs->l2gmap = graph->l2gmap;
1612: PetscObjectReference((PetscObject)BtoNmap);
1613: sub_schurs->BtoNmap = BtoNmap;
1614: sub_schurs->n_subs = n_all_cc;
1615: sub_schurs->is_subs = all_cc;
1616: if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1617: for (i=0;i<sub_schurs->n_subs;i++) {
1618: ISSort(sub_schurs->is_subs[i]);
1619: }
1620: }
1621: sub_schurs->S_Ej_all = NULL;
1622: sub_schurs->sum_S_Ej_all = NULL;
1623: sub_schurs->sum_S_Ej_inv_all = NULL;
1624: sub_schurs->sum_S_Ej_tilda_all = NULL;
1625: sub_schurs->is_Ej_all = NULL;
1626: return(0);
1627: }
1629: PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1630: {
1631: PCBDDCSubSchurs schurs_ctx;
1632: PetscErrorCode ierr;
1635: PetscNew(&schurs_ctx);
1636: schurs_ctx->n_subs = 0;
1637: *sub_schurs = schurs_ctx;
1638: return(0);
1639: }
1641: PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1642: {
1643: PetscInt i;
1647: if (!sub_schurs) return(0);
1648: MatDestroy(&sub_schurs->A);
1649: MatDestroy(&sub_schurs->S);
1650: ISDestroy(&sub_schurs->is_I);
1651: ISDestroy(&sub_schurs->is_B);
1652: ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);
1653: ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);
1654: MatDestroy(&sub_schurs->S_Ej_all);
1655: MatDestroy(&sub_schurs->sum_S_Ej_all);
1656: MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1657: MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);
1658: ISDestroy(&sub_schurs->is_Ej_all);
1659: ISDestroy(&sub_schurs->is_vertices);
1660: ISDestroy(&sub_schurs->is_dir);
1661: PetscBTDestroy(&sub_schurs->is_edge);
1662: for (i=0;i<sub_schurs->n_subs;i++) {
1663: ISDestroy(&sub_schurs->is_subs[i]);
1664: }
1665: if (sub_schurs->n_subs) {
1666: PetscFree(sub_schurs->is_subs);
1667: }
1668: if (sub_schurs->reuse_solver) {
1669: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1670: }
1671: PetscFree(sub_schurs->reuse_solver);
1672: if (sub_schurs->change) {
1673: for (i=0;i<sub_schurs->n_subs;i++) {
1674: KSPDestroy(&sub_schurs->change[i]);
1675: ISDestroy(&sub_schurs->change_primal_sub[i]);
1676: }
1677: }
1678: PetscFree(sub_schurs->change);
1679: PetscFree(sub_schurs->change_primal_sub);
1680: sub_schurs->n_subs = 0;
1681: return(0);
1682: }
1684: PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1685: {
1689: PCBDDCSubSchursReset(*sub_schurs);
1690: PetscFree(*sub_schurs);
1691: return(0);
1692: }
1694: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1695: {
1696: PetscInt i,j,n;
1700: n = 0;
1701: for (i=-n_prev;i<0;i++) {
1702: PetscInt start_dof = queue_tip[i];
1703: for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1704: PetscInt dof = adjncy[j];
1705: if (!PetscBTLookup(touched,dof)) {
1706: PetscBTSet(touched,dof);
1707: queue_tip[n] = dof;
1708: n++;
1709: }
1710: }
1711: }
1712: *n_added = n;
1713: return(0);
1714: }