Actual source code: bddcschurs.c
petsc-3.10.5 2019-03-28
1: #include <../src/ksp/pc/impls/bddc/bddc.h>
2: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3: #include <../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 PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer)
232: {
233: PCBDDCReuseSolvers ctx;
234: PetscBool iascii;
235: PetscErrorCode ierr;
238: PCShellGetContext(pc,(void **)&ctx);
239: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
240: if (iascii) {
241: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
242: }
243: MatView(ctx->F,viewer);
244: if (iascii) {
245: PetscViewerPopFormat(viewer);
246: }
247: return(0);
248: }
250: static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
251: {
252: PetscInt i;
256: MatDestroy(&reuse->F);
257: VecDestroy(&reuse->sol);
258: VecDestroy(&reuse->rhs);
259: PCDestroy(&reuse->interior_solver);
260: PCDestroy(&reuse->correction_solver);
261: ISDestroy(&reuse->is_R);
262: ISDestroy(&reuse->is_B);
263: VecScatterDestroy(&reuse->correction_scatter_B);
264: VecDestroy(&reuse->sol_B);
265: VecDestroy(&reuse->rhs_B);
266: for (i=0;i<reuse->benign_n;i++) {
267: ISDestroy(&reuse->benign_zerodiag_subs[i]);
268: }
269: PetscFree(reuse->benign_zerodiag_subs);
270: PetscFree(reuse->benign_save_vals);
271: MatDestroy(&reuse->benign_csAIB);
272: MatDestroy(&reuse->benign_AIIm1ones);
273: VecDestroy(&reuse->benign_corr_work);
274: VecDestroy(&reuse->benign_dummy_schur_vec);
275: return(0);
276: }
278: static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
279: {
280: Mat B, C, D, Bd, Cd, AinvBd;
281: KSP ksp;
282: PC pc;
283: PetscBool isLU, isILU, isCHOL, Bdense, Cdense;
284: PetscReal fill = 2.0;
285: PetscInt n_I;
286: PetscMPIInt size;
290: MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);
291: if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
292: if (reuse == MAT_REUSE_MATRIX) {
293: PetscBool Sdense;
295: PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);
296: if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
297: }
298: MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);
299: MatSchurComplementGetKSP(M, &ksp);
300: KSPGetPC(ksp, &pc);
301: PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);
302: PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);
303: PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);
304: PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);
305: PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);
306: MatGetSize(B,&n_I,NULL);
307: if (n_I) {
308: if (!Bdense) {
309: MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);
310: } else {
311: Bd = B;
312: }
314: if (isLU || isILU || isCHOL) {
315: Mat fact;
316: KSPSetUp(ksp);
317: PCFactorGetMatrix(pc, &fact);
318: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
319: MatMatSolve(fact, Bd, AinvBd);
320: } else {
321: PetscBool ex = PETSC_TRUE;
323: if (ex) {
324: Mat Ainvd;
326: PCComputeExplicitOperator(pc, &Ainvd);
327: MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);
328: MatDestroy(&Ainvd);
329: } else {
330: Vec sol,rhs;
331: PetscScalar *arrayrhs,*arraysol;
332: PetscInt i,nrhs,n;
334: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
335: MatGetSize(Bd,&n,&nrhs);
336: MatDenseGetArray(Bd,&arrayrhs);
337: MatDenseGetArray(AinvBd,&arraysol);
338: KSPGetSolution(ksp,&sol);
339: KSPGetRhs(ksp,&rhs);
340: for (i=0;i<nrhs;i++) {
341: VecPlaceArray(rhs,arrayrhs+i*n);
342: VecPlaceArray(sol,arraysol+i*n);
343: KSPSolve(ksp,rhs,sol);
344: VecResetArray(rhs);
345: VecResetArray(sol);
346: }
347: MatDenseRestoreArray(Bd,&arrayrhs);
348: MatDenseRestoreArray(AinvBd,&arrayrhs);
349: }
350: }
351: if (!Bdense & !issym) {
352: MatDestroy(&Bd);
353: }
355: if (!issym) {
356: if (!Cdense) {
357: MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);
358: } else {
359: Cd = C;
360: }
361: MatMatMult(Cd, AinvBd, reuse, fill, S);
362: if (!Cdense) {
363: MatDestroy(&Cd);
364: }
365: } else {
366: MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);
367: if (!Bdense) {
368: MatDestroy(&Bd);
369: }
370: }
371: MatDestroy(&AinvBd);
372: }
374: if (D) {
375: Mat Dd;
376: PetscBool Ddense;
378: PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);
379: if (!Ddense) {
380: MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);
381: } else {
382: Dd = D;
383: }
384: if (n_I) {
385: MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);
386: } else {
387: if (reuse == MAT_INITIAL_MATRIX) {
388: MatDuplicate(Dd,MAT_COPY_VALUES,S);
389: } else {
390: MatCopy(Dd,*S,SAME_NONZERO_PATTERN);
391: }
392: }
393: if (!Ddense) {
394: MatDestroy(&Dd);
395: }
396: } else {
397: MatScale(*S,-1.0);
398: }
399: return(0);
400: }
402: 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)
403: {
404: Mat F,A_II,A_IB,A_BI,A_BB,AE_II;
405: Mat S_all;
406: Mat global_schur_subsets,work_mat,*submats;
407: ISLocalToGlobalMapping l2gmap_subsets;
408: IS is_I,is_I_layer;
409: IS all_subsets,all_subsets_mult,all_subsets_n;
410: PetscInt *nnz,*all_local_idx_N;
411: PetscInt *auxnum1,*auxnum2;
412: PetscInt i,subset_size,max_subset_size;
413: PetscInt n_B,extra,local_size,global_size;
414: PetscBLASInt B_N,B_ierr,B_lwork,*pivots;
415: PetscScalar *Bwork;
416: MPI_Comm comm_n;
417: PetscBool deluxe = PETSC_TRUE;
418: PetscBool use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
419: PetscBool print_schurs = PETSC_FALSE;
420: PetscErrorCode ierr;
423: MatDestroy(&sub_schurs->A);
424: MatDestroy(&sub_schurs->S);
425: /* convert matrix if needed */
426: if (Ain) {
427: PetscBool isseqaij;
428: PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);
429: if (isseqaij) {
430: PetscObjectReference((PetscObject)Ain);
431: sub_schurs->A = Ain;
432: } else {
433: MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);
434: }
435: }
437: PetscObjectReference((PetscObject)Sin);
438: sub_schurs->S = Sin;
439: if (sub_schurs->schur_explicit) {
440: sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
441: }
443: /* preliminary checks */
444: 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");
446: if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
448: /* debug prints */
449: if (sub_schurs->debug) {
450: PetscMPIInt size,rank;
451: PetscInt nr,*print_schurs_ranks;
452: PetscBool flg;
454: MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);
455: MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
456: nr = size;
457: PetscMalloc1(nr,&print_schurs_ranks);
458: PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
459: PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);
460: if (!flg) print_schurs = PETSC_TRUE;
461: else {
462: for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; }
463: }
464: PetscOptionsEnd();
465: PetscFree(print_schurs_ranks);
466: }
469: /* restrict work on active processes */
470: if (sub_schurs->restrict_comm) {
471: PetscSubcomm subcomm;
472: PetscMPIInt color,rank;
474: color = 0;
475: if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
476: MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
477: PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);
478: PetscSubcommSetNumber(subcomm,2);
479: PetscSubcommSetTypeGeneral(subcomm,color,rank);
480: PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);
481: PetscSubcommDestroy(&subcomm);
482: if (!sub_schurs->n_subs) {
483: PetscCommDestroy(&comm_n);
484: return(0);
485: }
486: } else {
487: PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);
488: }
490: /* get Schur complement matrices */
491: if (!sub_schurs->schur_explicit) {
492: Mat tA_IB,tA_BI,tA_BB;
493: PetscBool isseqsbaij;
494: MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);
495: PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);
496: if (isseqsbaij) {
497: MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);
498: MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);
499: MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);
500: } else {
501: PetscObjectReference((PetscObject)tA_BB);
502: A_BB = tA_BB;
503: PetscObjectReference((PetscObject)tA_IB);
504: A_IB = tA_IB;
505: PetscObjectReference((PetscObject)tA_BI);
506: A_BI = tA_BI;
507: }
508: } else {
509: A_II = NULL;
510: A_IB = NULL;
511: A_BI = NULL;
512: A_BB = NULL;
513: }
514: S_all = NULL;
516: /* determine interior problems */
517: ISGetLocalSize(sub_schurs->is_I,&i);
518: if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
519: PetscBT touched;
520: const PetscInt* idx_B;
521: PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
523: if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
524: /* get sizes */
525: ISGetLocalSize(sub_schurs->is_I,&n_I);
526: ISGetLocalSize(sub_schurs->is_B,&n_B);
528: PetscMalloc1(n_I+n_B,&local_numbering);
529: PetscBTCreate(n_I+n_B,&touched);
530: PetscBTMemzero(n_I+n_B,touched);
532: /* all boundary dofs must be skipped when adding layers */
533: ISGetIndices(sub_schurs->is_B,&idx_B);
534: for (j=0;j<n_B;j++) {
535: PetscBTSet(touched,idx_B[j]);
536: }
537: PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));
538: ISRestoreIndices(sub_schurs->is_B,&idx_B);
540: /* add prescribed number of layers of dofs */
541: n_local_dofs = n_B;
542: n_prev_added = n_B;
543: for (layer=0;layer<nlayers;layer++) {
544: PetscInt n_added;
545: if (n_local_dofs == n_I+n_B) break;
546: 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);
547: PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);
548: n_prev_added = n_added;
549: n_local_dofs += n_added;
550: if (!n_added) break;
551: }
552: PetscBTDestroy(&touched);
554: /* IS for I layer dofs in original numbering */
555: ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);
556: PetscFree(local_numbering);
557: ISSort(is_I_layer);
558: /* IS for I layer dofs in I numbering */
559: if (!sub_schurs->schur_explicit) {
560: ISLocalToGlobalMapping ItoNmap;
561: ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);
562: ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);
563: ISLocalToGlobalMappingDestroy(&ItoNmap);
565: /* II block */
566: MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);
567: }
568: } else {
569: PetscInt n_I;
571: /* IS for I dofs in original numbering */
572: PetscObjectReference((PetscObject)sub_schurs->is_I);
573: is_I_layer = sub_schurs->is_I;
575: /* IS for I dofs in I numbering (strided 1) */
576: if (!sub_schurs->schur_explicit) {
577: ISGetSize(sub_schurs->is_I,&n_I);
578: ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);
580: /* II block is the same */
581: PetscObjectReference((PetscObject)A_II);
582: AE_II = A_II;
583: }
584: }
586: /* Get info on subset sizes and sum of all subsets sizes */
587: max_subset_size = 0;
588: local_size = 0;
589: for (i=0;i<sub_schurs->n_subs;i++) {
590: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
591: max_subset_size = PetscMax(subset_size,max_subset_size);
592: local_size += subset_size;
593: }
595: /* Work arrays for local indices */
596: extra = 0;
597: ISGetLocalSize(sub_schurs->is_B,&n_B);
598: if (sub_schurs->schur_explicit && is_I_layer) {
599: ISGetLocalSize(is_I_layer,&extra);
600: }
601: PetscMalloc1(n_B+extra,&all_local_idx_N);
602: if (extra) {
603: const PetscInt *idxs;
604: ISGetIndices(is_I_layer,&idxs);
605: PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));
606: ISRestoreIndices(is_I_layer,&idxs);
607: }
608: PetscMalloc1(local_size,&nnz);
609: PetscMalloc1(sub_schurs->n_subs,&auxnum1);
610: PetscMalloc1(sub_schurs->n_subs,&auxnum2);
612: /* Get local indices in local numbering */
613: local_size = 0;
614: for (i=0;i<sub_schurs->n_subs;i++) {
615: PetscInt j;
616: const PetscInt *idxs;
618: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
619: ISGetIndices(sub_schurs->is_subs[i],&idxs);
620: /* start (smallest in global ordering) and multiplicity */
621: auxnum1[i] = idxs[0];
622: auxnum2[i] = subset_size;
623: /* subset indices in local numbering */
624: PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));
625: ISRestoreIndices(sub_schurs->is_subs[i],&idxs);
626: for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
627: local_size += subset_size;
628: }
630: /* allocate extra workspace needed only for GETRI or SYTRF */
631: use_potr = use_sytr = PETSC_FALSE;
632: if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
633: use_potr = PETSC_TRUE;
634: } else if (sub_schurs->is_symmetric) {
635: use_sytr = PETSC_TRUE;
636: }
637: if (local_size && !use_potr) {
638: PetscScalar lwork,dummyscalar = 0.;
639: PetscBLASInt dummyint = 0;
641: B_lwork = -1;
642: PetscBLASIntCast(local_size,&B_N);
643: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
644: if (use_sytr) {
645: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
646: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
647: } else {
648: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
649: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
650: }
651: PetscFPTrapPop();
652: PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);
653: PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);
654: } else {
655: Bwork = NULL;
656: pivots = NULL;
657: }
659: /* prepare parallel matrices for summing up properly schurs on subsets */
660: ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);
661: ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);
662: ISDestroy(&all_subsets_n);
663: ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);
664: ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);
665: ISDestroy(&all_subsets);
666: ISDestroy(&all_subsets_mult);
667: ISGetLocalSize(all_subsets_n,&i);
668: if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
669: ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);
670: MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);
671: ISLocalToGlobalMappingDestroy(&l2gmap_subsets);
672: MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);
673: MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);
674: MatSetType(global_schur_subsets,MATMPIAIJ);
676: /* subset indices in local boundary numbering */
677: if (!sub_schurs->is_Ej_all) {
678: PetscInt *all_local_idx_B;
680: PetscMalloc1(local_size,&all_local_idx_B);
681: ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);
682: if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size);
683: ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);
684: }
686: if (change) {
687: ISLocalToGlobalMapping BtoS;
688: IS change_primal_B;
689: IS change_primal_all;
691: if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
692: if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
693: PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);
694: for (i=0;i<sub_schurs->n_subs;i++) {
695: ISLocalToGlobalMapping NtoS;
696: ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);
697: ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);
698: ISLocalToGlobalMappingDestroy(&NtoS);
699: }
700: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);
701: ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);
702: ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);
703: ISLocalToGlobalMappingDestroy(&BtoS);
704: ISDestroy(&change_primal_B);
705: PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);
706: for (i=0;i<sub_schurs->n_subs;i++) {
707: Mat change_sub;
709: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
710: KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);
711: KSPSetType(sub_schurs->change[i],KSPPREONLY);
712: if (!sub_schurs->change_with_qr) {
713: MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);
714: } else {
715: Mat change_subt;
716: MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);
717: MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);
718: MatDestroy(&change_subt);
719: }
720: KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);
721: MatDestroy(&change_sub);
722: KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);
723: KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");
724: }
725: ISDestroy(&change_primal_all);
726: }
728: /* Local matrix of all local Schur on subsets (transposed) */
729: if (!sub_schurs->S_Ej_all) {
730: MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);
731: MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);
732: MatSetType(sub_schurs->S_Ej_all,MATAIJ);
733: MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);
734: }
736: /* Compute Schur complements explicitly */
737: F = NULL;
738: if (!sub_schurs->schur_explicit) {
739: /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
740: it is not efficient, unless the economic version of the scaling is used */
741: Mat S_Ej_expl;
742: PetscScalar *work;
743: PetscInt j,*dummy_idx;
744: PetscBool Sdense;
746: PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
747: local_size = 0;
748: for (i=0;i<sub_schurs->n_subs;i++) {
749: IS is_subset_B;
750: Mat AE_EE,AE_IE,AE_EI,S_Ej;
752: /* subsets in original and boundary numbering */
753: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);
754: /* EE block */
755: MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);
756: /* IE block */
757: MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);
758: /* EI block */
759: if (sub_schurs->is_symmetric) {
760: MatCreateTranspose(AE_IE,&AE_EI);
761: } else if (sub_schurs->is_hermitian) {
762: MatCreateHermitianTranspose(AE_IE,&AE_EI);
763: } else {
764: MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);
765: }
766: ISDestroy(&is_subset_B);
767: MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);
768: MatDestroy(&AE_EE);
769: MatDestroy(&AE_IE);
770: MatDestroy(&AE_EI);
771: if (AE_II == A_II) { /* we can reuse the same ksp */
772: KSP ksp;
773: MatSchurComplementGetKSP(sub_schurs->S,&ksp);
774: MatSchurComplementSetKSP(S_Ej,ksp);
775: } else { /* build new ksp object which inherits ksp and pc types from the original one */
776: KSP origksp,schurksp;
777: PC origpc,schurpc;
778: KSPType ksp_type;
779: PetscInt n_internal;
780: PetscBool ispcnone;
782: MatSchurComplementGetKSP(sub_schurs->S,&origksp);
783: MatSchurComplementGetKSP(S_Ej,&schurksp);
784: KSPGetType(origksp,&ksp_type);
785: KSPSetType(schurksp,ksp_type);
786: KSPGetPC(schurksp,&schurpc);
787: KSPGetPC(origksp,&origpc);
788: PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);
789: if (!ispcnone) {
790: PCType pc_type;
791: PCGetType(origpc,&pc_type);
792: PCSetType(schurpc,pc_type);
793: } else {
794: PCSetType(schurpc,PCLU);
795: }
796: ISGetSize(is_I,&n_internal);
797: if (n_internal) { /* UMFPACK gives error with 0 sized problems */
798: MatSolverType solver = NULL;
799: PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);
800: if (solver) {
801: PCFactorSetMatSolverType(schurpc,solver);
802: }
803: }
804: KSPSetUp(schurksp);
805: }
806: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
807: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);
808: PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);
809: PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);
810: if (Sdense) {
811: for (j=0;j<subset_size;j++) {
812: dummy_idx[j]=local_size+j;
813: }
814: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
815: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
816: MatDestroy(&S_Ej);
817: MatDestroy(&S_Ej_expl);
818: local_size += subset_size;
819: }
820: PetscFree2(dummy_idx,work);
821: /* free */
822: ISDestroy(&is_I);
823: MatDestroy(&AE_II);
824: PetscFree(all_local_idx_N);
825: } else {
826: Mat A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
827: Vec Dall = NULL;
828: IS is_A_all,*is_p_r = NULL;
829: PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
830: PetscInt n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
831: PetscBool economic,solver_S,S_lower_triangular = PETSC_FALSE;
832: PetscBool schur_has_vertices,factor_workaround;
833: PetscBool use_cholesky;
835: /* get sizes */
836: n_I = 0;
837: if (is_I_layer) {
838: ISGetLocalSize(is_I_layer,&n_I);
839: }
840: economic = PETSC_FALSE;
841: ISGetLocalSize(sub_schurs->is_I,&cum);
842: if (cum != n_I) economic = PETSC_TRUE;
843: MatGetLocalSize(sub_schurs->A,&n,NULL);
844: size_active_schur = local_size;
846: /* import scaling vector (wrong formulation if we have 3D edges) */
847: if (scaling && compute_Stilda) {
848: const PetscScalar *array;
849: PetscScalar *array2;
850: const PetscInt *idxs;
851: PetscInt i;
853: ISGetIndices(sub_schurs->is_Ej_all,&idxs);
854: VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);
855: VecGetArrayRead(scaling,&array);
856: VecGetArray(Dall,&array2);
857: for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
858: VecRestoreArray(Dall,&array2);
859: VecRestoreArrayRead(scaling,&array);
860: ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);
861: deluxe = PETSC_FALSE;
862: }
864: /* size active schurs does not count any dirichlet or vertex dof on the interface */
865: factor_workaround = PETSC_FALSE;
866: schur_has_vertices = PETSC_FALSE;
867: cum = n_I+size_active_schur;
868: if (sub_schurs->is_dir) {
869: const PetscInt* idxs;
870: PetscInt n_dir;
872: ISGetLocalSize(sub_schurs->is_dir,&n_dir);
873: ISGetIndices(sub_schurs->is_dir,&idxs);
874: PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));
875: ISRestoreIndices(sub_schurs->is_dir,&idxs);
876: cum += n_dir;
877: factor_workaround = PETSC_TRUE;
878: }
879: /* include the primal vertices in the Schur complement */
880: if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
881: PetscInt n_v;
883: ISGetLocalSize(sub_schurs->is_vertices,&n_v);
884: if (n_v) {
885: const PetscInt* idxs;
887: ISGetIndices(sub_schurs->is_vertices,&idxs);
888: PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));
889: ISRestoreIndices(sub_schurs->is_vertices,&idxs);
890: cum += n_v;
891: factor_workaround = PETSC_TRUE;
892: schur_has_vertices = PETSC_TRUE;
893: }
894: }
895: size_schur = cum - n_I;
896: ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);
897: if (cum == n) {
898: ISSetPermutation(is_A_all);
899: MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);
900: } else {
901: MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);
902: }
903: MatSetOptionsPrefix(A,sub_schurs->prefix);
904: MatAppendOptionsPrefix(A,"sub_schurs_");
906: /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
907: n^2 instead of n^1.5 or something. This is a workaround */
908: if (benign_n) {
909: Vec v;
910: ISLocalToGlobalMapping N_to_reor;
911: IS is_p0,is_p0_p;
912: PetscScalar *cs_AIB,*AIIm1_data;
913: PetscInt sizeA;
915: ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);
916: ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);
917: ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);
918: ISDestroy(&is_p0);
919: MatCreateVecs(A,&v,NULL);
920: VecGetSize(v,&sizeA);
921: MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);
922: MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);
923: MatDenseGetArray(cs_AIB_mat,&cs_AIB);
924: MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
925: PetscMalloc1(benign_n,&is_p_r);
926: /* compute colsum of A_IB restricted to pressures */
927: for (i=0;i<benign_n;i++) {
928: Vec benign_AIIm1_ones;
929: PetscScalar *array;
930: const PetscInt *idxs;
931: PetscInt j,nz;
933: ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);
934: ISGetLocalSize(is_p_r[i],&nz);
935: ISGetIndices(is_p_r[i],&idxs);
936: for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
937: ISRestoreIndices(is_p_r[i],&idxs);
938: VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);
939: MatMult(A,benign_AIIm1_ones,v);
940: VecDestroy(&benign_AIIm1_ones);
941: VecGetArray(v,&array);
942: for (j=0;j<size_schur;j++) {
943: #if defined(PETSC_USE_COMPLEX)
944: cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
945: #else
946: cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
947: #endif
948: }
949: VecRestoreArray(v,&array);
950: }
951: MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
952: MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
953: VecDestroy(&v);
954: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);
955: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
956: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
957: MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);
958: ISDestroy(&is_p0_p);
959: ISLocalToGlobalMappingDestroy(&N_to_reor);
960: }
961: MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);
962: MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);
963: MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);
965: /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
966: use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
968: /* when using the benign subspace trick, the local Schur complements are SPD */
969: if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
971: if (n_I) {
972: IS is_schur;
974: if (use_cholesky) {
975: MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);
976: } else {
977: MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);
978: }
979: MatSetErrorIfFailure(A,PETSC_TRUE);
981: /* subsets ordered last */
982: ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);
983: MatFactorSetSchurIS(F,is_schur);
984: ISDestroy(&is_schur);
986: /* factorization step */
987: if (use_cholesky) {
988: MatCholeskyFactorSymbolic(F,A,NULL,NULL);
989: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
990: MatMumpsSetIcntl(F,19,2);
991: #endif
992: MatCholeskyFactorNumeric(F,A,NULL);
993: S_lower_triangular = PETSC_TRUE;
994: } else {
995: MatLUFactorSymbolic(F,A,NULL,NULL,NULL);
996: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
997: MatMumpsSetIcntl(F,19,3);
998: #endif
999: MatLUFactorNumeric(F,A,NULL);
1000: }
1001: MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");
1003: if (print_schurs) {
1004: PetscViewer viewer;
1005: char filename[256];
1006: Mat S;
1007: IS is;
1009: PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);
1010: PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);
1011: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
1012: PetscObjectSetName((PetscObject)A,"A");
1013: MatView(A,viewer);
1014: MatFactorCreateSchurComplement(F,&S,NULL);
1015: PetscObjectSetName((PetscObject)S,"S");
1016: MatView(S,viewer);
1017: MatDestroy(&S);
1018: ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);
1019: PetscObjectSetName((PetscObject)is,"I");
1020: ISView(is,viewer);
1021: ISDestroy(&is);
1022: ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);
1023: PetscObjectSetName((PetscObject)is,"B");
1024: ISView(is,viewer);
1025: ISDestroy(&is);
1026: PetscViewerDestroy(&viewer);
1027: }
1029: /* get explicit Schur Complement computed during numeric factorization */
1030: MatFactorGetSchurComplement(F,&S_all,NULL);
1031: MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);
1032: MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);
1034: /* we can reuse the solvers if we are not using the economic version */
1035: reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1036: factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1037: if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
1038: reuse_solvers = factor_workaround = PETSC_FALSE;
1040: solver_S = PETSC_TRUE;
1042: /* update the Schur complement with the change of basis on the pressures */
1043: if (benign_n) {
1044: PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1045: Vec v;
1046: PetscInt sizeA;
1048: MatDenseGetArray(S_all,&S_data);
1049: MatCreateVecs(A,&v,NULL);
1050: VecGetSize(v,&sizeA);
1051: #if defined(PETSC_HAVE_MUMPS)
1052: MatMumpsSetIcntl(F,26,0);
1053: #endif
1054: #if defined(PETSC_HAVE_MKL_PARDISO)
1055: MatMkl_PardisoSetCntl(F,70,1);
1056: #endif
1057: MatDenseGetArray(cs_AIB_mat,&cs_AIB);
1058: MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
1059: for (i=0;i<benign_n;i++) {
1060: Vec benign_AIIm1_ones;
1061: PetscScalar *array,sum = 0.,one = 1.;
1062: const PetscInt *idxs;
1063: PetscInt j,nz;
1064: PetscBLASInt B_k,B_n;
1066: VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);
1067: VecCopy(benign_AIIm1_ones,v);
1068: MatSolve(F,v,benign_AIIm1_ones);
1069: ISGetLocalSize(is_p_r[i],&nz);
1070: ISGetIndices(is_p_r[i],&idxs);
1071: /* p0 dof (eliminated) is excluded from the sum */
1072: for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1073: ISRestoreIndices(is_p_r[i],&idxs);
1074: MatMult(A,benign_AIIm1_ones,v);
1075: VecGetArray(v,&array);
1076: /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1077: /* cs_AIB already scaled by 1./nz */
1078: B_k = 1;
1079: PetscBLASIntCast(size_schur,&B_n);
1080: PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1081: sum = 1.;
1082: 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));
1083: VecRestoreArray(v,&array);
1084: /* set p0 entry of AIIm1_ones to zero */
1085: ISGetIndices(is_p_r[i],&idxs);
1086: for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1087: ISRestoreIndices(is_p_r[i],&idxs);
1088: VecDestroy(&benign_AIIm1_ones);
1089: }
1090: if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1091: PetscInt k,j;
1092: for (k=0;k<size_schur;k++) {
1093: for (j=k;j<size_schur;j++) {
1094: S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1095: }
1096: }
1097: }
1099: /* restore defaults */
1100: #if defined(PETSC_HAVE_MUMPS)
1101: MatMumpsSetIcntl(F,26,-1);
1102: #endif
1103: #if defined(PETSC_HAVE_MKL_PARDISO)
1104: MatMkl_PardisoSetCntl(F,70,0);
1105: #endif
1106: MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
1107: MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
1108: VecDestroy(&v);
1109: MatDenseRestoreArray(S_all,&S_data);
1110: }
1111: if (!reuse_solvers) {
1112: for (i=0;i<benign_n;i++) {
1113: ISDestroy(&is_p_r[i]);
1114: }
1115: PetscFree(is_p_r);
1116: MatDestroy(&cs_AIB_mat);
1117: MatDestroy(&benign_AIIm1_ones_mat);
1118: }
1119: } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1120: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);
1121: reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1122: factor_workaround = PETSC_FALSE;
1123: solver_S = PETSC_FALSE;
1124: }
1126: if (reuse_solvers) {
1127: Mat A_II,Afake;
1128: Vec vec1_B;
1129: PCBDDCReuseSolvers msolv_ctx;
1130: PetscInt n_R;
1132: if (sub_schurs->reuse_solver) {
1133: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1134: } else {
1135: PetscNew(&sub_schurs->reuse_solver);
1136: }
1137: msolv_ctx = sub_schurs->reuse_solver;
1138: MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);
1139: PetscObjectReference((PetscObject)F);
1140: msolv_ctx->F = F;
1141: MatCreateVecs(F,&msolv_ctx->sol,NULL);
1142: /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1143: {
1144: PetscScalar *array;
1145: PetscInt n;
1147: VecGetLocalSize(msolv_ctx->sol,&n);
1148: VecGetArray(msolv_ctx->sol,&array);
1149: VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);
1150: VecRestoreArray(msolv_ctx->sol,&array);
1151: }
1152: msolv_ctx->has_vertices = schur_has_vertices;
1154: /* interior solver */
1155: PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);
1156: PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);
1157: PCSetType(msolv_ctx->interior_solver,PCSHELL);
1158: PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");
1159: PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);
1160: PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);
1161: PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);
1162: PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);
1164: /* correction solver */
1165: PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);
1166: PCSetType(msolv_ctx->correction_solver,PCSHELL);
1167: PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");
1168: PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);
1169: PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);
1170: PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);
1171: PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);
1173: /* scatter and vecs for Schur complement solver */
1174: MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);
1175: MatCreateVecs(sub_schurs->S,&vec1_B,NULL);
1176: if (!schur_has_vertices) {
1177: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);
1178: VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);
1179: PetscObjectReference((PetscObject)is_A_all);
1180: msolv_ctx->is_R = is_A_all;
1181: } else {
1182: IS is_B_all;
1183: const PetscInt* idxs;
1184: PetscInt dual,n_v,n;
1186: ISGetLocalSize(sub_schurs->is_vertices,&n_v);
1187: dual = size_schur - n_v;
1188: ISGetLocalSize(is_A_all,&n);
1189: ISGetIndices(is_A_all,&idxs);
1190: ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);
1191: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);
1192: ISDestroy(&is_B_all);
1193: ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);
1194: VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);
1195: ISDestroy(&is_B_all);
1196: ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);
1197: ISRestoreIndices(is_A_all,&idxs);
1198: }
1199: ISGetLocalSize(msolv_ctx->is_R,&n_R);
1200: MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);
1201: MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);
1202: MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);
1203: PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);
1204: MatDestroy(&Afake);
1205: VecDestroy(&vec1_B);
1207: /* communicate benign info to solver context */
1208: if (benign_n) {
1209: PetscScalar *array;
1211: msolv_ctx->benign_n = benign_n;
1212: msolv_ctx->benign_zerodiag_subs = is_p_r;
1213: PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);
1214: msolv_ctx->benign_csAIB = cs_AIB_mat;
1215: MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);
1216: VecGetArray(msolv_ctx->benign_corr_work,&array);
1217: VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);
1218: VecRestoreArray(msolv_ctx->benign_corr_work,&array);
1219: msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1220: }
1221: } else {
1222: if (sub_schurs->reuse_solver) {
1223: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1224: }
1225: PetscFree(sub_schurs->reuse_solver);
1226: }
1227: MatDestroy(&A);
1228: ISDestroy(&is_A_all);
1230: /* Work arrays */
1231: PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
1233: /* matrices for deluxe scaling and adaptive selection */
1234: if (compute_Stilda) {
1235: if (!sub_schurs->sum_S_Ej_tilda_all) {
1236: MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);
1237: }
1238: if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1239: MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);
1240: }
1241: }
1243: /* S_Ej_all */
1244: cum = cum2 = 0;
1245: MatDenseGetArray(S_all,&S_data);
1246: for (i=0;i<sub_schurs->n_subs;i++) {
1247: PetscInt j;
1249: /* get S_E */
1250: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1251: if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1252: PetscInt k;
1253: for (k=0;k<subset_size;k++) {
1254: for (j=k;j<subset_size;j++) {
1255: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1256: work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1257: }
1258: }
1259: } else { /* just copy to workspace */
1260: PetscInt k;
1261: for (k=0;k<subset_size;k++) {
1262: for (j=0;j<subset_size;j++) {
1263: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1264: }
1265: }
1266: }
1267: /* insert S_E values */
1268: for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1269: if (sub_schurs->change) {
1270: Mat change_sub,SEj,T;
1272: /* change basis */
1273: KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1274: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1275: if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1276: Mat T2;
1277: MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1278: MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1279: MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1280: MatDestroy(&T2);
1281: } else {
1282: MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1283: }
1284: MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1285: MatDestroy(&T);
1286: MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);
1287: MatDestroy(&SEj);
1288: }
1289: if (deluxe) {
1290: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1291: /* if adaptivity is requested, invert S_E blocks */
1292: if (compute_Stilda) {
1293: PetscInt k;
1295: PetscBLASIntCast(subset_size,&B_N);
1296: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1297: if (use_potr) {
1298: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1299: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1300: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1301: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1302: for (k=0;k<subset_size;k++) {
1303: for (j=k;j<subset_size;j++) {
1304: work[j*subset_size+k] = work[k*subset_size+j];
1305: }
1306: }
1307: } else if (use_sytr) {
1308: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1309: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1310: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,work,&B_N,pivots,Bwork,&B_ierr));
1311: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1312: for (k=0;k<subset_size;k++) {
1313: for (j=k;j<subset_size;j++) {
1314: work[j*subset_size+k] = work[k*subset_size+j];
1315: }
1316: }
1317: } else {
1318: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1319: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1320: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1321: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1322: }
1323: PetscLogFlops(1.0*subset_size*subset_size*subset_size);
1324: PetscFPTrapPop();
1325: MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1326: }
1327: } else if (compute_Stilda) { /* not using deluxe */
1328: Mat SEj;
1329: Vec D;
1330: PetscScalar *array;
1332: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1333: VecGetArray(Dall,&array);
1334: VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);
1335: VecRestoreArray(Dall,&array);
1336: VecShift(D,-1.);
1337: MatDiagonalScale(SEj,D,D);
1338: MatDestroy(&SEj);
1339: VecDestroy(&D);
1340: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1341: }
1342: cum += subset_size;
1343: cum2 += subset_size*(size_schur + 1);
1344: }
1345: MatDenseRestoreArray(S_all,&S_data);
1347: if (solver_S) {
1348: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1349: }
1351: schur_factor = NULL;
1352: if (compute_Stilda && size_active_schur) {
1354: if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1355: PetscInt j;
1356: for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1357: MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);
1358: } else {
1359: Mat S_all_inv=NULL;
1360: if (solver_S) {
1361: /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1362: 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 */
1363: if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1364: PetscScalar *data;
1365: PetscInt nd = 0;
1367: if (!use_potr) {
1368: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1369: }
1370: MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1371: MatDenseGetArray(S_all_inv,&data);
1372: if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1373: ISGetLocalSize(sub_schurs->is_dir,&nd);
1374: }
1376: /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1377: if (schur_has_vertices) {
1378: Mat M;
1379: PetscScalar *tdata;
1380: PetscInt nv = 0, news;
1382: ISGetLocalSize(sub_schurs->is_vertices,&nv);
1383: news = size_active_schur + nv;
1384: PetscCalloc1(news*news,&tdata);
1385: for (i=0;i<size_active_schur;i++) {
1386: PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));
1387: PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));
1388: }
1389: for (i=0;i<nv;i++) {
1390: PetscInt k = i+size_active_schur;
1391: PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));
1392: }
1394: MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);
1395: MatSetOption(M,MAT_SPD,PETSC_TRUE);
1396: MatCholeskyFactor(M,NULL,NULL);
1397: /* save the factors */
1398: cum = 0;
1399: PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);
1400: for (i=0;i<size_active_schur;i++) {
1401: PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));
1402: cum += size_active_schur - i;
1403: }
1404: for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1405: MatSeqDenseInvertFactors_Private(M);
1406: /* move back just the active dofs to the Schur complement */
1407: for (i=0;i<size_active_schur;i++) {
1408: PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));
1409: }
1410: PetscFree(tdata);
1411: MatDestroy(&M);
1412: } else { /* we can factorize and invert just the activedofs */
1413: Mat M;
1414: PetscScalar *aux;
1416: PetscMalloc1(nd,&aux);
1417: for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1418: MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);
1419: MatSeqDenseSetLDA(M,size_schur);
1420: MatSetOption(M,MAT_SPD,PETSC_TRUE);
1421: MatCholeskyFactor(M,NULL,NULL);
1422: MatSeqDenseInvertFactors_Private(M);
1423: MatDestroy(&M);
1424: MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);
1425: MatZeroEntries(M);
1426: MatDestroy(&M);
1427: MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);
1428: MatSeqDenseSetLDA(M,size_schur);
1429: MatZeroEntries(M);
1430: MatDestroy(&M);
1431: for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
1432: PetscFree(aux);
1433: }
1434: MatDenseRestoreArray(S_all_inv,&data);
1435: } else { /* use MatFactor calls to invert S */
1436: MatFactorInvertSchurComplement(F);
1437: MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1438: }
1439: } else { /* we need to invert explicitly since we are not using MatFactor for S */
1440: PetscObjectReference((PetscObject)S_all);
1441: S_all_inv = S_all;
1442: MatDenseGetArray(S_all_inv,&S_data);
1443: PetscBLASIntCast(size_schur,&B_N);
1444: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1445: if (use_potr) {
1446: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1447: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1448: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1449: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1450: } else if (use_sytr) {
1451: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1452: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1453: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1454: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1455: } else {
1456: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1457: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1458: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1459: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1460: }
1461: PetscLogFlops(1.0*size_schur*size_schur*size_schur);
1462: PetscFPTrapPop();
1463: MatDenseRestoreArray(S_all_inv,&S_data);
1464: }
1465: /* S_Ej_tilda_all */
1466: cum = cum2 = 0;
1467: MatDenseGetArray(S_all_inv,&S_data);
1468: for (i=0;i<sub_schurs->n_subs;i++) {
1469: PetscInt j;
1471: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1472: /* get (St^-1)_E */
1473: /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1474: will be properly accessed later during adaptive selection */
1475: if (S_lower_triangular) {
1476: PetscInt k;
1477: if (sub_schurs->change) {
1478: for (k=0;k<subset_size;k++) {
1479: for (j=k;j<subset_size;j++) {
1480: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1481: work[j*subset_size+k] = work[k*subset_size+j];
1482: }
1483: }
1484: } else {
1485: for (k=0;k<subset_size;k++) {
1486: for (j=k;j<subset_size;j++) {
1487: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1488: }
1489: }
1490: }
1491: } else {
1492: PetscInt k;
1493: for (k=0;k<subset_size;k++) {
1494: for (j=0;j<subset_size;j++) {
1495: work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1496: }
1497: }
1498: }
1499: if (sub_schurs->change) {
1500: Mat change_sub,SEj,T;
1502: /* change basis */
1503: KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1504: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1505: if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1506: Mat T2;
1507: MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1508: MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1509: MatDestroy(&T2);
1510: MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1511: } else {
1512: MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1513: }
1514: MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1515: MatDestroy(&T);
1516: /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1517: MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);
1518: MatDestroy(&SEj);
1519: }
1520: for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1521: MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1522: cum += subset_size;
1523: cum2 += subset_size*(size_schur + 1);
1524: }
1525: MatDenseRestoreArray(S_all_inv,&S_data);
1526: if (solver_S) {
1527: if (schur_has_vertices) {
1528: MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);
1529: } else {
1530: MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);
1531: }
1532: }
1533: MatDestroy(&S_all_inv);
1534: }
1536: /* move back factors if needed */
1537: if (schur_has_vertices) {
1538: Mat S_tmp;
1539: PetscInt nd = 0;
1541: if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1542: MatFactorGetSchurComplement(F,&S_tmp,NULL);
1543: if (use_potr) {
1544: PetscScalar *data;
1546: MatDenseGetArray(S_tmp,&data);
1547: PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));
1549: if (S_lower_triangular) {
1550: cum = 0;
1551: for (i=0;i<size_active_schur;i++) {
1552: PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));
1553: cum += size_active_schur-i;
1554: }
1555: } else {
1556: PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));
1557: }
1558: if (sub_schurs->is_dir) {
1559: ISGetLocalSize(sub_schurs->is_dir,&nd);
1560: for (i=0;i<nd;i++) {
1561: data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1562: }
1563: }
1564: /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1565: set the diagonal entry of the Schur factor to a very large value */
1566: for (i=size_active_schur+nd;i<size_schur;i++) {
1567: data[i*(size_schur+1)] = infty;
1568: }
1569: MatDenseRestoreArray(S_tmp,&data);
1570: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1571: MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);
1572: }
1573: } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1574: PetscScalar *data;
1575: PetscInt nd = 0;
1577: if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1578: ISGetLocalSize(sub_schurs->is_dir,&nd);
1579: }
1580: MatFactorGetSchurComplement(F,&S_all,NULL);
1581: MatDenseGetArray(S_all,&data);
1582: for (i=0;i<size_active_schur;i++) {
1583: PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));
1584: }
1585: for (i=size_active_schur+nd;i<size_schur;i++) {
1586: PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));
1587: data[i*(size_schur+1)] = infty;
1588: }
1589: MatDenseRestoreArray(S_all,&data);
1590: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1591: }
1592: PetscFree2(dummy_idx,work);
1593: PetscFree(schur_factor);
1594: VecDestroy(&Dall);
1595: }
1596: PetscFree(nnz);
1597: MatDestroy(&F);
1598: ISDestroy(&is_I_layer);
1599: MatDestroy(&S_all);
1600: MatDestroy(&A_BB);
1601: MatDestroy(&A_IB);
1602: MatDestroy(&A_BI);
1604: /* speed up matrix assembly */
1605: PetscMalloc1(sub_schurs->n_subs,&nnz);
1606: for (i=0;i<sub_schurs->n_subs;i++) {
1607: ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);
1608: }
1609: ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);
1610: MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);
1611: MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1612: MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1613: if (compute_Stilda) {
1614: MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);
1615: MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1616: MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1617: if (deluxe) {
1618: MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);
1619: MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1620: MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1621: }
1622: }
1623: ISDestroy(&is_I_layer);
1625: /* Global matrix of all assembled Schur on subsets */
1626: MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);
1627: MatAssemblyBegin(work_mat,MAT_FINAL_ASSEMBLY);
1628: MatAssemblyEnd(work_mat,MAT_FINAL_ASSEMBLY);
1629: MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);
1630: MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);
1632: /* Get local part of (\sum_j S_Ej) */
1633: MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);
1634: if (!sub_schurs->sum_S_Ej_all) {
1635: MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);
1636: } else {
1637: MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);
1638: }
1640: /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1641: if (compute_Stilda) {
1642: MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);
1643: MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);
1644: MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);
1645: MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);
1646: if (deluxe) {
1647: MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);
1648: MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);
1649: MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);
1650: MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);
1651: } else {
1652: PetscScalar *array;
1653: PetscInt cum;
1655: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1656: cum = 0;
1657: for (i=0;i<sub_schurs->n_subs;i++) {
1658: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1659: PetscBLASIntCast(subset_size,&B_N);
1660: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1661: if (use_potr) {
1662: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1663: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1664: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1665: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1666: } else if (use_sytr) {
1667: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1668: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1669: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1670: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1671: } else {
1672: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1673: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1674: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1675: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1676: }
1677: PetscLogFlops(1.0*subset_size*subset_size*subset_size);
1678: PetscFPTrapPop();
1679: cum += subset_size*subset_size;
1680: }
1681: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1682: PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);
1683: MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1684: sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1685: }
1686: }
1687: MatDestroySubMatrices(1,&submats);
1688: if (print_schurs) {
1689: PetscViewer viewer;
1690: char filename[256];
1691: PetscInt cum;
1693: PetscSNPrintf(filename,sizeof(filename),"sub_schurs_mats_r%d.m",PetscGlobalRank);
1694: PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);
1695: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
1696: if (sub_schurs->S_Ej_all) {
1697: PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");
1698: MatView(sub_schurs->S_Ej_all,viewer);
1699: }
1700: if (sub_schurs->sum_S_Ej_all) {
1701: PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");
1702: MatView(sub_schurs->sum_S_Ej_all,viewer);
1703: }
1704: if (sub_schurs->sum_S_Ej_inv_all) {
1705: PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");
1706: MatView(sub_schurs->sum_S_Ej_inv_all,viewer);
1707: }
1708: if (sub_schurs->sum_S_Ej_tilda_all) {
1709: PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");
1710: MatView(sub_schurs->sum_S_Ej_tilda_all,viewer);
1711: }
1712: for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
1713: IS is;
1714: char name[16];
1716: PetscSNPrintf(name,sizeof(name),"IE%D",i);
1717: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1718: ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);
1719: PetscObjectSetName((PetscObject)is,name);
1720: ISView(is,viewer);
1721: ISDestroy(&is);
1722: cum += subset_size;
1723: }
1724: PetscViewerDestroy(&viewer);
1725: }
1727: /* free workspace */
1728: PetscFree(submats);
1729: PetscFree2(Bwork,pivots);
1730: MatDestroy(&global_schur_subsets);
1731: MatDestroy(&work_mat);
1732: ISDestroy(&all_subsets_n);
1733: PetscCommDestroy(&comm_n);
1734: return(0);
1735: }
1737: PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1738: {
1739: IS *faces,*edges,*all_cc,vertices;
1740: PetscInt i,n_faces,n_edges,n_all_cc;
1741: PetscBool is_sorted,ispetsc;
1742: PetscErrorCode ierr;
1745: ISSorted(is_I,&is_sorted);
1746: if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1747: ISSorted(is_B,&is_sorted);
1748: if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1750: /* reset any previous data */
1751: PCBDDCSubSchursReset(sub_schurs);
1753: /* get index sets for faces and edges (already sorted by global ordering) */
1754: PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1755: n_all_cc = n_faces+n_edges;
1756: PetscBTCreate(n_all_cc,&sub_schurs->is_edge);
1757: PetscMalloc1(n_all_cc,&all_cc);
1758: for (i=0;i<n_faces;i++) {
1759: if (copycc) {
1760: ISDuplicate(faces[i],&all_cc[i]);
1761: } else {
1762: PetscObjectReference((PetscObject)faces[i]);
1763: all_cc[i] = faces[i];
1764: }
1765: }
1766: for (i=0;i<n_edges;i++) {
1767: if (copycc) {
1768: ISDuplicate(edges[i],&all_cc[n_faces+i]);
1769: } else {
1770: PetscObjectReference((PetscObject)edges[i]);
1771: all_cc[n_faces+i] = edges[i];
1772: }
1773: PetscBTSet(sub_schurs->is_edge,n_faces+i);
1774: }
1775: PetscObjectReference((PetscObject)vertices);
1776: sub_schurs->is_vertices = vertices;
1777: PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1778: sub_schurs->is_dir = NULL;
1779: PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);
1781: /* Determine if MatFactor can be used */
1782: PetscStrallocpy(prefix,&sub_schurs->prefix);
1783: #if defined(PETSC_HAVE_MUMPS)
1784: PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);
1785: #elif defined(PETSC_HAVE_MKL_PARDISO)
1786: PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);
1787: #else
1788: PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);
1789: #endif
1790: #if defined(PETSC_USE_COMPLEX)
1791: sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1792: #else
1793: sub_schurs->is_hermitian = PETSC_TRUE;
1794: #endif
1795: sub_schurs->is_posdef = PETSC_TRUE;
1796: sub_schurs->is_symmetric = PETSC_TRUE;
1797: sub_schurs->debug = PETSC_FALSE;
1798: sub_schurs->restrict_comm = PETSC_FALSE;
1799: PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
1800: PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,64,NULL);
1801: PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);
1802: PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);
1803: PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);
1804: PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);
1805: PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);
1806: PetscOptionsEnd();
1807: PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);
1808: sub_schurs->schur_explicit = (PetscBool)!ispetsc;
1810: /* for reals, symmetric and hermitian are synonims */
1811: #if !defined(PETSC_USE_COMPLEX)
1812: sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1813: sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1814: #endif
1816: PetscObjectReference((PetscObject)is_I);
1817: sub_schurs->is_I = is_I;
1818: PetscObjectReference((PetscObject)is_B);
1819: sub_schurs->is_B = is_B;
1820: PetscObjectReference((PetscObject)graph->l2gmap);
1821: sub_schurs->l2gmap = graph->l2gmap;
1822: PetscObjectReference((PetscObject)BtoNmap);
1823: sub_schurs->BtoNmap = BtoNmap;
1824: sub_schurs->n_subs = n_all_cc;
1825: sub_schurs->is_subs = all_cc;
1826: sub_schurs->S_Ej_all = NULL;
1827: sub_schurs->sum_S_Ej_all = NULL;
1828: sub_schurs->sum_S_Ej_inv_all = NULL;
1829: sub_schurs->sum_S_Ej_tilda_all = NULL;
1830: sub_schurs->is_Ej_all = NULL;
1831: return(0);
1832: }
1834: PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1835: {
1836: PCBDDCSubSchurs schurs_ctx;
1837: PetscErrorCode ierr;
1840: PetscNew(&schurs_ctx);
1841: schurs_ctx->n_subs = 0;
1842: *sub_schurs = schurs_ctx;
1843: return(0);
1844: }
1846: PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1847: {
1848: PetscInt i;
1852: if (!sub_schurs) return(0);
1853: PetscFree(sub_schurs->prefix);
1854: MatDestroy(&sub_schurs->A);
1855: MatDestroy(&sub_schurs->S);
1856: ISDestroy(&sub_schurs->is_I);
1857: ISDestroy(&sub_schurs->is_B);
1858: ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);
1859: ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);
1860: MatDestroy(&sub_schurs->S_Ej_all);
1861: MatDestroy(&sub_schurs->sum_S_Ej_all);
1862: MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1863: MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);
1864: ISDestroy(&sub_schurs->is_Ej_all);
1865: ISDestroy(&sub_schurs->is_vertices);
1866: ISDestroy(&sub_schurs->is_dir);
1867: PetscBTDestroy(&sub_schurs->is_edge);
1868: for (i=0;i<sub_schurs->n_subs;i++) {
1869: ISDestroy(&sub_schurs->is_subs[i]);
1870: }
1871: if (sub_schurs->n_subs) {
1872: PetscFree(sub_schurs->is_subs);
1873: }
1874: if (sub_schurs->reuse_solver) {
1875: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1876: }
1877: PetscFree(sub_schurs->reuse_solver);
1878: if (sub_schurs->change) {
1879: for (i=0;i<sub_schurs->n_subs;i++) {
1880: KSPDestroy(&sub_schurs->change[i]);
1881: ISDestroy(&sub_schurs->change_primal_sub[i]);
1882: }
1883: }
1884: PetscFree(sub_schurs->change);
1885: PetscFree(sub_schurs->change_primal_sub);
1886: sub_schurs->n_subs = 0;
1887: return(0);
1888: }
1890: PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1891: {
1895: PCBDDCSubSchursReset(*sub_schurs);
1896: PetscFree(*sub_schurs);
1897: return(0);
1898: }
1900: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1901: {
1902: PetscInt i,j,n;
1906: n = 0;
1907: for (i=-n_prev;i<0;i++) {
1908: PetscInt start_dof = queue_tip[i];
1909: for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1910: PetscInt dof = adjncy[j];
1911: if (!PetscBTLookup(touched,dof)) {
1912: PetscBTSet(touched,dof);
1913: queue_tip[n] = dof;
1914: n++;
1915: }
1916: }
1917: }
1918: *n_added = n;
1919: return(0);
1920: }