Actual source code: bddcschurs.c

  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,&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,&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;
422:   PetscBool              flg;

425:   MatDestroy(&sub_schurs->A);
426:   MatDestroy(&sub_schurs->S);
427:   if (Ain) {
428:     PetscObjectReference((PetscObject)Ain);
429:     sub_schurs->A = Ain;
430:   }

432:   PetscObjectReference((PetscObject)Sin);
433:   sub_schurs->S = Sin;
434:   if (sub_schurs->schur_explicit) {
435:     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
436:   }

438:   /* preliminary checks */
439:   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");

441:   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;

443:   /* debug (MATLAB) */
444:   if (sub_schurs->debug) {
445:     PetscMPIInt size,rank;
446:     PetscInt    nr,*print_schurs_ranks,print_schurs = PETSC_FALSE;

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:   }

470:   /* restrict work on active processes */
471:   if (sub_schurs->restrict_comm) {
472:     PetscSubcomm subcomm;
473:     PetscMPIInt  color,rank;

475:     color = 0;
476:     if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
477:     MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
478:     PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);
479:     PetscSubcommSetNumber(subcomm,2);
480:     PetscSubcommSetTypeGeneral(subcomm,color,rank);
481:     PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);
482:     PetscSubcommDestroy(&subcomm);
483:     if (!sub_schurs->n_subs) {
484:       PetscCommDestroy(&comm_n);
485:       return(0);
486:     }
487:   } else {
488:     PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);
489:   }

491:   /* get Schur complement matrices */
492:   if (!sub_schurs->schur_explicit) {
493:     Mat       tA_IB,tA_BI,tA_BB;
494:     PetscBool isseqsbaij;
495:     MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);
496:     PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);
497:     if (isseqsbaij) {
498:       MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);
499:       MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);
500:       MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);
501:     } else {
502:       PetscObjectReference((PetscObject)tA_BB);
503:       A_BB = tA_BB;
504:       PetscObjectReference((PetscObject)tA_IB);
505:       A_IB = tA_IB;
506:       PetscObjectReference((PetscObject)tA_BI);
507:       A_BI = tA_BI;
508:     }
509:   } else {
510:     A_II = NULL;
511:     A_IB = NULL;
512:     A_BI = NULL;
513:     A_BB = NULL;
514:   }
515:   S_all = NULL;

517:   /* determine interior problems */
518:   ISGetLocalSize(sub_schurs->is_I,&i);
519:   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
520:     PetscBT                touched;
521:     const PetscInt*        idx_B;
522:     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;

524:     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
525:     /* get sizes */
526:     ISGetLocalSize(sub_schurs->is_I,&n_I);
527:     ISGetLocalSize(sub_schurs->is_B,&n_B);

529:     PetscMalloc1(n_I+n_B,&local_numbering);
530:     PetscBTCreate(n_I+n_B,&touched);
531:     PetscBTMemzero(n_I+n_B,touched);

533:     /* all boundary dofs must be skipped when adding layers */
534:     ISGetIndices(sub_schurs->is_B,&idx_B);
535:     for (j=0;j<n_B;j++) {
536:       PetscBTSet(touched,idx_B[j]);
537:     }
538:     PetscArraycpy(local_numbering,idx_B,n_B);
539:     ISRestoreIndices(sub_schurs->is_B,&idx_B);

541:     /* add prescribed number of layers of dofs */
542:     n_local_dofs = n_B;
543:     n_prev_added = n_B;
544:     for (layer=0;layer<nlayers;layer++) {
545:       PetscInt n_added;
546:       if (n_local_dofs == n_I+n_B) break;
547:       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);
548:       PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);
549:       n_prev_added = n_added;
550:       n_local_dofs += n_added;
551:       if (!n_added) break;
552:     }
553:     PetscBTDestroy(&touched);

555:     /* IS for I layer dofs in original numbering */
556:     ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);
557:     PetscFree(local_numbering);
558:     ISSort(is_I_layer);
559:     /* IS for I layer dofs in I numbering */
560:     if (!sub_schurs->schur_explicit) {
561:       ISLocalToGlobalMapping ItoNmap;
562:       ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);
563:       ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);
564:       ISLocalToGlobalMappingDestroy(&ItoNmap);

566:       /* II block */
567:       MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);
568:     }
569:   } else {
570:     PetscInt n_I;

572:     /* IS for I dofs in original numbering */
573:     PetscObjectReference((PetscObject)sub_schurs->is_I);
574:     is_I_layer = sub_schurs->is_I;

576:     /* IS for I dofs in I numbering (strided 1) */
577:     if (!sub_schurs->schur_explicit) {
578:       ISGetSize(sub_schurs->is_I,&n_I);
579:       ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);

581:       /* II block is the same */
582:       PetscObjectReference((PetscObject)A_II);
583:       AE_II = A_II;
584:     }
585:   }

587:   /* Get info on subset sizes and sum of all subsets sizes */
588:   max_subset_size = 0;
589:   local_size = 0;
590:   for (i=0;i<sub_schurs->n_subs;i++) {
591:     ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
592:     max_subset_size = PetscMax(subset_size,max_subset_size);
593:     local_size += subset_size;
594:   }

596:   /* Work arrays for local indices */
597:   extra = 0;
598:   ISGetLocalSize(sub_schurs->is_B,&n_B);
599:   if (sub_schurs->schur_explicit && is_I_layer) {
600:     ISGetLocalSize(is_I_layer,&extra);
601:   }
602:   PetscMalloc1(n_B+extra,&all_local_idx_N);
603:   if (extra) {
604:     const PetscInt *idxs;
605:     ISGetIndices(is_I_layer,&idxs);
606:     PetscArraycpy(all_local_idx_N,idxs,extra);
607:     ISRestoreIndices(is_I_layer,&idxs);
608:   }
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:   local_stash_size = 0;
615:   for (i=0;i<sub_schurs->n_subs;i++) {
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*subset_size;
623:     /* subset indices in local numbering */
624:     PetscArraycpy(all_local_idx_N+local_size+extra,idxs,subset_size);
625:     ISRestoreIndices(sub_schurs->is_subs[i],&idxs);
626:     local_size += subset_size;
627:     local_stash_size += subset_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 data 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_stash_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_stash_size);
669:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,local_stash_size,NULL,&lstash);
670:   VecCreateMPI(comm_n,PETSC_DECIDE,global_size,&gstash);
671:   VecScatterCreate(lstash,NULL,gstash,all_subsets_n,&sstash);
672:   ISDestroy(&all_subsets_n);

674:   /* subset indices in local boundary numbering */
675:   if (!sub_schurs->is_Ej_all) {
676:     PetscInt *all_local_idx_B;

678:     PetscMalloc1(local_size,&all_local_idx_B);
679:     ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);
680:     if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size);
681:     ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);
682:   }

684:   if (change) {
685:     ISLocalToGlobalMapping BtoS;
686:     IS                     change_primal_B;
687:     IS                     change_primal_all;

689:     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
690:     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
691:     PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);
692:     for (i=0;i<sub_schurs->n_subs;i++) {
693:       ISLocalToGlobalMapping NtoS;
694:       ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);
695:       ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);
696:       ISLocalToGlobalMappingDestroy(&NtoS);
697:     }
698:     ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);
699:     ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);
700:     ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);
701:     ISLocalToGlobalMappingDestroy(&BtoS);
702:     ISDestroy(&change_primal_B);
703:     PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);
704:     for (i=0;i<sub_schurs->n_subs;i++) {
705:       Mat change_sub;

707:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
708:       KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);
709:       KSPSetType(sub_schurs->change[i],KSPPREONLY);
710:       if (!sub_schurs->change_with_qr) {
711:         MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);
712:       } else {
713:         Mat change_subt;
714:         MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);
715:         MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);
716:         MatDestroy(&change_subt);
717:       }
718:       KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);
719:       MatDestroy(&change_sub);
720:       KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);
721:       KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");
722:     }
723:     ISDestroy(&change_primal_all);
724:   }

726:   /* Local matrix of all local Schur on subsets (transposed) */
727:   if (!sub_schurs->S_Ej_all) {
728:     Mat         T;
729:     PetscScalar *v;
730:     PetscInt    *ii,*jj;
731:     PetscInt    cum,i,j,k;

733:     /* MatSeqAIJSetPreallocation + MatSetValues is slow for these kind of matrices (may have large blocks)
734:        Allocate properly a representative matrix and duplicate */
735:     PetscMalloc3(local_size+1,&ii,local_stash_size,&jj,local_stash_size,&v);
736:     ii[0] = 0;
737:     cum   = 0;
738:     for (i=0;i<sub_schurs->n_subs;i++) {
739:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
740:       for (j=0;j<subset_size;j++) {
741:         const PetscInt row = cum+j;
742:         PetscInt col = cum;

744:         ii[row+1] = ii[row] + subset_size;
745:         for (k=ii[row];k<ii[row+1];k++) {
746:           jj[k] = col;
747:           col++;
748:         }
749:       }
750:       cum += subset_size;
751:     }
752:     MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,local_size,local_size,ii,jj,v,&T);
753:     MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&sub_schurs->S_Ej_all);
754:     MatDestroy(&T);
755:     PetscFree3(ii,jj,v);
756:   }
757:   /* matrices for deluxe scaling and adaptive selection */
758:   if (compute_Stilda) {
759:     if (!sub_schurs->sum_S_Ej_tilda_all) {
760:       MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_tilda_all);
761:     }
762:     if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
763:       MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_inv_all);
764:     }
765:   }

767:   /* Compute Schur complements explicitly */
768:   F = NULL;
769:   if (!sub_schurs->schur_explicit) {
770:     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
771:        it is not efficient, unless the economic version of the scaling is used */
772:     Mat         S_Ej_expl;
773:     PetscScalar *work;
774:     PetscInt    j,*dummy_idx;
775:     PetscBool   Sdense;

777:     PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
778:     local_size = 0;
779:     for (i=0;i<sub_schurs->n_subs;i++) {
780:       IS  is_subset_B;
781:       Mat AE_EE,AE_IE,AE_EI,S_Ej;

783:       /* subsets in original and boundary numbering */
784:       ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);
785:       /* EE block */
786:       MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);
787:       /* IE block */
788:       MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);
789:       /* EI block */
790:       if (sub_schurs->is_symmetric) {
791:         MatCreateTranspose(AE_IE,&AE_EI);
792:       } else if (sub_schurs->is_hermitian) {
793:         MatCreateHermitianTranspose(AE_IE,&AE_EI);
794:       } else {
795:         MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);
796:       }
797:       ISDestroy(&is_subset_B);
798:       MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);
799:       MatDestroy(&AE_EE);
800:       MatDestroy(&AE_IE);
801:       MatDestroy(&AE_EI);
802:       if (AE_II == A_II) { /* we can reuse the same ksp */
803:         KSP ksp;
804:         MatSchurComplementGetKSP(sub_schurs->S,&ksp);
805:         MatSchurComplementSetKSP(S_Ej,ksp);
806:       } else { /* build new ksp object which inherits ksp and pc types from the original one */
807:         KSP       origksp,schurksp;
808:         PC        origpc,schurpc;
809:         KSPType   ksp_type;
810:         PetscInt  n_internal;
811:         PetscBool ispcnone;

813:         MatSchurComplementGetKSP(sub_schurs->S,&origksp);
814:         MatSchurComplementGetKSP(S_Ej,&schurksp);
815:         KSPGetType(origksp,&ksp_type);
816:         KSPSetType(schurksp,ksp_type);
817:         KSPGetPC(schurksp,&schurpc);
818:         KSPGetPC(origksp,&origpc);
819:         PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);
820:         if (!ispcnone) {
821:           PCType pc_type;
822:           PCGetType(origpc,&pc_type);
823:           PCSetType(schurpc,pc_type);
824:         } else {
825:           PCSetType(schurpc,PCLU);
826:         }
827:         ISGetSize(is_I,&n_internal);
828:         if (!n_internal) { /* UMFPACK gives error with 0 sized problems */
829:           MatSolverType solver = NULL;
830:           PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);
831:           if (solver) {
832:             PCFactorSetMatSolverType(schurpc,solver);
833:           }
834:         }
835:         KSPSetUp(schurksp);
836:       }
837:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
838:       MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);
839:       PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);
840:       PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);
841:       if (Sdense) {
842:         for (j=0;j<subset_size;j++) {
843:           dummy_idx[j]=local_size+j;
844:         }
845:         MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
846:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
847:       MatDestroy(&S_Ej);
848:       MatDestroy(&S_Ej_expl);
849:       local_size += subset_size;
850:     }
851:     PetscFree2(dummy_idx,work);
852:     /* free */
853:     ISDestroy(&is_I);
854:     MatDestroy(&AE_II);
855:     PetscFree(all_local_idx_N);
856:   } else {
857:     Mat               A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
858:     Vec               Dall = NULL;
859:     IS                is_A_all,*is_p_r = NULL;
860:     MatType           Stype;
861:     PetscScalar       *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
862:     PetscScalar       *SEj_arr = NULL,*SEjinv_arr = NULL;
863:     const PetscScalar *rS_data;
864:     PetscInt          n,n_I,size_schur,size_active_schur,cum,cum2;
865:     PetscBool         economic,solver_S,S_lower_triangular = PETSC_FALSE;
866:     PetscBool         schur_has_vertices,factor_workaround;
867:     PetscBool         use_cholesky;
868: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
869:     PetscBool         oldpin;
870: #endif

872:     /* get sizes */
873:     n_I = 0;
874:     if (is_I_layer) {
875:       ISGetLocalSize(is_I_layer,&n_I);
876:     }
877:     economic = PETSC_FALSE;
878:     ISGetLocalSize(sub_schurs->is_I,&cum);
879:     if (cum != n_I) economic = PETSC_TRUE;
880:     MatGetLocalSize(sub_schurs->A,&n,NULL);
881:     size_active_schur = local_size;

883:     /* import scaling vector (wrong formulation if we have 3D edges) */
884:     if (scaling && compute_Stilda) {
885:       const PetscScalar *array;
886:       PetscScalar       *array2;
887:       const PetscInt    *idxs;
888:       PetscInt          i;

890:       ISGetIndices(sub_schurs->is_Ej_all,&idxs);
891:       VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);
892:       VecGetArrayRead(scaling,&array);
893:       VecGetArray(Dall,&array2);
894:       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
895:       VecRestoreArray(Dall,&array2);
896:       VecRestoreArrayRead(scaling,&array);
897:       ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);
898:       deluxe = PETSC_FALSE;
899:     }

901:     /* size active schurs does not count any dirichlet or vertex dof on the interface */
902:     factor_workaround = PETSC_FALSE;
903:     schur_has_vertices = PETSC_FALSE;
904:     cum = n_I+size_active_schur;
905:     if (sub_schurs->is_dir) {
906:       const PetscInt* idxs;
907:       PetscInt        n_dir;

909:       ISGetLocalSize(sub_schurs->is_dir,&n_dir);
910:       ISGetIndices(sub_schurs->is_dir,&idxs);
911:       PetscArraycpy(all_local_idx_N+cum,idxs,n_dir);
912:       ISRestoreIndices(sub_schurs->is_dir,&idxs);
913:       cum += n_dir;
914:       factor_workaround = PETSC_TRUE;
915:     }
916:     /* include the primal vertices in the Schur complement */
917:     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
918:       PetscInt n_v;

920:       ISGetLocalSize(sub_schurs->is_vertices,&n_v);
921:       if (n_v) {
922:         const PetscInt* idxs;

924:         ISGetIndices(sub_schurs->is_vertices,&idxs);
925:         PetscArraycpy(all_local_idx_N+cum,idxs,n_v);
926:         ISRestoreIndices(sub_schurs->is_vertices,&idxs);
927:         cum += n_v;
928:         factor_workaround = PETSC_TRUE;
929:         schur_has_vertices = PETSC_TRUE;
930:       }
931:     }
932:     size_schur = cum - n_I;
933:     ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);
934: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
935:     oldpin = sub_schurs->A->boundtocpu;
936:     MatBindToCPU(sub_schurs->A,PETSC_TRUE);
937: #endif
938:     if (cum == n) {
939:       ISSetPermutation(is_A_all);
940:       MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);
941:     } else {
942:       MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);
943:     }
944: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
945:     MatBindToCPU(sub_schurs->A,oldpin);
946: #endif
947:     MatSetOptionsPrefix(A,sub_schurs->prefix);
948:     MatAppendOptionsPrefix(A,"sub_schurs_");

950:     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
951:        this is a workaround */
952:     if (benign_n) {
953:       Vec                    v,benign_AIIm1_ones;
954:       ISLocalToGlobalMapping N_to_reor;
955:       IS                     is_p0,is_p0_p;
956:       PetscScalar            *cs_AIB,*AIIm1_data;
957:       PetscInt               sizeA;

959:       ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);
960:       ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);
961:       ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);
962:       ISDestroy(&is_p0);
963:       MatCreateVecs(A,&v,&benign_AIIm1_ones);
964:       VecGetSize(v,&sizeA);
965:       MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);
966:       MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);
967:       MatDenseGetArray(cs_AIB_mat,&cs_AIB);
968:       MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
969:       PetscMalloc1(benign_n,&is_p_r);
970:       /* compute colsum of A_IB restricted to pressures */
971:       for (i=0;i<benign_n;i++) {
972:         const PetscScalar *array;
973:         const PetscInt    *idxs;
974:         PetscInt          j,nz;

976:         ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);
977:         ISGetLocalSize(is_p_r[i],&nz);
978:         ISGetIndices(is_p_r[i],&idxs);
979:         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
980:         ISRestoreIndices(is_p_r[i],&idxs);
981:         VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);
982:         MatMult(A,benign_AIIm1_ones,v);
983:         VecResetArray(benign_AIIm1_ones);
984:         VecGetArrayRead(v,&array);
985:         for (j=0;j<size_schur;j++) {
986: #if defined(PETSC_USE_COMPLEX)
987:           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
988: #else
989:           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
990: #endif
991:         }
992:         VecRestoreArrayRead(v,&array);
993:       }
994:       MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
995:       MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
996:       VecDestroy(&v);
997:       VecDestroy(&benign_AIIm1_ones);
998:       MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);
999:       MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
1000:       MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
1001:       MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);
1002:       ISDestroy(&is_p0_p);
1003:       ISLocalToGlobalMappingDestroy(&N_to_reor);
1004:     }
1005:     MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);
1006:     MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);
1007:     MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);

1009:     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
1010:     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);

1012:     /* when using the benign subspace trick, the local Schur complements are SPD */
1013:     /* MKL_PARDISO does not handle well the computation of a Schur complement from a symmetric indefinite factorization
1014:        Use LU and adapt pivoting perturbation (still, solution is not as accurate as with using MUMPS) */
1015:     if (benign_trick) {
1016:       sub_schurs->is_posdef = PETSC_TRUE;
1017:       PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,&flg);
1018:       if (flg) use_cholesky = PETSC_FALSE;
1019:     }

1021:     if (n_I) {
1022:       IS        is_schur;
1023:       char      stype[64];
1024:       PetscBool gpu = PETSC_FALSE;

1026:       if (use_cholesky) {
1027:         MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);
1028:       } else {
1029:         MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);
1030:       }
1031:       MatSetErrorIfFailure(A,PETSC_TRUE);
1032: #if defined(PETSC_HAVE_MKL_PARDISO)
1033:       if (benign_trick) { MatMkl_PardisoSetCntl(F,10,10); }
1034: #endif
1035:       /* subsets ordered last */
1036:       ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);
1037:       MatFactorSetSchurIS(F,is_schur);
1038:       ISDestroy(&is_schur);

1040:       /* factorization step */
1041:       if (use_cholesky) {
1042:         MatCholeskyFactorSymbolic(F,A,NULL,NULL);
1043: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1044:         MatMumpsSetIcntl(F,19,2);
1045: #endif
1046:         MatCholeskyFactorNumeric(F,A,NULL);
1047:         S_lower_triangular = PETSC_TRUE;
1048:       } else {
1049:         MatLUFactorSymbolic(F,A,NULL,NULL,NULL);
1050: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1051:         MatMumpsSetIcntl(F,19,3);
1052: #endif
1053:         MatLUFactorNumeric(F,A,NULL);
1054:       }
1055:       MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");

1057:       if (matl_dbg_viewer) {
1058:         Mat S;
1059:         IS  is;

1061:         PetscObjectSetName((PetscObject)A,"A");
1062:         MatView(A,matl_dbg_viewer);
1063:         MatFactorCreateSchurComplement(F,&S,NULL);
1064:         PetscObjectSetName((PetscObject)S,"S");
1065:         MatView(S,matl_dbg_viewer);
1066:         MatDestroy(&S);
1067:         ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);
1068:         PetscObjectSetName((PetscObject)is,"I");
1069:         ISView(is,matl_dbg_viewer);
1070:         ISDestroy(&is);
1071:         ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);
1072:         PetscObjectSetName((PetscObject)is,"B");
1073:         ISView(is,matl_dbg_viewer);
1074:         ISDestroy(&is);
1075:         PetscObjectSetName((PetscObject)is_A_all,"IA");
1076:         ISView(is_A_all,matl_dbg_viewer);
1077:       }

1079:       /* get explicit Schur Complement computed during numeric factorization */
1080:       MatFactorGetSchurComplement(F,&S_all,NULL);
1081:       PetscStrncpy(stype,MATSEQDENSE,sizeof(stype));
1082: #if defined(PETSC_HAVE_CUDA)
1083:       PetscObjectTypeCompareAny((PetscObject)A,&gpu,MATSEQAIJVIENNACL,MATSEQAIJCUSPARSE,"");
1084: #endif
1085:       if (gpu) {
1086:         PetscStrncpy(stype,MATSEQDENSECUDA,sizeof(stype));
1087:       }
1088:       PetscOptionsGetString(NULL,sub_schurs->prefix,"-sub_schurs_schur_mat_type",stype,sizeof(stype),NULL);
1089:       MatConvert(S_all,stype,MAT_INPLACE_MATRIX,&S_all);
1090:       MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);
1091:       MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);
1092:       MatGetType(S_all,&Stype);

1094:       /* we can reuse the solvers if we are not using the economic version */
1095:       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1096:       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1097:       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
1098:         reuse_solvers = factor_workaround = PETSC_FALSE;

1100:       solver_S = PETSC_TRUE;

1102:       /* update the Schur complement with the change of basis on the pressures */
1103:       if (benign_n) {
1104:         const PetscScalar *cs_AIB;
1105:         PetscScalar       *S_data,*AIIm1_data;
1106:         Mat               S2 = NULL,S3 = NULL; /* dbg */
1107:         PetscScalar       *S2_data,*S3_data; /* dbg */
1108:         Vec               v,benign_AIIm1_ones;
1109:         PetscInt          sizeA;

1111:         MatDenseGetArray(S_all,&S_data);
1112:         MatCreateVecs(A,&v,&benign_AIIm1_ones);
1113:         VecGetSize(v,&sizeA);
1114: #if defined(PETSC_HAVE_MUMPS)
1115:         MatMumpsSetIcntl(F,26,0);
1116: #endif
1117: #if defined(PETSC_HAVE_MKL_PARDISO)
1118:         MatMkl_PardisoSetCntl(F,70,1);
1119: #endif
1120:         MatDenseGetArrayRead(cs_AIB_mat,&cs_AIB);
1121:         MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
1122:         if (matl_dbg_viewer) {
1123:           MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S2);
1124:           MatDuplicate(S_all,MAT_DO_NOT_COPY_VALUES,&S3);
1125:           MatDenseGetArray(S2,&S2_data);
1126:           MatDenseGetArray(S3,&S3_data);
1127:         }
1128:         for (i=0;i<benign_n;i++) {
1129:           PetscScalar    *array,sum = 0.,one = 1.,*sums;
1130:           const PetscInt *idxs;
1131:           PetscInt       k,j,nz;
1132:           PetscBLASInt   B_k,B_n;

1134:           PetscCalloc1(benign_n,&sums);
1135:           VecPlaceArray(benign_AIIm1_ones,AIIm1_data+sizeA*i);
1136:           VecCopy(benign_AIIm1_ones,v);
1137:           MatSolve(F,v,benign_AIIm1_ones);
1138:           MatMult(A,benign_AIIm1_ones,v);
1139:           VecResetArray(benign_AIIm1_ones);
1140:           /* p0 dofs (eliminated) are excluded from the sums */
1141:           for (k=0;k<benign_n;k++) {
1142:             ISGetLocalSize(is_p_r[k],&nz);
1143:             ISGetIndices(is_p_r[k],&idxs);
1144:             for (j=0;j<nz-1;j++) sums[k] -= AIIm1_data[idxs[j]+sizeA*i];
1145:             ISRestoreIndices(is_p_r[k],&idxs);
1146:           }
1147:           VecGetArrayRead(v,(const PetscScalar**)&array);
1148:           if (matl_dbg_viewer) {
1149:             Vec  vv;
1150:             char name[16];

1152:             VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array+n_I,&vv);
1153:             PetscSNPrintf(name,sizeof(name),"Pvs%D",i);
1154:             PetscObjectSetName((PetscObject)vv,name);
1155:             VecView(vv,matl_dbg_viewer);
1156:           }
1157:           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1158:           /* cs_AIB already scaled by 1./nz */
1159:           B_k = 1;
1160:           for (k=0;k<benign_n;k++) {
1161:             sum  = sums[k];
1162:             PetscBLASIntCast(size_schur,&B_n);

1164:             if (PetscAbsScalar(sum) == 0.0) continue;
1165:             if (k == i) {
1166:               PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1167:               if (matl_dbg_viewer) {
1168:                 PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S3_data,&B_n));
1169:               }
1170:             } else { /* XXX Is it correct to use symmetric rank-2 update with half of the sum? */
1171:               sum /= 2.0;
1172:               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));
1173:               if (matl_dbg_viewer) {
1174:                 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));
1175:               }
1176:             }
1177:           }
1178:           sum  = 1.;
1179:           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));
1180:           if (matl_dbg_viewer) {
1181:             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));
1182:           }
1183:           VecRestoreArrayRead(v,(const PetscScalar**)&array);
1184:           /* set p0 entry of AIIm1_ones to zero */
1185:           ISGetLocalSize(is_p_r[i],&nz);
1186:           ISGetIndices(is_p_r[i],&idxs);
1187:           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1188:           ISRestoreIndices(is_p_r[i],&idxs);
1189:           PetscFree(sums);
1190:         }
1191:         VecDestroy(&benign_AIIm1_ones);
1192:         if (matl_dbg_viewer) {
1193:           MatDenseRestoreArray(S2,&S2_data);
1194:           MatDenseRestoreArray(S3,&S3_data);
1195:         }
1196:         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1197:           PetscInt k,j;
1198:           for (k=0;k<size_schur;k++) {
1199:             for (j=k;j<size_schur;j++) {
1200:               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1201:             }
1202:           }
1203:         }

1205:         /* restore defaults */
1206: #if defined(PETSC_HAVE_MUMPS)
1207:         MatMumpsSetIcntl(F,26,-1);
1208: #endif
1209: #if defined(PETSC_HAVE_MKL_PARDISO)
1210:         MatMkl_PardisoSetCntl(F,70,0);
1211: #endif
1212:         MatDenseRestoreArrayRead(cs_AIB_mat,&cs_AIB);
1213:         MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
1214:         VecDestroy(&v);
1215:         MatDenseRestoreArray(S_all,&S_data);
1216:         if (matl_dbg_viewer) {
1217:           Mat S;

1219:           MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1220:           MatFactorCreateSchurComplement(F,&S,NULL);
1221:           PetscObjectSetName((PetscObject)S,"Sb");
1222:           MatView(S,matl_dbg_viewer);
1223:           MatDestroy(&S);
1224:           PetscObjectSetName((PetscObject)S2,"S2P");
1225:           MatView(S2,matl_dbg_viewer);
1226:           PetscObjectSetName((PetscObject)S3,"S3P");
1227:           MatView(S3,matl_dbg_viewer);
1228:           PetscObjectSetName((PetscObject)cs_AIB_mat,"cs");
1229:           MatView(cs_AIB_mat,matl_dbg_viewer);
1230:           MatFactorGetSchurComplement(F,&S_all,NULL);
1231:         }
1232:         MatDestroy(&S2);
1233:         MatDestroy(&S3);
1234:       }
1235:       if (!reuse_solvers) {
1236:         for (i=0;i<benign_n;i++) {
1237:           ISDestroy(&is_p_r[i]);
1238:         }
1239:         PetscFree(is_p_r);
1240:         MatDestroy(&cs_AIB_mat);
1241:         MatDestroy(&benign_AIIm1_ones_mat);
1242:       }
1243:     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1244:       MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);
1245:       MatGetType(S_all,&Stype);
1246:       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1247:       factor_workaround = PETSC_FALSE;
1248:       solver_S = PETSC_FALSE;
1249:     }

1251:     if (reuse_solvers) {
1252:       Mat                A_II,Afake;
1253:       Vec                vec1_B;
1254:       PCBDDCReuseSolvers msolv_ctx;
1255:       PetscInt           n_R;

1257:       if (sub_schurs->reuse_solver) {
1258:         PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1259:       } else {
1260:         PetscNew(&sub_schurs->reuse_solver);
1261:       }
1262:       msolv_ctx = sub_schurs->reuse_solver;
1263:       MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);
1264:       PetscObjectReference((PetscObject)F);
1265:       msolv_ctx->F = F;
1266:       MatCreateVecs(F,&msolv_ctx->sol,NULL);
1267:       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1268:       {
1269:         PetscScalar *array;
1270:         PetscInt    n;

1272:         VecGetLocalSize(msolv_ctx->sol,&n);
1273:         VecGetArray(msolv_ctx->sol,&array);
1274:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);
1275:         VecRestoreArray(msolv_ctx->sol,&array);
1276:       }
1277:       msolv_ctx->has_vertices = schur_has_vertices;

1279:       /* interior solver */
1280:       PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);
1281:       PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);
1282:       PCSetType(msolv_ctx->interior_solver,PCSHELL);
1283:       PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");
1284:       PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);
1285:       PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);
1286:       PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);
1287:       PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);

1289:       /* correction solver */
1290:       PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);
1291:       PCSetType(msolv_ctx->correction_solver,PCSHELL);
1292:       PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");
1293:       PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);
1294:       PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);
1295:       PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);
1296:       PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);

1298:       /* scatter and vecs for Schur complement solver */
1299:       MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);
1300:       MatCreateVecs(sub_schurs->S,&vec1_B,NULL);
1301:       if (!schur_has_vertices) {
1302:         ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);
1303:         VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);
1304:         PetscObjectReference((PetscObject)is_A_all);
1305:         msolv_ctx->is_R = is_A_all;
1306:       } else {
1307:         IS              is_B_all;
1308:         const PetscInt* idxs;
1309:         PetscInt        dual,n_v,n;

1311:         ISGetLocalSize(sub_schurs->is_vertices,&n_v);
1312:         dual = size_schur - n_v;
1313:         ISGetLocalSize(is_A_all,&n);
1314:         ISGetIndices(is_A_all,&idxs);
1315:         ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);
1316:         ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);
1317:         ISDestroy(&is_B_all);
1318:         ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);
1319:         VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);
1320:         ISDestroy(&is_B_all);
1321:         ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);
1322:         ISRestoreIndices(is_A_all,&idxs);
1323:       }
1324:       ISGetLocalSize(msolv_ctx->is_R,&n_R);
1325:       MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);
1326:       MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);
1327:       MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);
1328:       PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);
1329:       MatDestroy(&Afake);
1330:       VecDestroy(&vec1_B);

1332:       /* communicate benign info to solver context */
1333:       if (benign_n) {
1334:         PetscScalar *array;

1336:         msolv_ctx->benign_n = benign_n;
1337:         msolv_ctx->benign_zerodiag_subs = is_p_r;
1338:         PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);
1339:         msolv_ctx->benign_csAIB = cs_AIB_mat;
1340:         MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);
1341:         VecGetArray(msolv_ctx->benign_corr_work,&array);
1342:         VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);
1343:         VecRestoreArray(msolv_ctx->benign_corr_work,&array);
1344:         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1345:       }
1346:     } else {
1347:       if (sub_schurs->reuse_solver) {
1348:         PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1349:       }
1350:       PetscFree(sub_schurs->reuse_solver);
1351:     }
1352:     MatDestroy(&A);
1353:     ISDestroy(&is_A_all);

1355:     /* Work arrays */
1356:     PetscMalloc1(max_subset_size*max_subset_size,&work);

1358:     /* S_Ej_all */
1359:     cum = cum2 = 0;
1360:     MatDenseGetArrayRead(S_all,&rS_data);
1361:     MatSeqAIJGetArray(sub_schurs->S_Ej_all,&SEj_arr);
1362:     if (compute_Stilda) {
1363:       MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&SEjinv_arr);
1364:     }
1365:     for (i=0;i<sub_schurs->n_subs;i++) {
1366:       PetscInt j;

1368:       /* get S_E */
1369:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1370:       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1371:         PetscInt k;
1372:         for (k=0;k<subset_size;k++) {
1373:           for (j=k;j<subset_size;j++) {
1374:             work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1375:             work[j*subset_size+k] = PetscConj(rS_data[cum2+k*size_schur+j]);
1376:           }
1377:         }
1378:       } else { /* just copy to workspace */
1379:         PetscInt k;
1380:         for (k=0;k<subset_size;k++) {
1381:           for (j=0;j<subset_size;j++) {
1382:             work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1383:           }
1384:         }
1385:       }
1386:       /* insert S_E values */
1387:       if (sub_schurs->change) {
1388:         Mat change_sub,SEj,T;

1390:         /* change basis */
1391:         KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1392:         MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1393:         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1394:           Mat T2;
1395:           MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1396:           MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1397:           MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1398:           MatDestroy(&T2);
1399:         } else {
1400:           MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1401:         }
1402:         MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1403:         MatDestroy(&T);
1404:         MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);
1405:         MatDestroy(&SEj);
1406:       }
1407:       if (deluxe) {
1408:         PetscArraycpy(SEj_arr,work,subset_size*subset_size);
1409:         /* if adaptivity is requested, invert S_E blocks */
1410:         if (compute_Stilda) {
1411:           Mat               M;
1412:           const PetscScalar *vals;
1413:           PetscBool         isdense,isdensecuda;

1415:           MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&M);
1416:           MatSetOption(M,MAT_SPD,sub_schurs->is_posdef);
1417:           MatSetOption(M,MAT_HERMITIAN,sub_schurs->is_hermitian);
1418:           if (!PetscBTLookup(sub_schurs->is_edge,i)) {
1419:             MatSetType(M,Stype);
1420:           }
1421:           PetscObjectTypeCompare((PetscObject)M,MATSEQDENSE,&isdense);
1422:           PetscObjectTypeCompare((PetscObject)M,MATSEQDENSECUDA,&isdensecuda);
1423:           if (use_cholesky) {
1424:             MatCholeskyFactor(M,NULL,NULL);
1425:           } else {
1426:             MatLUFactor(M,NULL,NULL,NULL);
1427:           }
1428:           if (isdense) {
1429:             MatSeqDenseInvertFactors_Private(M);
1430: #if defined(PETSC_HAVE_CUDA)
1431:           } else if (isdensecuda) {
1432:             MatSeqDenseCUDAInvertFactors_Private(M);
1433: #endif
1434:           } else SETERRQ1(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"Not implemented for type %s",Stype);
1435:           MatDenseGetArrayRead(M,&vals);
1436:           PetscArraycpy(SEjinv_arr,vals,subset_size*subset_size);
1437:           MatDenseRestoreArrayRead(M,&vals);
1438:           MatDestroy(&M);
1439:         }
1440:       } else if (compute_Stilda) { /* not using deluxe */
1441:         Mat         SEj;
1442:         Vec         D;
1443:         PetscScalar *array;

1445:         MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1446:         VecGetArray(Dall,&array);
1447:         VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);
1448:         VecRestoreArray(Dall,&array);
1449:         VecShift(D,-1.);
1450:         MatDiagonalScale(SEj,D,D);
1451:         MatDestroy(&SEj);
1452:         VecDestroy(&D);
1453:         PetscArraycpy(SEj_arr,work,subset_size*subset_size);
1454:       }
1455:       cum += subset_size;
1456:       cum2 += subset_size*(size_schur + 1);
1457:       SEj_arr += subset_size*subset_size;
1458:       if (SEjinv_arr) SEjinv_arr += subset_size*subset_size;
1459:     }
1460:     MatDenseRestoreArrayRead(S_all,&rS_data);
1461:     MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&SEj_arr);
1462:     if (compute_Stilda) {
1463:       MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&SEjinv_arr);
1464:     }
1465:     if (solver_S) {
1466:       MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1467:     }

1469:     /* may prevent from unneeded copies, since MUMPS or MKL_Pardiso always use CPU memory
1470:        however, preliminary tests indicate using GPUs is still faster in the solve phase */
1471: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1472:     if (reuse_solvers) {
1473:       Mat                  St;
1474:       MatFactorSchurStatus st;

1476:       flg  = PETSC_FALSE;
1477:       PetscOptionsGetBool(NULL,sub_schurs->prefix,"-sub_schurs_schur_pin_to_cpu",&flg,NULL);
1478:       MatFactorGetSchurComplement(F,&St,&st);
1479:       MatBindToCPU(St,flg);
1480:       MatFactorRestoreSchurComplement(F,&St,st);
1481:     }
1482: #endif

1484:     schur_factor = NULL;
1485:     if (compute_Stilda && size_active_schur) {

1487:       MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&SEjinv_arr);
1488:       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1489:         PetscArraycpy(SEjinv_arr,work,size_schur*size_schur);
1490:       } else {
1491:         Mat S_all_inv=NULL;

1493:         if (solver_S) {
1494:           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1495:              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 */
1496:           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1497:             PetscScalar *data;
1498:             PetscInt     nd = 0;

1500:             if (!use_potr) {
1501:               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1502:             }
1503:             MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1504:             MatDenseGetArray(S_all_inv,&data);
1505:             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1506:               ISGetLocalSize(sub_schurs->is_dir,&nd);
1507:             }

1509:             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1510:             if (schur_has_vertices) {
1511:               Mat          M;
1512:               PetscScalar *tdata;
1513:               PetscInt     nv = 0, news;

1515:               ISGetLocalSize(sub_schurs->is_vertices,&nv);
1516:               news = size_active_schur + nv;
1517:               PetscCalloc1(news*news,&tdata);
1518:               for (i=0;i<size_active_schur;i++) {
1519:                 PetscArraycpy(tdata+i*(news+1),data+i*(size_schur+1),size_active_schur-i);
1520:                 PetscArraycpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv);
1521:               }
1522:               for (i=0;i<nv;i++) {
1523:                 PetscInt k = i+size_active_schur;
1524:                 PetscArraycpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),nv-i);
1525:               }

1527:               MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);
1528:               MatSetOption(M,MAT_SPD,PETSC_TRUE);
1529:               MatCholeskyFactor(M,NULL,NULL);
1530:               /* save the factors */
1531:               cum = 0;
1532:               PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);
1533:               for (i=0;i<size_active_schur;i++) {
1534:                 PetscArraycpy(schur_factor+cum,tdata+i*(news+1),size_active_schur-i);
1535:                 cum += size_active_schur - i;
1536:               }
1537:               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1538:               MatSeqDenseInvertFactors_Private(M);
1539:               /* move back just the active dofs to the Schur complement */
1540:               for (i=0;i<size_active_schur;i++) {
1541:                 PetscArraycpy(data+i*size_schur,tdata+i*news,size_active_schur);
1542:               }
1543:               PetscFree(tdata);
1544:               MatDestroy(&M);
1545:             } else { /* we can factorize and invert just the activedofs */
1546:               Mat         M;
1547:               PetscScalar *aux;

1549:               PetscMalloc1(nd,&aux);
1550:               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1551:               MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);
1552:               MatDenseSetLDA(M,size_schur);
1553:               MatSetOption(M,MAT_SPD,PETSC_TRUE);
1554:               MatCholeskyFactor(M,NULL,NULL);
1555:               MatSeqDenseInvertFactors_Private(M);
1556:               MatDestroy(&M);
1557:               MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);
1558:               MatZeroEntries(M);
1559:               MatDestroy(&M);
1560:               MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);
1561:               MatDenseSetLDA(M,size_schur);
1562:               MatZeroEntries(M);
1563:               MatDestroy(&M);
1564:               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
1565:               PetscFree(aux);
1566:             }
1567:             MatDenseRestoreArray(S_all_inv,&data);
1568:           } else { /* use MatFactor calls to invert S */
1569:             MatFactorInvertSchurComplement(F);
1570:             MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1571:           }
1572:         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1573:           PetscObjectReference((PetscObject)S_all);
1574:           S_all_inv = S_all;
1575:           MatDenseGetArray(S_all_inv,&S_data);
1576:           PetscBLASIntCast(size_schur,&B_N);
1577:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1578:           if (use_potr) {
1579:             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1580:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1581:             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1582:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1583:           } else if (use_sytr) {
1584:             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1585:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1586:             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1587:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1588:           } else {
1589:             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1590:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1591:             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1592:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1593:           }
1594:           PetscLogFlops(1.0*size_schur*size_schur*size_schur);
1595:           PetscFPTrapPop();
1596:           MatDenseRestoreArray(S_all_inv,&S_data);
1597:         }
1598:         /* S_Ej_tilda_all */
1599:         cum = cum2 = 0;
1600:         MatDenseGetArrayRead(S_all_inv,&rS_data);
1601:         for (i=0;i<sub_schurs->n_subs;i++) {
1602:           PetscInt j;

1604:           ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1605:           /* get (St^-1)_E */
1606:           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1607:              will be properly accessed later during adaptive selection */
1608:           if (S_lower_triangular) {
1609:             PetscInt k;
1610:             if (sub_schurs->change) {
1611:               for (k=0;k<subset_size;k++) {
1612:                 for (j=k;j<subset_size;j++) {
1613:                   work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1614:                   work[j*subset_size+k] = work[k*subset_size+j];
1615:                 }
1616:               }
1617:             } else {
1618:               for (k=0;k<subset_size;k++) {
1619:                 for (j=k;j<subset_size;j++) {
1620:                   work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1621:                 }
1622:               }
1623:             }
1624:           } else {
1625:             PetscInt k;
1626:             for (k=0;k<subset_size;k++) {
1627:               for (j=0;j<subset_size;j++) {
1628:                 work[k*subset_size+j] = rS_data[cum2+k*size_schur+j];
1629:               }
1630:             }
1631:           }
1632:           if (sub_schurs->change) {
1633:             Mat change_sub,SEj,T;

1635:             /* change basis */
1636:             KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1637:             MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1638:             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1639:               Mat T2;
1640:               MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1641:               MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1642:               MatDestroy(&T2);
1643:               MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1644:             } else {
1645:               MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1646:             }
1647:             MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1648:             MatDestroy(&T);
1649:             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1650:             MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);
1651:             MatDestroy(&SEj);
1652:           }
1653:           PetscArraycpy(SEjinv_arr,work,subset_size*subset_size);
1654:           cum += subset_size;
1655:           cum2 += subset_size*(size_schur + 1);
1656:           SEjinv_arr += subset_size*subset_size;
1657:         }
1658:         MatDenseRestoreArrayRead(S_all_inv,&rS_data);
1659:         if (solver_S) {
1660:           if (schur_has_vertices) {
1661:             MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);
1662:           } else {
1663:             MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);
1664:           }
1665:         }
1666:         MatDestroy(&S_all_inv);
1667:       }
1668:       MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&SEjinv_arr);

1670:       /* move back factors if needed */
1671:       if (schur_has_vertices) {
1672:         Mat      S_tmp;
1673:         PetscInt nd = 0;

1675:         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1676:         MatFactorGetSchurComplement(F,&S_tmp,NULL);
1677:         if (use_potr) {
1678:           PetscScalar *data;

1680:           MatDenseGetArray(S_tmp,&data);
1681:           PetscArrayzero(data,size_schur*size_schur);

1683:           if (S_lower_triangular) {
1684:             cum = 0;
1685:             for (i=0;i<size_active_schur;i++) {
1686:               PetscArraycpy(data+i*(size_schur+1),schur_factor+cum,size_active_schur-i);
1687:               cum += size_active_schur-i;
1688:             }
1689:           } else {
1690:             PetscArraycpy(data,schur_factor,size_schur*size_schur);
1691:           }
1692:           if (sub_schurs->is_dir) {
1693:             ISGetLocalSize(sub_schurs->is_dir,&nd);
1694:             for (i=0;i<nd;i++) {
1695:               data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1696:             }
1697:           }
1698:           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1699:              set the diagonal entry of the Schur factor to a very large value */
1700:           for (i=size_active_schur+nd;i<size_schur;i++) {
1701:             data[i*(size_schur+1)] = infty;
1702:           }
1703:           MatDenseRestoreArray(S_tmp,&data);
1704:         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1705:         MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);
1706:       }
1707:     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1708:       PetscScalar *data;
1709:       PetscInt    nd = 0;

1711:       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1712:         ISGetLocalSize(sub_schurs->is_dir,&nd);
1713:       }
1714:       MatFactorGetSchurComplement(F,&S_all,NULL);
1715:       MatDenseGetArray(S_all,&data);
1716:       for (i=0;i<size_active_schur;i++) {
1717:         PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);
1718:       }
1719:       for (i=size_active_schur+nd;i<size_schur;i++) {
1720:         PetscArrayzero(data+i*size_schur+size_active_schur,size_schur-size_active_schur);
1721:         data[i*(size_schur+1)] = infty;
1722:       }
1723:       MatDenseRestoreArray(S_all,&data);
1724:       MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1725:     }
1726:     PetscFree(work);
1727:     PetscFree(schur_factor);
1728:     VecDestroy(&Dall);
1729:   }
1730:   ISDestroy(&is_I_layer);
1731:   MatDestroy(&S_all);
1732:   MatDestroy(&A_BB);
1733:   MatDestroy(&A_IB);
1734:   MatDestroy(&A_BI);
1735:   MatDestroy(&F);

1737:   PetscMalloc1(sub_schurs->n_subs,&nnz);
1738:   for (i=0;i<sub_schurs->n_subs;i++) {
1739:     ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);
1740:   }
1741:   ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);
1742:   MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);
1743:   MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1744:   MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1745:   if (compute_Stilda) {
1746:     MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);
1747:     MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1748:     MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1749:     if (deluxe) {
1750:       MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);
1751:       MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1752:       MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1753:     }
1754:   }
1755:   ISDestroy(&is_I_layer);

1757:   /* Get local part of (\sum_j S_Ej) */
1758:   if (!sub_schurs->sum_S_Ej_all) {
1759:     MatDuplicate(sub_schurs->S_Ej_all,MAT_DO_NOT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);
1760:   }
1761:   VecSet(gstash,0.0);
1762:   MatSeqAIJGetArray(sub_schurs->S_Ej_all,&stasharray);
1763:   VecPlaceArray(lstash,stasharray);
1764:   VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1765:   VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1766:   MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&stasharray);
1767:   VecResetArray(lstash);
1768:   MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&stasharray);
1769:   VecPlaceArray(lstash,stasharray);
1770:   VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1771:   VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1772:   MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&stasharray);
1773:   VecResetArray(lstash);

1775:   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1776:   if (compute_Stilda) {
1777:     VecSet(gstash,0.0);
1778:     MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);
1779:     VecPlaceArray(lstash,stasharray);
1780:     VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1781:     VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1782:     VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1783:     VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1784:     MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&stasharray);
1785:     VecResetArray(lstash);
1786:     if (deluxe) {
1787:       VecSet(gstash,0.0);
1788:       MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);
1789:       VecPlaceArray(lstash,stasharray);
1790:       VecScatterBegin(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1791:       VecScatterEnd(sstash,lstash,gstash,ADD_VALUES,SCATTER_FORWARD);
1792:       VecScatterBegin(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1793:       VecScatterEnd(sstash,gstash,lstash,INSERT_VALUES,SCATTER_REVERSE);
1794:       MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&stasharray);
1795:       VecResetArray(lstash);
1796:     } else {
1797:       PetscScalar *array;
1798:       PetscInt    cum;

1800:       MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1801:       cum = 0;
1802:       for (i=0;i<sub_schurs->n_subs;i++) {
1803:         ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1804:         PetscBLASIntCast(subset_size,&B_N);
1805:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1806:         if (use_potr) {
1807:           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1808:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1809:           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1810:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1811:         } else if (use_sytr) {
1812:           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1813:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1814:           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1815:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1816:         } else {
1817:           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1818:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1819:           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1820:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1821:         }
1822:         PetscLogFlops(1.0*subset_size*subset_size*subset_size);
1823:         PetscFPTrapPop();
1824:         cum += subset_size*subset_size;
1825:       }
1826:       MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1827:       PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);
1828:       MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1829:       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1830:     }
1831:   }
1832:   VecDestroy(&lstash);
1833:   VecDestroy(&gstash);
1834:   VecScatterDestroy(&sstash);

1836:   if (matl_dbg_viewer) {
1837:     PetscInt cum;

1839:     if (sub_schurs->S_Ej_all) {
1840:       PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");
1841:       MatView(sub_schurs->S_Ej_all,matl_dbg_viewer);
1842:     }
1843:     if (sub_schurs->sum_S_Ej_all) {
1844:       PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");
1845:       MatView(sub_schurs->sum_S_Ej_all,matl_dbg_viewer);
1846:     }
1847:     if (sub_schurs->sum_S_Ej_inv_all) {
1848:       PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");
1849:       MatView(sub_schurs->sum_S_Ej_inv_all,matl_dbg_viewer);
1850:     }
1851:     if (sub_schurs->sum_S_Ej_tilda_all) {
1852:       PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");
1853:       MatView(sub_schurs->sum_S_Ej_tilda_all,matl_dbg_viewer);
1854:     }
1855:     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
1856:       IS   is;
1857:       char name[16];

1859:       PetscSNPrintf(name,sizeof(name),"IE%D",i);
1860:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1861:       ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);
1862:       PetscObjectSetName((PetscObject)is,name);
1863:       ISView(is,matl_dbg_viewer);
1864:       ISDestroy(&is);
1865:       cum += subset_size;
1866:     }
1867:   }

1869:   /* free workspace */
1870:   PetscViewerDestroy(&matl_dbg_viewer);
1871:   PetscFree2(Bwork,pivots);
1872:   PetscCommDestroy(&comm_n);
1873:   return(0);
1874: }

1876: PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1877: {
1878:   IS              *faces,*edges,*all_cc,vertices;
1879:   PetscInt        i,n_faces,n_edges,n_all_cc;
1880:   PetscBool       is_sorted,ispardiso,ismumps;
1881:   PetscErrorCode  ierr;

1884:   ISSorted(is_I,&is_sorted);
1885:   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1886:   ISSorted(is_B,&is_sorted);
1887:   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");

1889:   /* reset any previous data */
1890:   PCBDDCSubSchursReset(sub_schurs);

1892:   /* get index sets for faces and edges (already sorted by global ordering) */
1893:   PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1894:   n_all_cc = n_faces+n_edges;
1895:   PetscBTCreate(n_all_cc,&sub_schurs->is_edge);
1896:   PetscMalloc1(n_all_cc,&all_cc);
1897:   for (i=0;i<n_faces;i++) {
1898:     if (copycc) {
1899:       ISDuplicate(faces[i],&all_cc[i]);
1900:     } else {
1901:       PetscObjectReference((PetscObject)faces[i]);
1902:       all_cc[i] = faces[i];
1903:     }
1904:   }
1905:   for (i=0;i<n_edges;i++) {
1906:     if (copycc) {
1907:       ISDuplicate(edges[i],&all_cc[n_faces+i]);
1908:     } else {
1909:       PetscObjectReference((PetscObject)edges[i]);
1910:       all_cc[n_faces+i] = edges[i];
1911:     }
1912:     PetscBTSet(sub_schurs->is_edge,n_faces+i);
1913:   }
1914:   PetscObjectReference((PetscObject)vertices);
1915:   sub_schurs->is_vertices = vertices;
1916:   PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1917:   sub_schurs->is_dir = NULL;
1918:   PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);

1920:   /* Determine if MatFactor can be used */
1921:   PetscStrallocpy(prefix,&sub_schurs->prefix);
1922: #if defined(PETSC_HAVE_MUMPS)
1923:   PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,sizeof(sub_schurs->mat_solver_type));
1924: #elif defined(PETSC_HAVE_MKL_PARDISO)
1925:   PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,sizeof(sub_schurs->mat_solver_type));
1926: #else
1927:   PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,sizeof(sub_schurs->mat_solver_type));
1928: #endif
1929: #if defined(PETSC_USE_COMPLEX)
1930:   sub_schurs->is_hermitian  = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1931: #else
1932:   sub_schurs->is_hermitian  = PETSC_TRUE;
1933: #endif
1934:   sub_schurs->is_posdef     = PETSC_TRUE;
1935:   sub_schurs->is_symmetric  = PETSC_TRUE;
1936:   sub_schurs->debug         = PETSC_FALSE;
1937:   sub_schurs->restrict_comm = PETSC_FALSE;
1938:   PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
1939:   PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,sizeof(sub_schurs->mat_solver_type),NULL);
1940:   PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);
1941:   PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);
1942:   PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);
1943:   PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);
1944:   PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);
1945:   PetscOptionsEnd();
1946:   PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMUMPS,&ismumps);
1947:   PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,&ispardiso);
1948:   sub_schurs->schur_explicit = (PetscBool)(ispardiso || ismumps);

1950:   /* for reals, symmetric and hermitian are synonims */
1951: #if !defined(PETSC_USE_COMPLEX)
1952:   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1953:   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1954: #endif

1956:   PetscObjectReference((PetscObject)is_I);
1957:   sub_schurs->is_I = is_I;
1958:   PetscObjectReference((PetscObject)is_B);
1959:   sub_schurs->is_B = is_B;
1960:   PetscObjectReference((PetscObject)graph->l2gmap);
1961:   sub_schurs->l2gmap = graph->l2gmap;
1962:   PetscObjectReference((PetscObject)BtoNmap);
1963:   sub_schurs->BtoNmap = BtoNmap;
1964:   sub_schurs->n_subs = n_all_cc;
1965:   sub_schurs->is_subs = all_cc;
1966:   sub_schurs->S_Ej_all = NULL;
1967:   sub_schurs->sum_S_Ej_all = NULL;
1968:   sub_schurs->sum_S_Ej_inv_all = NULL;
1969:   sub_schurs->sum_S_Ej_tilda_all = NULL;
1970:   sub_schurs->is_Ej_all = NULL;
1971:   return(0);
1972: }

1974: PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1975: {
1976:   PCBDDCSubSchurs schurs_ctx;
1977:   PetscErrorCode  ierr;

1980:   PetscNew(&schurs_ctx);
1981:   schurs_ctx->n_subs = 0;
1982:   *sub_schurs = schurs_ctx;
1983:   return(0);
1984: }

1986: PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1987: {
1988:   PetscInt       i;

1992:   if (!sub_schurs) return(0);
1993:   PetscFree(sub_schurs->prefix);
1994:   MatDestroy(&sub_schurs->A);
1995:   MatDestroy(&sub_schurs->S);
1996:   ISDestroy(&sub_schurs->is_I);
1997:   ISDestroy(&sub_schurs->is_B);
1998:   ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);
1999:   ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);
2000:   MatDestroy(&sub_schurs->S_Ej_all);
2001:   MatDestroy(&sub_schurs->sum_S_Ej_all);
2002:   MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
2003:   MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);
2004:   ISDestroy(&sub_schurs->is_Ej_all);
2005:   ISDestroy(&sub_schurs->is_vertices);
2006:   ISDestroy(&sub_schurs->is_dir);
2007:   PetscBTDestroy(&sub_schurs->is_edge);
2008:   for (i=0;i<sub_schurs->n_subs;i++) {
2009:     ISDestroy(&sub_schurs->is_subs[i]);
2010:   }
2011:   if (sub_schurs->n_subs) {
2012:     PetscFree(sub_schurs->is_subs);
2013:   }
2014:   if (sub_schurs->reuse_solver) {
2015:     PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
2016:   }
2017:   PetscFree(sub_schurs->reuse_solver);
2018:   if (sub_schurs->change) {
2019:     for (i=0;i<sub_schurs->n_subs;i++) {
2020:       KSPDestroy(&sub_schurs->change[i]);
2021:       ISDestroy(&sub_schurs->change_primal_sub[i]);
2022:     }
2023:   }
2024:   PetscFree(sub_schurs->change);
2025:   PetscFree(sub_schurs->change_primal_sub);
2026:   sub_schurs->n_subs = 0;
2027:   return(0);
2028: }

2030: PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
2031: {

2035:   PCBDDCSubSchursReset(*sub_schurs);
2036:   PetscFree(*sub_schurs);
2037:   return(0);
2038: }

2040: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
2041: {
2042:   PetscInt       i,j,n;

2046:   n = 0;
2047:   for (i=-n_prev;i<0;i++) {
2048:     PetscInt start_dof = queue_tip[i];
2049:     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
2050:       PetscInt dof = adjncy[j];
2051:       if (!PetscBTLookup(touched,dof)) {
2052:         PetscBTSet(touched,dof);
2053:         queue_tip[n] = dof;
2054:         n++;
2055:       }
2056:     }
2057:   }
2058:   *n_added = n;
2059:   return(0);
2060: }