Actual source code: bddcschurs.c
petsc-3.13.1 2020-05-02
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: PetscArraycpy(array2,array,nl);
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: PetscArraycpy(array_solver,array,n);
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: PetscArraycpy(array,array_solver,n);
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: PCComputeOperator(pc, MATDENSE, &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: Vec gstash,lstash;
407: VecScatter sstash;
408: IS is_I,is_I_layer;
409: IS all_subsets,all_subsets_mult,all_subsets_n;
410: PetscScalar *stasharray,*Bwork;
411: PetscInt *nnz,*all_local_idx_N;
412: PetscInt *auxnum1,*auxnum2;
413: PetscInt i,subset_size,max_subset_size;
414: PetscInt n_B,extra,local_size,global_size;
415: PetscInt local_stash_size;
416: PetscBLASInt B_N,B_ierr,B_lwork,*pivots;
417: MPI_Comm comm_n;
418: PetscBool deluxe = PETSC_TRUE;
419: PetscBool use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
420: PetscViewer matl_dbg_viewer = NULL;
421: PetscErrorCode ierr;
424: MatDestroy(&sub_schurs->A);
425: MatDestroy(&sub_schurs->S);
426: if (Ain) {
427: PetscObjectReference((PetscObject)Ain);
428: sub_schurs->A = Ain;
429: }
431: PetscObjectReference((PetscObject)Sin);
432: sub_schurs->S = Sin;
433: if (sub_schurs->schur_explicit) {
434: sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
435: }
437: /* preliminary checks */
438: 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");
440: if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
442: /* debug (MATLAB) */
443: if (sub_schurs->debug) {
444: PetscMPIInt size,rank;
445: PetscInt nr,*print_schurs_ranks,print_schurs = PETSC_FALSE;
446: PetscBool flg;
448: MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);
449: MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
450: nr = size;
451: PetscMalloc1(nr,&print_schurs_ranks);
452: PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
453: PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);
454: if (!flg) print_schurs = PETSC_TRUE;
455: else {
456: print_schurs = PETSC_FALSE;
457: for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; }
458: }
459: PetscOptionsEnd();
460: PetscFree(print_schurs_ranks);
461: if (print_schurs) {
462: char filename[256];
464: PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);
465: PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&matl_dbg_viewer);
466: PetscViewerPushFormat(matl_dbg_viewer,PETSC_VIEWER_ASCII_MATLAB);
467: }
468: }
471: /* restrict work on active processes */
472: if (sub_schurs->restrict_comm) {
473: PetscSubcomm subcomm;
474: PetscMPIInt color,rank;
476: color = 0;
477: if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
478: MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
479: PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);
480: PetscSubcommSetNumber(subcomm,2);
481: PetscSubcommSetTypeGeneral(subcomm,color,rank);
482: PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);
483: PetscSubcommDestroy(&subcomm);
484: if (!sub_schurs->n_subs) {
485: PetscCommDestroy(&comm_n);
486: return(0);
487: }
488: } else {
489: PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);
490: }
492: /* get Schur complement matrices */
493: if (!sub_schurs->schur_explicit) {
494: Mat tA_IB,tA_BI,tA_BB;
495: PetscBool isseqsbaij;
496: MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);
497: PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);
498: if (isseqsbaij) {
499: MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);
500: MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);
501: MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);
502: } else {
503: PetscObjectReference((PetscObject)tA_BB);
504: A_BB = tA_BB;
505: PetscObjectReference((PetscObject)tA_IB);
506: A_IB = tA_IB;
507: PetscObjectReference((PetscObject)tA_BI);
508: A_BI = tA_BI;
509: }
510: } else {
511: A_II = NULL;
512: A_IB = NULL;
513: A_BI = NULL;
514: A_BB = NULL;
515: }
516: S_all = NULL;
518: /* determine interior problems */
519: ISGetLocalSize(sub_schurs->is_I,&i);
520: if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
521: PetscBT touched;
522: const PetscInt* idx_B;
523: PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
525: if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
526: /* get sizes */
527: ISGetLocalSize(sub_schurs->is_I,&n_I);
528: ISGetLocalSize(sub_schurs->is_B,&n_B);
530: PetscMalloc1(n_I+n_B,&local_numbering);
531: PetscBTCreate(n_I+n_B,&touched);
532: PetscBTMemzero(n_I+n_B,touched);
534: /* all boundary dofs must be skipped when adding layers */
535: ISGetIndices(sub_schurs->is_B,&idx_B);
536: for (j=0;j<n_B;j++) {
537: PetscBTSet(touched,idx_B[j]);
538: }
539: PetscArraycpy(local_numbering,idx_B,n_B);
540: ISRestoreIndices(sub_schurs->is_B,&idx_B);
542: /* add prescribed number of layers of dofs */
543: n_local_dofs = n_B;
544: n_prev_added = n_B;
545: for (layer=0;layer<nlayers;layer++) {
546: PetscInt n_added;
547: if (n_local_dofs == n_I+n_B) break;
548: 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);
549: PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);
550: n_prev_added = n_added;
551: n_local_dofs += n_added;
552: if (!n_added) break;
553: }
554: PetscBTDestroy(&touched);
556: /* IS for I layer dofs in original numbering */
557: ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);
558: PetscFree(local_numbering);
559: ISSort(is_I_layer);
560: /* IS for I layer dofs in I numbering */
561: if (!sub_schurs->schur_explicit) {
562: ISLocalToGlobalMapping ItoNmap;
563: ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);
564: ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);
565: ISLocalToGlobalMappingDestroy(&ItoNmap);
567: /* II block */
568: MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);
569: }
570: } else {
571: PetscInt n_I;
573: /* IS for I dofs in original numbering */
574: PetscObjectReference((PetscObject)sub_schurs->is_I);
575: is_I_layer = sub_schurs->is_I;
577: /* IS for I dofs in I numbering (strided 1) */
578: if (!sub_schurs->schur_explicit) {
579: ISGetSize(sub_schurs->is_I,&n_I);
580: ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);
582: /* II block is the same */
583: PetscObjectReference((PetscObject)A_II);
584: AE_II = A_II;
585: }
586: }
588: /* Get info on subset sizes and sum of all subsets sizes */
589: max_subset_size = 0;
590: local_size = 0;
591: for (i=0;i<sub_schurs->n_subs;i++) {
592: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
593: max_subset_size = PetscMax(subset_size,max_subset_size);
594: local_size += subset_size;
595: }
597: /* Work arrays for local indices */
598: extra = 0;
599: ISGetLocalSize(sub_schurs->is_B,&n_B);
600: if (sub_schurs->schur_explicit && is_I_layer) {
601: ISGetLocalSize(is_I_layer,&extra);
602: }
603: PetscMalloc1(n_B+extra,&all_local_idx_N);
604: if (extra) {
605: const PetscInt *idxs;
606: ISGetIndices(is_I_layer,&idxs);
607: PetscArraycpy(all_local_idx_N,idxs,extra);
608: ISRestoreIndices(is_I_layer,&idxs);
609: }
610: PetscMalloc1(sub_schurs->n_subs,&auxnum1);
611: PetscMalloc1(sub_schurs->n_subs,&auxnum2);
613: /* Get local indices in local numbering */
614: local_size = 0;
615: local_stash_size = 0;
616: for (i=0;i<sub_schurs->n_subs;i++) {
617: const PetscInt *idxs;
619: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
620: ISGetIndices(sub_schurs->is_subs[i],&idxs);
621: /* start (smallest in global ordering) and multiplicity */
622: auxnum1[i] = idxs[0];
623: auxnum2[i] = subset_size*subset_size;
624: /* subset indices in local numbering */
625: PetscArraycpy(all_local_idx_N+local_size+extra,idxs,subset_size);
626: ISRestoreIndices(sub_schurs->is_subs[i],&idxs);
627: local_size += subset_size;
628: local_stash_size += subset_size*subset_size;
629: }
631: /* allocate extra workspace needed only for GETRI or SYTRF */
632: use_potr = use_sytr = PETSC_FALSE;
633: if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
634: use_potr = PETSC_TRUE;
635: } else if (sub_schurs->is_symmetric) {
636: use_sytr = PETSC_TRUE;
637: }
638: if (local_size && !use_potr) {
639: PetscScalar lwork,dummyscalar = 0.;
640: PetscBLASInt dummyint = 0;
642: B_lwork = -1;
643: PetscBLASIntCast(local_size,&B_N);
644: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
645: if (use_sytr) {
646: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
647: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
648: } else {
649: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
650: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
651: }
652: PetscFPTrapPop();
653: PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);
654: PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);
655: } else {
656: Bwork = NULL;
657: pivots = NULL;
658: }
660: /* prepare data for summing up properly schurs on subsets */
661: ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);
662: ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);
663: ISDestroy(&all_subsets_n);
664: ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);
665: ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);
666: ISDestroy(&all_subsets);
667: ISDestroy(&all_subsets_mult);
668: ISGetLocalSize(all_subsets_n,&i);
669: if (i != local_stash_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_stash_size);
670: VecCreateSeqWithArray(PETSC_COMM_SELF,1,local_stash_size,NULL,&lstash);
671: VecCreateMPI(comm_n,PETSC_DECIDE,global_size,&gstash);
672: VecScatterCreate(lstash,NULL,gstash,all_subsets_n,&sstash);
673: ISDestroy(&all_subsets_n);
675: /* subset indices in local boundary numbering */
676: if (!sub_schurs->is_Ej_all) {
677: PetscInt *all_local_idx_B;
679: PetscMalloc1(local_size,&all_local_idx_B);
680: ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);
681: if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size);
682: ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);
683: }
685: if (change) {
686: ISLocalToGlobalMapping BtoS;
687: IS change_primal_B;
688: IS change_primal_all;
690: if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
691: if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
692: PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);
693: for (i=0;i<sub_schurs->n_subs;i++) {
694: ISLocalToGlobalMapping NtoS;
695: ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);
696: ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);
697: ISLocalToGlobalMappingDestroy(&NtoS);
698: }
699: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);
700: ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);
701: ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);
702: ISLocalToGlobalMappingDestroy(&BtoS);
703: ISDestroy(&change_primal_B);
704: PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);
705: for (i=0;i<sub_schurs->n_subs;i++) {
706: Mat change_sub;
708: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
709: KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);
710: KSPSetType(sub_schurs->change[i],KSPPREONLY);
711: if (!sub_schurs->change_with_qr) {
712: MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);
713: } else {
714: Mat change_subt;
715: MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);
716: MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);
717: MatDestroy(&change_subt);
718: }
719: KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);
720: MatDestroy(&change_sub);
721: KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);
722: KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");
723: }
724: ISDestroy(&change_primal_all);
725: }
727: /* Local matrix of all local Schur on subsets (transposed) */
728: if (!sub_schurs->S_Ej_all) {
729: Mat T;
730: PetscScalar *v;
731: PetscInt *ii,*jj;
732: PetscInt cum,i,j,k;
734: /* MatSeqAIJSetPreallocation + MatSetValues is slow for these kind of matrices (may have large blocks)
735: Allocate properly a representative matrix and duplicate */
736: PetscMalloc3(local_size+1,&ii,local_stash_size,&jj,local_stash_size,&v);
737: ii[0] = 0;
738: cum = 0;
739: for (i=0;i<sub_schurs->n_subs;i++) {
740: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
741: for (j=0;j<subset_size;j++) {
742: const PetscInt row = cum+j;
743: PetscInt col = cum;
745: ii[row+1] = ii[row] + subset_size;
746: for (k=ii[row];k<ii[row+1];k++) {
747: jj[k] = col;
748: col++;
749: }
750: }
751: cum += subset_size;
752: }
753: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,local_size,local_size,ii,jj,v,&T);
754: MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&sub_schurs->S_Ej_all);
755: MatDestroy(&T);
756: PetscFree3(ii,jj,v);
757: }
758: /* matrices for deluxe scaling and adaptive selection */
759: if (compute_Stilda) {
760: if (!sub_schurs->sum_S_Ej_tilda_all) {
761: MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_tilda_all);
762: }
763: if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
764: MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_inv_all);
765: }
766: }
768: /* Compute Schur complements explicitly */
769: F = NULL;
770: if (!sub_schurs->schur_explicit) {
771: /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
772: it is not efficient, unless the economic version of the scaling is used */
773: Mat S_Ej_expl;
774: PetscScalar *work;
775: PetscInt j,*dummy_idx;
776: PetscBool Sdense;
778: PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
779: local_size = 0;
780: for (i=0;i<sub_schurs->n_subs;i++) {
781: IS is_subset_B;
782: Mat AE_EE,AE_IE,AE_EI,S_Ej;
784: /* subsets in original and boundary numbering */
785: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);
786: /* EE block */
787: MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);
788: /* IE block */
789: MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);
790: /* EI block */
791: if (sub_schurs->is_symmetric) {
792: MatCreateTranspose(AE_IE,&AE_EI);
793: } else if (sub_schurs->is_hermitian) {
794: MatCreateHermitianTranspose(AE_IE,&AE_EI);
795: } else {
796: MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);
797: }
798: ISDestroy(&is_subset_B);
799: MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);
800: MatDestroy(&AE_EE);
801: MatDestroy(&AE_IE);
802: MatDestroy(&AE_EI);
803: if (AE_II == A_II) { /* we can reuse the same ksp */
804: KSP ksp;
805: MatSchurComplementGetKSP(sub_schurs->S,&ksp);
806: MatSchurComplementSetKSP(S_Ej,ksp);
807: } else { /* build new ksp object which inherits ksp and pc types from the original one */
808: KSP origksp,schurksp;
809: PC origpc,schurpc;
810: KSPType ksp_type;
811: PetscInt n_internal;
812: PetscBool ispcnone;
814: MatSchurComplementGetKSP(sub_schurs->S,&origksp);
815: MatSchurComplementGetKSP(S_Ej,&schurksp);
816: KSPGetType(origksp,&ksp_type);
817: KSPSetType(schurksp,ksp_type);
818: KSPGetPC(schurksp,&schurpc);
819: KSPGetPC(origksp,&origpc);
820: PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);
821: if (!ispcnone) {
822: PCType pc_type;
823: PCGetType(origpc,&pc_type);
824: PCSetType(schurpc,pc_type);
825: } else {
826: PCSetType(schurpc,PCLU);
827: }
828: ISGetSize(is_I,&n_internal);
829: if (!n_internal) { /* UMFPACK gives error with 0 sized problems */
830: MatSolverType solver = NULL;
831: PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);
832: if (solver) {
833: PCFactorSetMatSolverType(schurpc,solver);
834: }
835: }
836: KSPSetUp(schurksp);
837: }
838: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
839: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);
840: PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);
841: PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);
842: if (Sdense) {
843: for (j=0;j<subset_size;j++) {
844: dummy_idx[j]=local_size+j;
845: }
846: MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
847: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
848: MatDestroy(&S_Ej);
849: MatDestroy(&S_Ej_expl);
850: local_size += subset_size;
851: }
852: PetscFree2(dummy_idx,work);
853: /* free */
854: ISDestroy(&is_I);
855: MatDestroy(&AE_II);
856: PetscFree(all_local_idx_N);
857: } else {
858: Mat A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
859: Vec Dall = NULL;
860: IS is_A_all,*is_p_r = NULL;
861: MatType Stype;
862: PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
863: PetscScalar *SEj_arr = NULL,*SEjinv_arr = NULL;
864: const PetscScalar *rS_data;
865: PetscInt n,n_I,size_schur,size_active_schur,cum,cum2;
866: PetscBool economic,solver_S,S_lower_triangular = PETSC_FALSE;
867: PetscBool schur_has_vertices,factor_workaround;
868: PetscBool use_cholesky;
869: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
870: PetscBool oldpin;
871: #endif
873: /* get sizes */
874: n_I = 0;
875: if (is_I_layer) {
876: ISGetLocalSize(is_I_layer,&n_I);
877: }
878: economic = PETSC_FALSE;
879: ISGetLocalSize(sub_schurs->is_I,&cum);
880: if (cum != n_I) economic = PETSC_TRUE;
881: MatGetLocalSize(sub_schurs->A,&n,NULL);
882: size_active_schur = local_size;
884: /* import scaling vector (wrong formulation if we have 3D edges) */
885: if (scaling && compute_Stilda) {
886: const PetscScalar *array;
887: PetscScalar *array2;
888: const PetscInt *idxs;
889: PetscInt i;
891: ISGetIndices(sub_schurs->is_Ej_all,&idxs);
892: VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);
893: VecGetArrayRead(scaling,&array);
894: VecGetArray(Dall,&array2);
895: for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
896: VecRestoreArray(Dall,&array2);
897: VecRestoreArrayRead(scaling,&array);
898: ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);
899: deluxe = PETSC_FALSE;
900: }
902: /* size active schurs does not count any dirichlet or vertex dof on the interface */
903: factor_workaround = PETSC_FALSE;
904: schur_has_vertices = PETSC_FALSE;
905: cum = n_I+size_active_schur;
906: if (sub_schurs->is_dir) {
907: const PetscInt* idxs;
908: PetscInt n_dir;
910: ISGetLocalSize(sub_schurs->is_dir,&n_dir);
911: ISGetIndices(sub_schurs->is_dir,&idxs);
912: PetscArraycpy(all_local_idx_N+cum,idxs,n_dir);
913: ISRestoreIndices(sub_schurs->is_dir,&idxs);
914: cum += n_dir;
915: factor_workaround = PETSC_TRUE;
916: }
917: /* include the primal vertices in the Schur complement */
918: if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
919: PetscInt n_v;
921: ISGetLocalSize(sub_schurs->is_vertices,&n_v);
922: if (n_v) {
923: const PetscInt* idxs;
925: ISGetIndices(sub_schurs->is_vertices,&idxs);
926: PetscArraycpy(all_local_idx_N+cum,idxs,n_v);
927: ISRestoreIndices(sub_schurs->is_vertices,&idxs);
928: cum += n_v;
929: factor_workaround = PETSC_TRUE;
930: schur_has_vertices = PETSC_TRUE;
931: }
932: }
933: size_schur = cum - n_I;
934: ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);
935: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
936: oldpin = sub_schurs->A->boundtocpu;
937: MatBindToCPU(sub_schurs->A,PETSC_TRUE);
938: #endif
939: if (cum == n) {
940: ISSetPermutation(is_A_all);
941: MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);
942: } else {
943: MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);
944: }
945: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
946: MatBindToCPU(sub_schurs->A,oldpin);
947: #endif
948: MatSetOptionsPrefix(A,sub_schurs->prefix);
949: MatAppendOptionsPrefix(A,"sub_schurs_");
951: /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
952: this is a workaround */
953: if (benign_n) {
954: Vec v,benign_AIIm1_ones;
955: ISLocalToGlobalMapping N_to_reor;
956: IS is_p0,is_p0_p;
957: PetscScalar *cs_AIB,*AIIm1_data;
958: PetscInt sizeA;
960: ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);
961: ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);
962: ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);
963: ISDestroy(&is_p0);
964: MatCreateVecs(A,&v,&benign_AIIm1_ones);
965: VecGetSize(v,&sizeA);
966: MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);
967: MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);
968: MatDenseGetArray(cs_AIB_mat,&cs_AIB);
969: MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
970: PetscMalloc1(benign_n,&is_p_r);
971: /* compute colsum of A_IB restricted to pressures */
972: for (i=0;i<benign_n;i++) {
973: const PetscScalar *array;
974: const PetscInt *idxs;
975: PetscInt j,nz;
977: ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);
978: ISGetLocalSize(is_p_r[i],&nz);
979: ISGetIndices(is_p_r[i],&idxs);
980: for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
981: ISRestoreIndices(is_p_r[i],&idxs);
982: VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);
983: MatMult(A,benign_AIIm1_ones,v);
984: VecResetArray(benign_AIIm1_ones);
985: VecGetArrayRead(v,&array);
986: for (j=0;j<size_schur;j++) {
987: #if defined(PETSC_USE_COMPLEX)
988: cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
989: #else
990: cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
991: #endif
992: }
993: VecRestoreArrayRead(v,&array);
994: }
995: MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
996: MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
997: VecDestroy(&v);
998: VecDestroy(&benign_AIIm1_ones);
999: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);
1000: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
1001: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
1002: MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);
1003: ISDestroy(&is_p0_p);
1004: ISLocalToGlobalMappingDestroy(&N_to_reor);
1005: }
1006: MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);
1007: MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);
1008: MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);
1010: /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
1011: use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
1013: /* when using the benign subspace trick, the local Schur complements are SPD */
1014: if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
1016: if (n_I) {
1017: IS is_schur;
1018: char stype[64];
1019: PetscBool gpu;
1021: if (use_cholesky) {
1022: MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);
1023: } else {
1024: MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);
1025: }
1026: MatSetErrorIfFailure(A,PETSC_TRUE);
1028: /* subsets ordered last */
1029: ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);
1030: MatFactorSetSchurIS(F,is_schur);
1031: ISDestroy(&is_schur);
1033: /* factorization step */
1034: if (use_cholesky) {
1035: MatCholeskyFactorSymbolic(F,A,NULL,NULL);
1036: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1037: MatMumpsSetIcntl(F,19,2);
1038: #endif
1039: MatCholeskyFactorNumeric(F,A,NULL);
1040: S_lower_triangular = PETSC_TRUE;
1041: } else {
1042: MatLUFactorSymbolic(F,A,NULL,NULL,NULL);
1043: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1044: MatMumpsSetIcntl(F,19,3);
1045: #endif
1046: MatLUFactorNumeric(F,A,NULL);
1047: }
1048: MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");
1050: if (matl_dbg_viewer) {
1051: Mat S;
1052: IS is;
1054: PetscObjectSetName((PetscObject)A,"A");
1055: MatView(A,matl_dbg_viewer);
1056: MatFactorCreateSchurComplement(F,&S,NULL);
1057: PetscObjectSetName((PetscObject)S,"S");
1058: MatView(S,matl_dbg_viewer);
1059: MatDestroy(&S);
1060: ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);
1061: PetscObjectSetName((PetscObject)is,"I");
1062: ISView(is,matl_dbg_viewer);
1063: ISDestroy(&is);
1064: ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);
1065: PetscObjectSetName((PetscObject)is,"B");
1066: ISView(is,matl_dbg_viewer);
1067: ISDestroy(&is);
1068: PetscObjectSetName((PetscObject)is_A_all,"IA");
1069: ISView(is_A_all,matl_dbg_viewer);
1070: }
1072: /* get explicit Schur Complement computed during numeric factorization */
1073: MatFactorGetSchurComplement(F,&S_all,NULL);
1074: PetscStrncpy(stype,MATSEQDENSE,sizeof(stype));
1075: PetscObjectTypeCompareAny((PetscObject)A,&gpu,MATSEQAIJVIENNACL,MATSEQAIJCUSPARSE,"");
1076: if (gpu) {
1077: PetscStrncpy(stype,MATSEQDENSECUDA,sizeof(stype));
1078: }
1079: PetscOptionsGetString(NULL,sub_schurs->prefix,"-sub_schurs_schur_mat_type",stype,sizeof(stype),NULL);
1080: MatConvert(S_all,stype,MAT_INPLACE_MATRIX,&S_all);
1081: MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);
1082: MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);
1083: MatGetType(S_all,&Stype);
1085: /* we can reuse the solvers if we are not using the economic version */
1086: reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1087: factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1088: if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
1089: reuse_solvers = factor_workaround = PETSC_FALSE;
1091: solver_S = PETSC_TRUE;
1093: /* update the Schur complement with the change of basis on the pressures */
1094: if (benign_n) {
1095: const PetscScalar *cs_AIB;
1096: PetscScalar *S_data,*AIIm1_data;
1097: Mat S2 = NULL,S3 = NULL; /* dbg */
1098: PetscScalar *S2_data,*S3_data; /* dbg */
1099: Vec v,benign_AIIm1_ones;
1100: PetscInt sizeA;
1102: MatDenseGetArray(S_all,&S_data);
1103: MatCreateVecs(A,&v,&benign_AIIm1_ones);
1104: VecGetSize(v,&sizeA);
1105: #if defined(PETSC_HAVE_MUMPS)
1106: MatMumpsSetIcntl(F,26,0);
1107: #endif
1108: #if defined(PETSC_HAVE_MKL_PARDISO)
1109: MatMkl_PardisoSetCntl(F,70,1);
1110: #endif
1111: MatDenseGetArrayRead(cs_AIB_mat,&cs_AIB);
1112: MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
1113: if (matl_dbg_viewer) {
1114: MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S2);
1115: MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S3);
1116: MatDenseGetArray(S2,&S2_data);
1117: MatDenseGetArray(S3,&S3_data);
1118: }
1119: for (i=0;i<benign_n;i++) {
1120: PetscScalar *array,sum = 0.,one = 1.,*sums;
1121: const PetscInt *idxs;
1122: PetscInt k,j,nz;
1123: PetscBLASInt B_k,B_n;
1125: PetscCalloc1(benign_n,&sums);
1126: VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);
1127: VecCopy(benign_AIIm1_ones,v);
1128: MatSolve(F,v,benign_AIIm1_ones);
1129: MatMult(A,benign_AIIm1_ones,v);
1130: VecResetArray(benign_AIIm1_ones);
1131: /* p0 dofs (eliminated) are excluded from the sums */
1132: for (k=0;k<benign_n;k++) {
1133: ISGetLocalSize(is_p_r[k],&nz);
1134: ISGetIndices(is_p_r[k],&idxs);
1135: for (j=0;j<nz-1;j++) sums[k] -= AIIm1_data[idxs[j]+sizeA*i];
1136: ISRestoreIndices(is_p_r[k],&idxs);
1137: }
1138: VecGetArrayRead(v,(const PetscScalar**)&array);
1139: if (matl_dbg_viewer) {
1140: Vec vv;
1141: char name[16];
1143: VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array+n_I,&vv);
1144: PetscSNPrintf(name,sizeof(name),"Pvs%D",i);
1145: PetscObjectSetName((PetscObject)vv,name);
1146: VecView(vv,matl_dbg_viewer);
1147: }
1148: /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1149: /* cs_AIB already scaled by 1./nz */
1150: B_k = 1;
1151: for (k=0;k<benign_n;k++) {
1152: sum = sums[k];
1153: PetscBLASIntCast(size_schur,&B_n);
1155: if (PetscAbsScalar(sum) == 0.0) continue;
1156: if (k == i) {
1157: PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1158: if (matl_dbg_viewer) {
1159: PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n));
1160: }
1161: } else { /* XXX Is it correct to use symmetric rank-2 update with half of the sum? */
1162: sum /= 2.0;
1163: PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,cs_AIB+k*size_schur,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1164: if (matl_dbg_viewer) {
1165: PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,cs_AIB+k*size_schur,&B_n,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n));
1166: }
1167: }
1168: }
1169: sum = 1.;
1170: 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));
1171: if (matl_dbg_viewer) {
1172: PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S2_data,&B_n));
1173: }
1174: VecRestoreArrayRead(v,(const PetscScalar**)&array);
1175: /* set p0 entry of AIIm1_ones to zero */
1176: ISGetLocalSize(is_p_r[i],&nz);
1177: ISGetIndices(is_p_r[i],&idxs);
1178: for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1179: ISRestoreIndices(is_p_r[i],&idxs);
1180: PetscFree(sums);
1181: }
1182: VecDestroy(&benign_AIIm1_ones);
1183: if (matl_dbg_viewer) {
1184: MatDenseRestoreArray(S2,&S2_data);
1185: MatDenseRestoreArray(S3,&S3_data);
1186: }
1187: if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1188: PetscInt k,j;
1189: for (k=0;k<size_schur;k++) {
1190: for (j=k;j<size_schur;j++) {
1191: S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1192: }
1193: }
1194: }
1196: /* restore defaults */
1197: #if defined(PETSC_HAVE_MUMPS)
1198: MatMumpsSetIcntl(F,26,-1);
1199: #endif
1200: #if defined(PETSC_HAVE_MKL_PARDISO)
1201: MatMkl_PardisoSetCntl(F,70,0);
1202: #endif
1203: MatDenseRestoreArrayRead(cs_AIB_mat,&cs_AIB);
1204: MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
1205: VecDestroy(&v);
1206: MatDenseRestoreArray(S_all,&S_data);
1207: if (matl_dbg_viewer) {
1208: Mat S;
1210: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1211: MatFactorCreateSchurComplement(F,&S,NULL);
1212: PetscObjectSetName((PetscObject)S,"Sb");
1213: MatView(S,matl_dbg_viewer);
1214: MatDestroy(&S);
1215: PetscObjectSetName((PetscObject)S2,"S2P");
1216: MatView(S2,matl_dbg_viewer);
1217: PetscObjectSetName((PetscObject)S3,"S3P");
1218: MatView(S3,matl_dbg_viewer);
1219: PetscObjectSetName((PetscObject)cs_AIB_mat,"cs");
1220: MatView(cs_AIB_mat,matl_dbg_viewer);
1221: MatFactorGetSchurComplement(F,&S_all,NULL);
1222: }
1223: MatDestroy(&S2);
1224: MatDestroy(&S3);
1225: }
1226: if (!reuse_solvers) {
1227: for (i=0;i<benign_n;i++) {
1228: ISDestroy(&is_p_r[i]);
1229: }
1230: PetscFree(is_p_r);
1231: MatDestroy(&cs_AIB_mat);
1232: MatDestroy(&benign_AIIm1_ones_mat);
1233: }
1234: } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1235: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);
1236: MatGetType(S_all,&Stype);
1237: reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1238: factor_workaround = PETSC_FALSE;
1239: solver_S = PETSC_FALSE;
1240: }
1242: if (reuse_solvers) {
1243: Mat A_II,Afake;
1244: Vec vec1_B;
1245: PCBDDCReuseSolvers msolv_ctx;
1246: PetscInt n_R;
1248: if (sub_schurs->reuse_solver) {
1249: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1250: } else {
1251: PetscNew(&sub_schurs->reuse_solver);
1252: }
1253: msolv_ctx = sub_schurs->reuse_solver;
1254: MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);
1255: PetscObjectReference((PetscObject)F);
1256: msolv_ctx->F = F;
1257: MatCreateVecs(F,&msolv_ctx->sol,NULL);
1258: /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1259: {
1260: PetscScalar *array;
1261: PetscInt n;
1263: VecGetLocalSize(msolv_ctx->sol,&n);
1264: VecGetArray(msolv_ctx->sol,&array);
1265: VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);
1266: VecRestoreArray(msolv_ctx->sol,&array);
1267: }
1268: msolv_ctx->has_vertices = schur_has_vertices;
1270: /* interior solver */
1271: PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);
1272: PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);
1273: PCSetType(msolv_ctx->interior_solver,PCSHELL);
1274: PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");
1275: PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);
1276: PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);
1277: PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);
1278: PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);
1280: /* correction solver */
1281: PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);
1282: PCSetType(msolv_ctx->correction_solver,PCSHELL);
1283: PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");
1284: PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);
1285: PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);
1286: PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);
1287: PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);
1289: /* scatter and vecs for Schur complement solver */
1290: MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);
1291: MatCreateVecs(sub_schurs->S,&vec1_B,NULL);
1292: if (!schur_has_vertices) {
1293: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);
1294: VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);
1295: PetscObjectReference((PetscObject)is_A_all);
1296: msolv_ctx->is_R = is_A_all;
1297: } else {
1298: IS is_B_all;
1299: const PetscInt* idxs;
1300: PetscInt dual,n_v,n;
1302: ISGetLocalSize(sub_schurs->is_vertices,&n_v);
1303: dual = size_schur - n_v;
1304: ISGetLocalSize(is_A_all,&n);
1305: ISGetIndices(is_A_all,&idxs);
1306: ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);
1307: ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);
1308: ISDestroy(&is_B_all);
1309: ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);
1310: VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);
1311: ISDestroy(&is_B_all);
1312: ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);
1313: ISRestoreIndices(is_A_all,&idxs);
1314: }
1315: ISGetLocalSize(msolv_ctx->is_R,&n_R);
1316: MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);
1317: MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);
1318: MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);
1319: PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);
1320: MatDestroy(&Afake);
1321: VecDestroy(&vec1_B);
1323: /* communicate benign info to solver context */
1324: if (benign_n) {
1325: PetscScalar *array;
1327: msolv_ctx->benign_n = benign_n;
1328: msolv_ctx->benign_zerodiag_subs = is_p_r;
1329: PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);
1330: msolv_ctx->benign_csAIB = cs_AIB_mat;
1331: MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);
1332: VecGetArray(msolv_ctx->benign_corr_work,&array);
1333: VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);
1334: VecRestoreArray(msolv_ctx->benign_corr_work,&array);
1335: msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1336: }
1337: } else {
1338: if (sub_schurs->reuse_solver) {
1339: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1340: }
1341: PetscFree(sub_schurs->reuse_solver);
1342: }
1343: MatDestroy(&A);
1344: ISDestroy(&is_A_all);
1346: /* Work arrays */
1347: PetscMalloc1(max_subset_size*max_subset_size,&work);
1349: /* S_Ej_all */
1350: cum = cum2 = 0;
1351: MatDenseGetArrayRead(S_all,&rS_data);
1352: MatSeqAIJGetArray(sub_schurs->S_Ej_all,&SEj_arr);
1353: if (compute_Stilda) {
1354: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&SEjinv_arr);
1355: }
1356: for (i=0;i<sub_schurs->n_subs;i++) {
1357: PetscInt j;
1359: /* get S_E */
1360: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1361: if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1362: PetscInt k;
1363: for (k=0;k<subset_size;k++) {
1364: for (j=k;j<subset_size;j++) {
1365: work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1366: work[j*subset_size+k] = PetscConj(rS_data[cum2+k*size_schur+j]);
1367: }
1368: }
1369: } else { /* just copy to workspace */
1370: PetscInt k;
1371: for (k=0;k<subset_size;k++) {
1372: for (j=0;j<subset_size;j++) {
1373: work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1374: }
1375: }
1376: }
1377: /* insert S_E values */
1378: if (sub_schurs->change) {
1379: Mat change_sub,SEj,T;
1381: /* change basis */
1382: KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1383: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1384: if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1385: Mat T2;
1386: MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1387: MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1388: MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1389: MatDestroy(&T2);
1390: } else {
1391: MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1392: }
1393: MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1394: MatDestroy(&T);
1395: MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);
1396: MatDestroy(&SEj);
1397: }
1398: if (deluxe) {
1399: PetscArraycpy(SEj_arr,work,subset_size*subset_size);
1400: /* if adaptivity is requested, invert S_E blocks */
1401: if (compute_Stilda) {
1402: Mat M;
1403: const PetscScalar *vals;
1404: PetscBool isdense,isdensecuda;
1406: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&M);
1407: MatSetOption(M,MAT_SPD,sub_schurs->is_posdef);
1408: MatSetOption(M,MAT_HERMITIAN,sub_schurs->is_hermitian);
1409: if (!PetscBTLookup(sub_schurs->is_edge,i)) {
1410: MatSetType(M,Stype);
1411: }
1412: PetscObjectTypeCompare((PetscObject)M,MATSEQDENSE,&isdense);
1413: PetscObjectTypeCompare((PetscObject)M,MATSEQDENSECUDA,&isdensecuda);
1414: if (use_cholesky) {
1415: MatCholeskyFactor(M,NULL,NULL);
1416: } else {
1417: MatLUFactor(M,NULL,NULL,NULL);
1418: }
1419: if (isdense) {
1420: MatSeqDenseInvertFactors_Private(M);
1421: #if defined(PETSC_HAVE_CUDA)
1422: } else if (isdensecuda) {
1423: MatSeqDenseCUDAInvertFactors_Private(M);
1424: #endif
1425: } else SETERRQ1(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"Not implemented for type %s",Stype);
1426: MatDenseGetArrayRead(M,&vals);
1427: PetscArraycpy(SEjinv_arr,vals,subset_size*subset_size);
1428: MatDenseRestoreArrayRead(M,&vals);
1429: MatDestroy(&M);
1430: }
1431: } else if (compute_Stilda) { /* not using deluxe */
1432: Mat SEj;
1433: Vec D;
1434: PetscScalar *array;
1436: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1437: VecGetArray(Dall,&array);
1438: VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);
1439: VecRestoreArray(Dall,&array);
1440: VecShift(D,-1.);
1441: MatDiagonalScale(SEj,D,D);
1442: MatDestroy(&SEj);
1443: VecDestroy(&D);
1444: PetscArraycpy(SEj_arr,work,subset_size*subset_size);
1445: }
1446: cum += subset_size;
1447: cum2 += subset_size*(size_schur + 1);
1448: SEj_arr += subset_size*subset_size;
1449: if (SEjinv_arr) SEjinv_arr += subset_size*subset_size;
1450: }
1451: MatDenseRestoreArrayRead(S_all,&rS_data);
1452: MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&SEj_arr);
1453: if (compute_Stilda) {
1454: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&SEjinv_arr);
1455: }
1456: if (solver_S) {
1457: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1458: }
1460: /* may prevent from unneeded copies, since MUMPS or MKL_Pardiso always use CPU memory
1461: however, preliminary tests indicate using GPUs is still faster in the solve phase */
1462: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1463: if (reuse_solvers) {
1464: Mat St;
1465: MatFactorSchurStatus st;
1466: PetscBool flg = PETSC_FALSE;
1468: PetscOptionsGetBool(NULL,sub_schurs->prefix,"-sub_schurs_schur_pin_to_cpu",&flg,NULL);
1469: MatFactorGetSchurComplement(F,&St,&st);
1470: MatBindToCPU(St,flg);
1471: MatFactorRestoreSchurComplement(F,&St,st);
1472: }
1473: #endif
1475: schur_factor = NULL;
1476: if (compute_Stilda && size_active_schur) {
1478: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&SEjinv_arr);
1479: if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1480: PetscArraycpy(SEjinv_arr,work,size_schur*size_schur);
1481: } else {
1482: Mat S_all_inv=NULL;
1484: if (solver_S) {
1485: /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1486: 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 */
1487: if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1488: PetscScalar *data;
1489: PetscInt nd = 0;
1491: if (!use_potr) {
1492: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1493: }
1494: MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1495: MatDenseGetArray(S_all_inv,&data);
1496: if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1497: ISGetLocalSize(sub_schurs->is_dir,&nd);
1498: }
1500: /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1501: if (schur_has_vertices) {
1502: Mat M;
1503: PetscScalar *tdata;
1504: PetscInt nv = 0, news;
1506: ISGetLocalSize(sub_schurs->is_vertices,&nv);
1507: news = size_active_schur + nv;
1508: PetscCalloc1(news*news,&tdata);
1509: for (i=0;i<size_active_schur;i++) {
1510: PetscArraycpy(tdata+i*(news+1),data+i*(size_schur+1),size_active_schur-i);
1511: PetscArraycpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv);
1512: }
1513: for (i=0;i<nv;i++) {
1514: PetscInt k = i+size_active_schur;
1515: PetscArraycpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),nv-i);
1516: }
1518: MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);
1519: MatSetOption(M,MAT_SPD,PETSC_TRUE);
1520: MatCholeskyFactor(M,NULL,NULL);
1521: /* save the factors */
1522: cum = 0;
1523: PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);
1524: for (i=0;i<size_active_schur;i++) {
1525: PetscArraycpy(schur_factor+cum,tdata+i*(news+1),size_active_schur-i);
1526: cum += size_active_schur - i;
1527: }
1528: for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1529: MatSeqDenseInvertFactors_Private(M);
1530: /* move back just the active dofs to the Schur complement */
1531: for (i=0;i<size_active_schur;i++) {
1532: PetscArraycpy(data+i*size_schur,tdata+i*news,size_active_schur);
1533: }
1534: PetscFree(tdata);
1535: MatDestroy(&M);
1536: } else { /* we can factorize and invert just the activedofs */
1537: Mat M;
1538: PetscScalar *aux;
1540: PetscMalloc1(nd,&aux);
1541: for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1542: MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);
1543: MatSeqDenseSetLDA(M,size_schur);
1544: MatSetOption(M,MAT_SPD,PETSC_TRUE);
1545: MatCholeskyFactor(M,NULL,NULL);
1546: MatSeqDenseInvertFactors_Private(M);
1547: MatDestroy(&M);
1548: MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);
1549: MatZeroEntries(M);
1550: MatDestroy(&M);
1551: MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);
1552: MatSeqDenseSetLDA(M,size_schur);
1553: MatZeroEntries(M);
1554: MatDestroy(&M);
1555: for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
1556: PetscFree(aux);
1557: }
1558: MatDenseRestoreArray(S_all_inv,&data);
1559: } else { /* use MatFactor calls to invert S */
1560: MatFactorInvertSchurComplement(F);
1561: MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1562: }
1563: } else { /* we need to invert explicitly since we are not using MatFactor for S */
1564: PetscObjectReference((PetscObject)S_all);
1565: S_all_inv = S_all;
1566: MatDenseGetArray(S_all_inv,&S_data);
1567: PetscBLASIntCast(size_schur,&B_N);
1568: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1569: if (use_potr) {
1570: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1571: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1572: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1573: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1574: } else if (use_sytr) {
1575: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1576: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1577: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1578: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1579: } else {
1580: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1581: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1582: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1583: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1584: }
1585: PetscLogFlops(1.0*size_schur*size_schur*size_schur);
1586: PetscFPTrapPop();
1587: MatDenseRestoreArray(S_all_inv,&S_data);
1588: }
1589: /* S_Ej_tilda_all */
1590: cum = cum2 = 0;
1591: MatDenseGetArrayRead(S_all_inv,&rS_data);
1592: for (i=0;i<sub_schurs->n_subs;i++) {
1593: PetscInt j;
1595: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1596: /* get (St^-1)_E */
1597: /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1598: will be properly accessed later during adaptive selection */
1599: if (S_lower_triangular) {
1600: PetscInt k;
1601: if (sub_schurs->change) {
1602: for (k=0;k<subset_size;k++) {
1603: for (j=k;j<subset_size;j++) {
1604: work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1605: work[j*subset_size+k] = work[k*subset_size+j];
1606: }
1607: }
1608: } else {
1609: for (k=0;k<subset_size;k++) {
1610: for (j=k;j<subset_size;j++) {
1611: work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1612: }
1613: }
1614: }
1615: } else {
1616: PetscInt k;
1617: for (k=0;k<subset_size;k++) {
1618: for (j=0;j<subset_size;j++) {
1619: work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1620: }
1621: }
1622: }
1623: if (sub_schurs->change) {
1624: Mat change_sub,SEj,T;
1626: /* change basis */
1627: KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1628: MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1629: if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1630: Mat T2;
1631: MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1632: MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1633: MatDestroy(&T2);
1634: MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1635: } else {
1636: MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1637: }
1638: MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1639: MatDestroy(&T);
1640: /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1641: MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);
1642: MatDestroy(&SEj);
1643: }
1644: PetscArraycpy(SEjinv_arr,work,subset_size*subset_size);
1645: cum += subset_size;
1646: cum2 += subset_size*(size_schur + 1);
1647: SEjinv_arr += subset_size*subset_size;
1648: }
1649: MatDenseRestoreArrayRead(S_all_inv,&rS_data);
1650: if (solver_S) {
1651: if (schur_has_vertices) {
1652: MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);
1653: } else {
1654: MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);
1655: }
1656: }
1657: MatDestroy(&S_all_inv);
1658: }
1659: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&SEjinv_arr);
1661: /* move back factors if needed */
1662: if (schur_has_vertices) {
1663: Mat S_tmp;
1664: PetscInt nd = 0;
1666: if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1667: MatFactorGetSchurComplement(F,&S_tmp,NULL);
1668: if (use_potr) {
1669: PetscScalar *data;
1671: MatDenseGetArray(S_tmp,&data);
1672: PetscArrayzero(data,size_schur*size_schur);
1674: if (S_lower_triangular) {
1675: cum = 0;
1676: for (i=0;i<size_active_schur;i++) {
1677: PetscArraycpy(data+i*(size_schur+1),schur_factor+cum,size_active_schur-i);
1678: cum += size_active_schur-i;
1679: }
1680: } else {
1681: PetscArraycpy(data,schur_factor,size_schur*size_schur);
1682: }
1683: if (sub_schurs->is_dir) {
1684: ISGetLocalSize(sub_schurs->is_dir,&nd);
1685: for (i=0;i<nd;i++) {
1686: data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1687: }
1688: }
1689: /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1690: set the diagonal entry of the Schur factor to a very large value */
1691: for (i=size_active_schur+nd;i<size_schur;i++) {
1692: data[i*(size_schur+1)] = infty;
1693: }
1694: MatDenseRestoreArray(S_tmp,&data);
1695: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1696: MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);
1697: }
1698: } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1699: PetscScalar *data;
1700: PetscInt nd = 0;
1702: if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1703: ISGetLocalSize(sub_schurs->is_dir,&nd);
1704: }
1705: MatFactorGetSchurComplement(F,&S_all,NULL);
1706: MatDenseGetArray(S_all,&data);
1707: for (i=0;i<size_active_schur;i++) {
1708: PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);
1709: }
1710: for (i=size_active_schur+nd;i<size_schur;i++) {
1711: PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);
1712: data[i*(size_schur+1)] = infty;
1713: }
1714: MatDenseRestoreArray(S_all,&data);
1715: MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1716: }
1717: PetscFree(work);
1718: PetscFree(schur_factor);
1719: VecDestroy(&Dall);
1720: }
1721: ISDestroy(&is_I_layer);
1722: MatDestroy(&S_all);
1723: MatDestroy(&A_BB);
1724: MatDestroy(&A_IB);
1725: MatDestroy(&A_BI);
1726: MatDestroy(&F);
1728: PetscMalloc1(sub_schurs->n_subs,&nnz);
1729: for (i=0;i<sub_schurs->n_subs;i++) {
1730: ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);
1731: }
1732: ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);
1733: MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);
1734: MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1735: MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1736: if (compute_Stilda) {
1737: MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);
1738: MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1739: MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1740: if (deluxe) {
1741: MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);
1742: MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1743: MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1744: }
1745: }
1746: ISDestroy(&is_I_layer);
1748: /* Get local part of (\sum_j S_Ej) */
1749: if (!sub_schurs->sum_S_Ej_all) {
1750: MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);
1751: }
1752: VecSet(gstash,0.0);
1753: MatSeqAIJGetArray(sub_schurs->S_Ej_all,&stasharray);
1754: VecPlaceArray(lstash,stasharray);
1755: VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1756: VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1757: MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&stasharray);
1758: VecResetArray(lstash);
1759: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&stasharray);
1760: VecPlaceArray(lstash,stasharray);
1761: VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1762: VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1763: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&stasharray);
1764: VecResetArray(lstash);
1766: /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1767: if (compute_Stilda) {
1768: VecSet(gstash,0.0);
1769: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);
1770: VecPlaceArray(lstash,stasharray);
1771: VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1772: VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1773: VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1774: VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1775: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);
1776: VecResetArray(lstash);
1777: if (deluxe) {
1778: VecSet(gstash,0.0);
1779: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);
1780: VecPlaceArray(lstash,stasharray);
1781: VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1782: VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1783: VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1784: VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1785: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);
1786: VecResetArray(lstash);
1787: } else {
1788: PetscScalar *array;
1789: PetscInt cum;
1791: MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1792: cum = 0;
1793: for (i=0;i<sub_schurs->n_subs;i++) {
1794: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1795: PetscBLASIntCast(subset_size,&B_N);
1796: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1797: if (use_potr) {
1798: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1799: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1800: PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1801: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1802: } else if (use_sytr) {
1803: PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1804: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1805: PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1806: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1807: } else {
1808: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1809: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1810: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1811: if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1812: }
1813: PetscLogFlops(1.0*subset_size*subset_size*subset_size);
1814: PetscFPTrapPop();
1815: cum += subset_size*subset_size;
1816: }
1817: MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1818: PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);
1819: MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1820: sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1821: }
1822: }
1823: VecDestroy(&lstash);
1824: VecDestroy(&gstash);
1825: VecScatterDestroy(&sstash);
1827: if (matl_dbg_viewer) {
1828: PetscInt cum;
1830: if (sub_schurs->S_Ej_all) {
1831: PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");
1832: MatView(sub_schurs->S_Ej_all,matl_dbg_viewer);
1833: }
1834: if (sub_schurs->sum_S_Ej_all) {
1835: PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");
1836: MatView(sub_schurs->sum_S_Ej_all,matl_dbg_viewer);
1837: }
1838: if (sub_schurs->sum_S_Ej_inv_all) {
1839: PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");
1840: MatView(sub_schurs->sum_S_Ej_inv_all,matl_dbg_viewer);
1841: }
1842: if (sub_schurs->sum_S_Ej_tilda_all) {
1843: PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");
1844: MatView(sub_schurs->sum_S_Ej_tilda_all,matl_dbg_viewer);
1845: }
1846: for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
1847: IS is;
1848: char name[16];
1850: PetscSNPrintf(name,sizeof(name),"IE%D",i);
1851: ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1852: ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);
1853: PetscObjectSetName((PetscObject)is,name);
1854: ISView(is,matl_dbg_viewer);
1855: ISDestroy(&is);
1856: cum += subset_size;
1857: }
1858: }
1860: /* free workspace */
1861: PetscViewerDestroy(&matl_dbg_viewer);
1862: PetscFree2(Bwork,pivots);
1863: PetscCommDestroy(&comm_n);
1864: return(0);
1865: }
1867: PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1868: {
1869: IS *faces,*edges,*all_cc,vertices;
1870: PetscInt i,n_faces,n_edges,n_all_cc;
1871: PetscBool is_sorted,ispardiso,ismumps;
1872: PetscErrorCode ierr;
1875: ISSorted(is_I,&is_sorted);
1876: if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1877: ISSorted(is_B,&is_sorted);
1878: if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1880: /* reset any previous data */
1881: PCBDDCSubSchursReset(sub_schurs);
1883: /* get index sets for faces and edges (already sorted by global ordering) */
1884: PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1885: n_all_cc = n_faces+n_edges;
1886: PetscBTCreate(n_all_cc,&sub_schurs->is_edge);
1887: PetscMalloc1(n_all_cc,&all_cc);
1888: for (i=0;i<n_faces;i++) {
1889: if (copycc) {
1890: ISDuplicate(faces[i],&all_cc[i]);
1891: } else {
1892: PetscObjectReference((PetscObject)faces[i]);
1893: all_cc[i] = faces[i];
1894: }
1895: }
1896: for (i=0;i<n_edges;i++) {
1897: if (copycc) {
1898: ISDuplicate(edges[i],&all_cc[n_faces+i]);
1899: } else {
1900: PetscObjectReference((PetscObject)edges[i]);
1901: all_cc[n_faces+i] = edges[i];
1902: }
1903: PetscBTSet(sub_schurs->is_edge,n_faces+i);
1904: }
1905: PetscObjectReference((PetscObject)vertices);
1906: sub_schurs->is_vertices = vertices;
1907: PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1908: sub_schurs->is_dir = NULL;
1909: PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);
1911: /* Determine if MatFactor can be used */
1912: PetscStrallocpy(prefix,&sub_schurs->prefix);
1913: #if defined(PETSC_HAVE_MUMPS)
1914: PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);
1915: #elif defined(PETSC_HAVE_MKL_PARDISO)
1916: PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);
1917: #else
1918: PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);
1919: #endif
1920: #if defined(PETSC_USE_COMPLEX)
1921: sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1922: #else
1923: sub_schurs->is_hermitian = PETSC_TRUE;
1924: #endif
1925: sub_schurs->is_posdef = PETSC_TRUE;
1926: sub_schurs->is_symmetric = PETSC_TRUE;
1927: sub_schurs->debug = PETSC_FALSE;
1928: sub_schurs->restrict_comm = PETSC_FALSE;
1929: PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
1930: PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,64,NULL);
1931: PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);
1932: PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);
1933: PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);
1934: PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);
1935: PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);
1936: PetscOptionsEnd();
1937: PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMUMPS,&ismumps);
1938: PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,&ispardiso);
1939: sub_schurs->schur_explicit = (PetscBool)(ispardiso || ismumps);
1941: /* for reals, symmetric and hermitian are synonims */
1942: #if !defined(PETSC_USE_COMPLEX)
1943: sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1944: sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1945: #endif
1947: PetscObjectReference((PetscObject)is_I);
1948: sub_schurs->is_I = is_I;
1949: PetscObjectReference((PetscObject)is_B);
1950: sub_schurs->is_B = is_B;
1951: PetscObjectReference((PetscObject)graph->l2gmap);
1952: sub_schurs->l2gmap = graph->l2gmap;
1953: PetscObjectReference((PetscObject)BtoNmap);
1954: sub_schurs->BtoNmap = BtoNmap;
1955: sub_schurs->n_subs = n_all_cc;
1956: sub_schurs->is_subs = all_cc;
1957: sub_schurs->S_Ej_all = NULL;
1958: sub_schurs->sum_S_Ej_all = NULL;
1959: sub_schurs->sum_S_Ej_inv_all = NULL;
1960: sub_schurs->sum_S_Ej_tilda_all = NULL;
1961: sub_schurs->is_Ej_all = NULL;
1962: return(0);
1963: }
1965: PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1966: {
1967: PCBDDCSubSchurs schurs_ctx;
1968: PetscErrorCode ierr;
1971: PetscNew(&schurs_ctx);
1972: schurs_ctx->n_subs = 0;
1973: *sub_schurs = schurs_ctx;
1974: return(0);
1975: }
1977: PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1978: {
1979: PetscInt i;
1983: if (!sub_schurs) return(0);
1984: PetscFree(sub_schurs->prefix);
1985: MatDestroy(&sub_schurs->A);
1986: MatDestroy(&sub_schurs->S);
1987: ISDestroy(&sub_schurs->is_I);
1988: ISDestroy(&sub_schurs->is_B);
1989: ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);
1990: ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);
1991: MatDestroy(&sub_schurs->S_Ej_all);
1992: MatDestroy(&sub_schurs->sum_S_Ej_all);
1993: MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1994: MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);
1995: ISDestroy(&sub_schurs->is_Ej_all);
1996: ISDestroy(&sub_schurs->is_vertices);
1997: ISDestroy(&sub_schurs->is_dir);
1998: PetscBTDestroy(&sub_schurs->is_edge);
1999: for (i=0;i<sub_schurs->n_subs;i++) {
2000: ISDestroy(&sub_schurs->is_subs[i]);
2001: }
2002: if (sub_schurs->n_subs) {
2003: PetscFree(sub_schurs->is_subs);
2004: }
2005: if (sub_schurs->reuse_solver) {
2006: PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
2007: }
2008: PetscFree(sub_schurs->reuse_solver);
2009: if (sub_schurs->change) {
2010: for (i=0;i<sub_schurs->n_subs;i++) {
2011: KSPDestroy(&sub_schurs->change[i]);
2012: ISDestroy(&sub_schurs->change_primal_sub[i]);
2013: }
2014: }
2015: PetscFree(sub_schurs->change);
2016: PetscFree(sub_schurs->change_primal_sub);
2017: sub_schurs->n_subs = 0;
2018: return(0);
2019: }
2021: PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
2022: {
2026: PCBDDCSubSchursReset(*sub_schurs);
2027: PetscFree(*sub_schurs);
2028: return(0);
2029: }
2031: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
2032: {
2033: PetscInt i,j,n;
2037: n = 0;
2038: for (i=-n_prev;i<0;i++) {
2039: PetscInt start_dof = queue_tip[i];
2040: for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
2041: PetscInt dof = adjncy[j];
2042: if (!PetscBTLookup(touched,dof)) {
2043: PetscBTSet(touched,dof);
2044: queue_tip[n] = dof;
2045: n++;
2046: }
2047: }
2048: }
2049: *n_added = n;
2050: return(0);
2051: }