Actual source code: bddcschurs.c

petsc-3.9.3 2018-07-02
Report Typos and Errors
  1:  #include <../src/ksp/pc/impls/bddc/bddc.h>
  2:  #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
  3:  #include <../src/mat/impls/dense/seq/dense.h>
  4:  #include <petscblaslapack.h>

  6: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
  7: static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
  8: static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
  9: static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);

 11: /* if v2 is not present, correction is done in-place */
 12: PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
 13: {
 14:   PetscScalar    *array;
 15:   PetscScalar    *array2;

 19:   if (!ctx->benign_n) return(0);
 20:   if (sol && full) {
 21:     PetscInt n_I,size_schur;

 23:     /* get sizes */
 24:     MatGetSize(ctx->benign_csAIB,&size_schur,NULL);
 25:     VecGetSize(v,&n_I);
 26:     n_I = n_I - size_schur;
 27:     /* get schur sol from array */
 28:     VecGetArray(v,&array);
 29:     VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);
 30:     VecRestoreArray(v,&array);
 31:     /* apply interior sol correction */
 32:     MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);
 33:     VecResetArray(ctx->benign_dummy_schur_vec);
 34:     MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);
 35:   }
 36:   if (v2) {
 37:     PetscInt nl;

 39:     VecGetArrayRead(v,(const PetscScalar**)&array);
 40:     VecGetLocalSize(v2,&nl);
 41:     VecGetArray(v2,&array2);
 42:     PetscMemcpy(array2,array,nl*sizeof(PetscScalar));
 43:   } else {
 44:     VecGetArray(v,&array);
 45:     array2 = array;
 46:   }
 47:   if (!sol) { /* change rhs */
 48:     PetscInt n;
 49:     for (n=0;n<ctx->benign_n;n++) {
 50:       PetscScalar    sum = 0.;
 51:       const PetscInt *cols;
 52:       PetscInt       nz,i;

 54:       ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);
 55:       ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);
 56:       for (i=0;i<nz-1;i++) sum += array[cols[i]];
 57: #if defined(PETSC_USE_COMPLEX)
 58:       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
 59: #else
 60:       sum = -sum/nz;
 61: #endif
 62:       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
 63:       ctx->benign_save_vals[n] = array2[cols[nz-1]];
 64:       array2[cols[nz-1]] = sum;
 65:       ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);
 66:     }
 67:   } else {
 68:     PetscInt n;
 69:     for (n=0;n<ctx->benign_n;n++) {
 70:       PetscScalar    sum = 0.;
 71:       const PetscInt *cols;
 72:       PetscInt       nz,i;
 73:       ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);
 74:       ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);
 75:       for (i=0;i<nz-1;i++) sum += array[cols[i]];
 76: #if defined(PETSC_USE_COMPLEX)
 77:       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
 78: #else
 79:       sum = -sum/nz;
 80: #endif
 81:       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
 82:       array2[cols[nz-1]] = ctx->benign_save_vals[n];
 83:       ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);
 84:     }
 85:   }
 86:   if (v2) {
 87:     VecRestoreArrayRead(v,(const PetscScalar**)&array);
 88:     VecRestoreArray(v2,&array2);
 89:   } else {
 90:     VecRestoreArray(v,&array);
 91:   }
 92:   if (!sol && full) {
 93:     Vec      usedv;
 94:     PetscInt n_I,size_schur;

 96:     /* get sizes */
 97:     MatGetSize(ctx->benign_csAIB,&size_schur,NULL);
 98:     VecGetSize(v,&n_I);
 99:     n_I = n_I - size_schur;
100:     /* compute schur rhs correction */
101:     if (v2) {
102:       usedv = v2;
103:     } else {
104:       usedv = v;
105:     }
106:     /* apply schur rhs correction */
107:     MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);
108:     VecGetArrayRead(usedv,(const PetscScalar**)&array);
109:     VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);
110:     VecRestoreArrayRead(usedv,(const PetscScalar**)&array);
111:     MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);
112:     VecResetArray(ctx->benign_dummy_schur_vec);
113:   }
114:   return(0);
115: }

117: static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
118: {
119:   PCBDDCReuseSolvers ctx;
120:   PetscBool          copy = PETSC_FALSE;
121:   PetscErrorCode     ierr;

124:   PCShellGetContext(pc,(void **)&ctx);
125:   if (full) {
126: #if defined(PETSC_HAVE_MUMPS)
127:     MatMumpsSetIcntl(ctx->F,26,-1);
128: #endif
129: #if defined(PETSC_HAVE_MKL_PARDISO)
130:     MatMkl_PardisoSetCntl(ctx->F,70,0);
131: #endif
132:     copy = ctx->has_vertices;
133:   } else { /* interior solver */
134: #if defined(PETSC_HAVE_MUMPS)
135:     MatMumpsSetIcntl(ctx->F,26,0);
136: #endif
137: #if defined(PETSC_HAVE_MKL_PARDISO)
138:     MatMkl_PardisoSetCntl(ctx->F,70,1);
139: #endif
140:     copy = PETSC_TRUE;
141:   }
142:   /* copy rhs into factored matrix workspace */
143:   if (copy) {
144:     PetscInt    n;
145:     PetscScalar *array,*array_solver;

147:     VecGetLocalSize(rhs,&n);
148:     VecGetArrayRead(rhs,(const PetscScalar**)&array);
149:     VecGetArray(ctx->rhs,&array_solver);
150:     PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));
151:     VecRestoreArray(ctx->rhs,&array_solver);
152:     VecRestoreArrayRead(rhs,(const PetscScalar**)&array);

154:     PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);
155:     if (transpose) {
156:       MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);
157:     } else {
158:       MatSolve(ctx->F,ctx->rhs,ctx->sol);
159:     }
160:     PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);

162:     /* get back data to caller worskpace */
163:     VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);
164:     VecGetArray(sol,&array);
165:     PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));
166:     VecRestoreArray(sol,&array);
167:     VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);
168:   } else {
169:     if (ctx->benign_n) {
170:       PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);
171:       if (transpose) {
172:         MatSolveTranspose(ctx->F,ctx->rhs,sol);
173:       } else {
174:         MatSolve(ctx->F,ctx->rhs,sol);
175:       }
176:       PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);
177:     } else {
178:       if (transpose) {
179:         MatSolveTranspose(ctx->F,rhs,sol);
180:       } else {
181:         MatSolve(ctx->F,rhs,sol);
182:       }
183:     }
184:   }
185:   /* restore defaults */
186: #if defined(PETSC_HAVE_MUMPS)
187:   MatMumpsSetIcntl(ctx->F,26,-1);
188: #endif
189: #if defined(PETSC_HAVE_MKL_PARDISO)
190:   MatMkl_PardisoSetCntl(ctx->F,70,0);
191: #endif
192:   return(0);
193: }

195: static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
196: {
197:   PetscErrorCode   ierr;

200:   PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);
201:   return(0);
202: }

204: static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
205: {
206:   PetscErrorCode   ierr;

209:   PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);
210:   return(0);
211: }

213: static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
214: {
215:   PetscErrorCode   ierr;

218:   PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);
219:   return(0);
220: }

222: static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
223: {
224:   PetscErrorCode   ierr;

227:   PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);
228:   return(0);
229: }

231: static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
232: {
233:   PetscInt       i;

237:   MatDestroy(&reuse->F);
238:   VecDestroy(&reuse->sol);
239:   VecDestroy(&reuse->rhs);
240:   PCDestroy(&reuse->interior_solver);
241:   PCDestroy(&reuse->correction_solver);
242:   ISDestroy(&reuse->is_R);
243:   ISDestroy(&reuse->is_B);
244:   VecScatterDestroy(&reuse->correction_scatter_B);
245:   VecDestroy(&reuse->sol_B);
246:   VecDestroy(&reuse->rhs_B);
247:   for (i=0;i<reuse->benign_n;i++) {
248:     ISDestroy(&reuse->benign_zerodiag_subs[i]);
249:   }
250:   PetscFree(reuse->benign_zerodiag_subs);
251:   PetscFree(reuse->benign_save_vals);
252:   MatDestroy(&reuse->benign_csAIB);
253:   MatDestroy(&reuse->benign_AIIm1ones);
254:   VecDestroy(&reuse->benign_corr_work);
255:   VecDestroy(&reuse->benign_dummy_schur_vec);
256:   return(0);
257: }

259: static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
260: {
261:   Mat            B, C, D, Bd, Cd, AinvBd;
262:   KSP            ksp;
263:   PC             pc;
264:   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
265:   PetscReal      fill = 2.0;
266:   PetscInt       n_I;
267:   PetscMPIInt    size;

271:   MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);
272:   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
273:   if (reuse == MAT_REUSE_MATRIX) {
274:     PetscBool Sdense;

276:     PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);
277:     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
278:   }
279:   MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);
280:   MatSchurComplementGetKSP(M, &ksp);
281:   KSPGetPC(ksp, &pc);
282:   PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);
283:   PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);
284:   PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);
285:   PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);
286:   PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);
287:   MatGetSize(B,&n_I,NULL);
288:   if (n_I) {
289:     if (!Bdense) {
290:       MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);
291:     } else {
292:       Bd = B;
293:     }

295:     if (isLU || isILU || isCHOL) {
296:       Mat fact;
297:       KSPSetUp(ksp);
298:       PCFactorGetMatrix(pc, &fact);
299:       MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
300:       MatMatSolve(fact, Bd, AinvBd);
301:     } else {
302:       PetscBool ex = PETSC_TRUE;

304:       if (ex) {
305:         Mat Ainvd;

307:         PCComputeExplicitOperator(pc, &Ainvd);
308:         MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);
309:         MatDestroy(&Ainvd);
310:       } else {
311:         Vec         sol,rhs;
312:         PetscScalar *arrayrhs,*arraysol;
313:         PetscInt    i,nrhs,n;

315:         MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
316:         MatGetSize(Bd,&n,&nrhs);
317:         MatDenseGetArray(Bd,&arrayrhs);
318:         MatDenseGetArray(AinvBd,&arraysol);
319:         KSPGetSolution(ksp,&sol);
320:         KSPGetRhs(ksp,&rhs);
321:         for (i=0;i<nrhs;i++) {
322:           VecPlaceArray(rhs,arrayrhs+i*n);
323:           VecPlaceArray(sol,arraysol+i*n);
324:           KSPSolve(ksp,rhs,sol);
325:           VecResetArray(rhs);
326:           VecResetArray(sol);
327:         }
328:         MatDenseRestoreArray(Bd,&arrayrhs);
329:         MatDenseRestoreArray(AinvBd,&arrayrhs);
330:       }
331:     }
332:     if (!Bdense & !issym) {
333:       MatDestroy(&Bd);
334:     }

336:     if (!issym) {
337:       if (!Cdense) {
338:         MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);
339:       } else {
340:         Cd = C;
341:       }
342:       MatMatMult(Cd, AinvBd, reuse, fill, S);
343:       if (!Cdense) {
344:         MatDestroy(&Cd);
345:       }
346:     } else {
347:       MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);
348:       if (!Bdense) {
349:         MatDestroy(&Bd);
350:       }
351:     }
352:     MatDestroy(&AinvBd);
353:   }

355:   if (D) {
356:     Mat       Dd;
357:     PetscBool Ddense;

359:     PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);
360:     if (!Ddense) {
361:       MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);
362:     } else {
363:       Dd = D;
364:     }
365:     if (n_I) {
366:       MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);
367:     } else {
368:       if (reuse == MAT_INITIAL_MATRIX) {
369:         MatDuplicate(Dd,MAT_COPY_VALUES,S);
370:       } else {
371:         MatCopy(Dd,*S,SAME_NONZERO_PATTERN);
372:       }
373:     }
374:     if (!Ddense) {
375:       MatDestroy(&Dd);
376:     }
377:   } else {
378:     MatScale(*S,-1.0);
379:   }
380:   return(0);
381: }

383: PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal)
384: {
385:   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
386:   Mat                    S_all;
387:   Mat                    global_schur_subsets,work_mat,*submats;
388:   ISLocalToGlobalMapping l2gmap_subsets;
389:   IS                     is_I,is_I_layer;
390:   IS                     all_subsets,all_subsets_mult,all_subsets_n;
391:   PetscInt               *nnz,*all_local_idx_N;
392:   PetscInt               *auxnum1,*auxnum2;
393:   PetscInt               i,subset_size,max_subset_size;
394:   PetscInt               n_B,extra,local_size,global_size;
395:   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
396:   PetscScalar            *Bwork;
397:   PetscSubcomm           subcomm;
398:   PetscMPIInt            color,rank;
399:   MPI_Comm               comm_n;
400:   PetscBool              deluxe = PETSC_TRUE;
401:   PetscBool              use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
402:   PetscErrorCode         ierr;

405:   MatDestroy(&sub_schurs->A);
406:   MatDestroy(&sub_schurs->S);
407:   /* convert matrix if needed */
408:   if (Ain) {
409:     PetscBool isseqaij;
410:     PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);
411:     if (isseqaij) {
412:       PetscObjectReference((PetscObject)Ain);
413:       sub_schurs->A = Ain;
414:     } else {
415:       MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);
416:     }
417:   }

419:   PetscObjectReference((PetscObject)Sin);
420:   sub_schurs->S = Sin;
421:   if (sub_schurs->schur_explicit) {
422:     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
423:   }

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

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

430:   /* restrict work on active processes */
431:   color = 0;
432:   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
433:   MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);
434:   PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);
435:   PetscSubcommSetNumber(subcomm,2);
436:   PetscSubcommSetTypeGeneral(subcomm,color,rank);
437:   PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);
438:   PetscSubcommDestroy(&subcomm);
439:   if (!sub_schurs->n_subs) {
440:     PetscCommDestroy(&comm_n);
441:     return(0);
442:   }
443:   /* PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL); */

445:   /* get Schur complement matrices */
446:   if (!sub_schurs->schur_explicit) {
447:     Mat       tA_IB,tA_BI,tA_BB;
448:     PetscBool isseqsbaij;
449:     MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);
450:     PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);
451:     if (isseqsbaij) {
452:       MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);
453:       MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);
454:       MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);
455:     } else {
456:       PetscObjectReference((PetscObject)tA_BB);
457:       A_BB = tA_BB;
458:       PetscObjectReference((PetscObject)tA_IB);
459:       A_IB = tA_IB;
460:       PetscObjectReference((PetscObject)tA_BI);
461:       A_BI = tA_BI;
462:     }
463:   } else {
464:     A_II = NULL;
465:     A_IB = NULL;
466:     A_BI = NULL;
467:     A_BB = NULL;
468:   }
469:   S_all = NULL;

471:   /* determine interior problems */
472:   ISGetLocalSize(sub_schurs->is_I,&i);
473:   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
474:     PetscBT                touched;
475:     const PetscInt*        idx_B;
476:     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;

478:     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
479:     /* get sizes */
480:     ISGetLocalSize(sub_schurs->is_I,&n_I);
481:     ISGetLocalSize(sub_schurs->is_B,&n_B);

483:     PetscMalloc1(n_I+n_B,&local_numbering);
484:     PetscBTCreate(n_I+n_B,&touched);
485:     PetscBTMemzero(n_I+n_B,touched);

487:     /* all boundary dofs must be skipped when adding layers */
488:     ISGetIndices(sub_schurs->is_B,&idx_B);
489:     for (j=0;j<n_B;j++) {
490:       PetscBTSet(touched,idx_B[j]);
491:     }
492:     PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));
493:     ISRestoreIndices(sub_schurs->is_B,&idx_B);

495:     /* add prescribed number of layers of dofs */
496:     n_local_dofs = n_B;
497:     n_prev_added = n_B;
498:     for (layer=0;layer<nlayers;layer++) {
499:       PetscInt n_added;
500:       if (n_local_dofs == n_I+n_B) break;
501:       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);
502:       PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);
503:       n_prev_added = n_added;
504:       n_local_dofs += n_added;
505:       if (!n_added) break;
506:     }
507:     PetscBTDestroy(&touched);

509:     /* IS for I layer dofs in original numbering */
510:     ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);
511:     PetscFree(local_numbering);
512:     ISSort(is_I_layer);
513:     /* IS for I layer dofs in I numbering */
514:     if (!sub_schurs->schur_explicit) {
515:       ISLocalToGlobalMapping ItoNmap;
516:       ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);
517:       ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);
518:       ISLocalToGlobalMappingDestroy(&ItoNmap);

520:       /* II block */
521:       MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);
522:     }
523:   } else {
524:     PetscInt n_I;

526:     /* IS for I dofs in original numbering */
527:     PetscObjectReference((PetscObject)sub_schurs->is_I);
528:     is_I_layer = sub_schurs->is_I;

530:     /* IS for I dofs in I numbering (strided 1) */
531:     if (!sub_schurs->schur_explicit) {
532:       ISGetSize(sub_schurs->is_I,&n_I);
533:       ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);

535:       /* II block is the same */
536:       PetscObjectReference((PetscObject)A_II);
537:       AE_II = A_II;
538:     }
539:   }

541:   /* Get info on subset sizes and sum of all subsets sizes */
542:   max_subset_size = 0;
543:   local_size = 0;
544:   for (i=0;i<sub_schurs->n_subs;i++) {
545:     ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
546:     max_subset_size = PetscMax(subset_size,max_subset_size);
547:     local_size += subset_size;
548:   }

550:   /* Work arrays for local indices */
551:   extra = 0;
552:   ISGetLocalSize(sub_schurs->is_B,&n_B);
553:   if (sub_schurs->schur_explicit && is_I_layer) {
554:     ISGetLocalSize(is_I_layer,&extra);
555:   }
556:   PetscMalloc1(n_B+extra,&all_local_idx_N);
557:   if (extra) {
558:     const PetscInt *idxs;
559:     ISGetIndices(is_I_layer,&idxs);
560:     PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));
561:     ISRestoreIndices(is_I_layer,&idxs);
562:   }
563:   PetscMalloc1(local_size,&nnz);
564:   PetscMalloc1(sub_schurs->n_subs,&auxnum1);
565:   PetscMalloc1(sub_schurs->n_subs,&auxnum2);

567:   /* Get local indices in local numbering */
568:   local_size = 0;
569:   for (i=0;i<sub_schurs->n_subs;i++) {
570:     PetscInt j;
571:     const    PetscInt *idxs;

573:     ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
574:     ISGetIndices(sub_schurs->is_subs[i],&idxs);
575:     /* start (smallest in global ordering) and multiplicity */
576:     auxnum1[i] = idxs[0];
577:     auxnum2[i] = subset_size;
578:     /* subset indices in local numbering */
579:     PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));
580:     ISRestoreIndices(sub_schurs->is_subs[i],&idxs);
581:     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
582:     local_size += subset_size;
583:   }

585:   /* allocate extra workspace needed only for GETRI or SYTRF */
586:   use_potr = use_sytr = PETSC_FALSE;
587:   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
588:     use_potr = PETSC_TRUE;
589:   } else if (sub_schurs->is_symmetric) {
590:     use_sytr = PETSC_TRUE;
591:   }
592:   if (local_size && !use_potr) {
593:     PetscScalar  lwork,dummyscalar = 0.;
594:     PetscBLASInt dummyint = 0;

596:     B_lwork = -1;
597:     PetscBLASIntCast(local_size,&B_N);
598:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
599:     if (use_sytr) {
600:       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
601:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
602:     } else {
603:       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
604:       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
605:     }
606:     PetscFPTrapPop();
607:     PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);
608:     PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);
609:   } else {
610:     Bwork = NULL;
611:     pivots = NULL;
612:   }

614:   /* prepare parallel matrices for summing up properly schurs on subsets */
615:   ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);
616:   ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);
617:   ISDestroy(&all_subsets_n);
618:   ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);
619:   ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);
620:   ISDestroy(&all_subsets);
621:   ISDestroy(&all_subsets_mult);
622:   ISGetLocalSize(all_subsets_n,&i);
623:   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
624:   ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);
625:   MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);
626:   ISLocalToGlobalMappingDestroy(&l2gmap_subsets);
627:   MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);
628:   MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);
629:   MatSetType(global_schur_subsets,MATMPIAIJ);

631:   /* subset indices in local boundary numbering */
632:   if (!sub_schurs->is_Ej_all) {
633:     PetscInt *all_local_idx_B;

635:     PetscMalloc1(local_size,&all_local_idx_B);
636:     ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);
637:     if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D\n",subset_size,local_size);
638:     ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);
639:   }

641:   if (change) {
642:     ISLocalToGlobalMapping BtoS;
643:     IS                     change_primal_B;
644:     IS                     change_primal_all;

646:     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
647:     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
648:     PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);
649:     for (i=0;i<sub_schurs->n_subs;i++) {
650:       ISLocalToGlobalMapping NtoS;
651:       ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);
652:       ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);
653:       ISLocalToGlobalMappingDestroy(&NtoS);
654:     }
655:     ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);
656:     ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);
657:     ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);
658:     ISLocalToGlobalMappingDestroy(&BtoS);
659:     ISDestroy(&change_primal_B);
660:     PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);
661:     for (i=0;i<sub_schurs->n_subs;i++) {
662:       Mat change_sub;

664:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
665:       KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);
666:       KSPSetType(sub_schurs->change[i],KSPPREONLY);
667:       if (!sub_schurs->change_with_qr) {
668:         MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);
669:       } else {
670:         Mat change_subt;
671:         MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);
672:         MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);
673:         MatDestroy(&change_subt);
674:       }
675:       KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);
676:       MatDestroy(&change_sub);
677:       KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);
678:       KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");
679:     }
680:     ISDestroy(&change_primal_all);
681:   }

683:   /* Local matrix of all local Schur on subsets (transposed) */
684:   if (!sub_schurs->S_Ej_all) {
685:     MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);
686:     MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);
687:     MatSetType(sub_schurs->S_Ej_all,MATAIJ);
688:     MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);
689:   }

691:   /* Compute Schur complements explicitly */
692:   F = NULL;
693:   if (!sub_schurs->schur_explicit) {
694:     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
695:        it is not efficient, unless the economic version of the scaling is used */
696:     Mat         S_Ej_expl;
697:     PetscScalar *work;
698:     PetscInt    j,*dummy_idx;
699:     PetscBool   Sdense;

701:     PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);
702:     local_size = 0;
703:     for (i=0;i<sub_schurs->n_subs;i++) {
704:       IS  is_subset_B;
705:       Mat AE_EE,AE_IE,AE_EI,S_Ej;

707:       /* subsets in original and boundary numbering */
708:       ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);
709:       /* EE block */
710:       MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);
711:       /* IE block */
712:       MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);
713:       /* EI block */
714:       if (sub_schurs->is_symmetric) {
715:         MatCreateTranspose(AE_IE,&AE_EI);
716:       } else if (sub_schurs->is_hermitian) {
717:         MatCreateHermitianTranspose(AE_IE,&AE_EI);
718:       } else {
719:         MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);
720:       }
721:       ISDestroy(&is_subset_B);
722:       MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);
723:       MatDestroy(&AE_EE);
724:       MatDestroy(&AE_IE);
725:       MatDestroy(&AE_EI);
726:       if (AE_II == A_II) { /* we can reuse the same ksp */
727:         KSP ksp;
728:         MatSchurComplementGetKSP(sub_schurs->S,&ksp);
729:         MatSchurComplementSetKSP(S_Ej,ksp);
730:       } else { /* build new ksp object which inherits ksp and pc types from the original one */
731:         KSP       origksp,schurksp;
732:         PC        origpc,schurpc;
733:         KSPType   ksp_type;
734:         PetscInt  n_internal;
735:         PetscBool ispcnone;

737:         MatSchurComplementGetKSP(sub_schurs->S,&origksp);
738:         MatSchurComplementGetKSP(S_Ej,&schurksp);
739:         KSPGetType(origksp,&ksp_type);
740:         KSPSetType(schurksp,ksp_type);
741:         KSPGetPC(schurksp,&schurpc);
742:         KSPGetPC(origksp,&origpc);
743:         PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);
744:         if (!ispcnone) {
745:           PCType pc_type;
746:           PCGetType(origpc,&pc_type);
747:           PCSetType(schurpc,pc_type);
748:         } else {
749:           PCSetType(schurpc,PCLU);
750:         }
751:         ISGetSize(is_I,&n_internal);
752:         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
753:           MatSolverType solver = NULL;
754:           PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);
755:           if (solver) {
756:             PCFactorSetMatSolverType(schurpc,solver);
757:           }
758:         }
759:         KSPSetUp(schurksp);
760:       }
761:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
762:       MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);
763:       PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);
764:       PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);
765:       if (Sdense) {
766:         for (j=0;j<subset_size;j++) {
767:           dummy_idx[j]=local_size+j;
768:         }
769:         MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
770:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
771:       MatDestroy(&S_Ej);
772:       MatDestroy(&S_Ej_expl);
773:       local_size += subset_size;
774:     }
775:     PetscFree2(dummy_idx,work);
776:     /* free */
777:     ISDestroy(&is_I);
778:     MatDestroy(&AE_II);
779:     PetscFree(all_local_idx_N);
780:   } else {
781:     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
782:     Vec         Dall = NULL;
783:     IS          is_A_all,*is_p_r = NULL;
784:     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
785:     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
786:     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE;
787:     PetscBool   schur_has_vertices,factor_workaround;
788:     PetscBool   use_cholesky;

790:     /* get sizes */
791:     n_I = 0;
792:     if (is_I_layer) {
793:       ISGetLocalSize(is_I_layer,&n_I);
794:     }
795:     economic = PETSC_FALSE;
796:     ISGetLocalSize(sub_schurs->is_I,&cum);
797:     if (cum != n_I) economic = PETSC_TRUE;
798:     MatGetLocalSize(sub_schurs->A,&n,NULL);
799:     size_active_schur = local_size;

801:     /* import scaling vector (wrong formulation if we have 3D edges) */
802:     if (scaling && compute_Stilda) {
803:       const PetscScalar *array;
804:       PetscScalar       *array2;
805:       const PetscInt    *idxs;
806:       PetscInt          i;

808:       ISGetIndices(sub_schurs->is_Ej_all,&idxs);
809:       VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);
810:       VecGetArrayRead(scaling,&array);
811:       VecGetArray(Dall,&array2);
812:       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
813:       VecRestoreArray(Dall,&array2);
814:       VecRestoreArrayRead(scaling,&array);
815:       ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);
816:       deluxe = PETSC_FALSE;
817:     }

819:     /* size active schurs does not count any dirichlet or vertex dof on the interface */
820:     factor_workaround = PETSC_FALSE;
821:     schur_has_vertices = PETSC_FALSE;
822:     cum = n_I+size_active_schur;
823:     if (sub_schurs->is_dir) {
824:       const PetscInt* idxs;
825:       PetscInt        n_dir;

827:       ISGetLocalSize(sub_schurs->is_dir,&n_dir);
828:       ISGetIndices(sub_schurs->is_dir,&idxs);
829:       PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));
830:       ISRestoreIndices(sub_schurs->is_dir,&idxs);
831:       cum += n_dir;
832:       factor_workaround = PETSC_TRUE;
833:     }
834:     /* include the primal vertices in the Schur complement */
835:     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
836:       PetscInt n_v;

838:       ISGetLocalSize(sub_schurs->is_vertices,&n_v);
839:       if (n_v) {
840:         const PetscInt* idxs;

842:         ISGetIndices(sub_schurs->is_vertices,&idxs);
843:         PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));
844:         ISRestoreIndices(sub_schurs->is_vertices,&idxs);
845:         cum += n_v;
846:         factor_workaround = PETSC_TRUE;
847:         schur_has_vertices = PETSC_TRUE;
848:       }
849:     }
850:     size_schur = cum - n_I;
851:     ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);
852:     if (cum == n) {
853:       ISSetPermutation(is_A_all);
854:       MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);
855:     } else {
856:       MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);
857:     }
858:     MatSetOptionsPrefix(A,sub_schurs->prefix);
859:     MatAppendOptionsPrefix(A,"sub_schurs_");

861:     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
862:        n^2 instead of n^1.5 or something. This is a workaround */
863:     if (benign_n) {
864:       Vec                    v;
865:       ISLocalToGlobalMapping N_to_reor;
866:       IS                     is_p0,is_p0_p;
867:       PetscScalar            *cs_AIB,*AIIm1_data;
868:       PetscInt               sizeA;

870:       ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);
871:       ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);
872:       ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);
873:       ISDestroy(&is_p0);
874:       MatCreateVecs(A,&v,NULL);
875:       VecGetSize(v,&sizeA);
876:       MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);
877:       MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);
878:       MatDenseGetArray(cs_AIB_mat,&cs_AIB);
879:       MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
880:       PetscMalloc1(benign_n,&is_p_r);
881:       /* compute colsum of A_IB restricted to pressures */
882:       for (i=0;i<benign_n;i++) {
883:         Vec            benign_AIIm1_ones;
884:         PetscScalar    *array;
885:         const PetscInt *idxs;
886:         PetscInt       j,nz;

888:         ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);
889:         ISGetLocalSize(is_p_r[i],&nz);
890:         ISGetIndices(is_p_r[i],&idxs);
891:         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
892:         ISRestoreIndices(is_p_r[i],&idxs);
893:         VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);
894:         MatMult(A,benign_AIIm1_ones,v);
895:         VecDestroy(&benign_AIIm1_ones);
896:         VecGetArray(v,&array);
897:         for (j=0;j<size_schur;j++) {
898: #if defined(PETSC_USE_COMPLEX)
899:           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
900: #else
901:           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
902: #endif
903:         }
904:         VecRestoreArray(v,&array);
905:       }
906:       MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
907:       MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
908:       VecDestroy(&v);
909:       MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);
910:       MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
911:       MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
912:       MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);
913:       ISDestroy(&is_p0_p);
914:       ISLocalToGlobalMappingDestroy(&N_to_reor);
915:     }
916:     MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);
917:     MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);
918:     MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);

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

923:     /* when using the benign subspace trick, the local Schur complements are SPD */
924:     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;

926:     if (n_I) {
927:       IS is_schur;

929:       if (use_cholesky) {
930:         MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);
931:       } else {
932:         MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);
933:       }
934:       MatSetErrorIfFailure(A,PETSC_TRUE);

936:       /* subsets ordered last */
937:       ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);
938:       MatFactorSetSchurIS(F,is_schur);
939:       ISDestroy(&is_schur);

941:       /* factorization step */
942:       if (use_cholesky) {
943:         MatCholeskyFactorSymbolic(F,A,NULL,NULL);
944: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
945:         MatMumpsSetIcntl(F,19,2);
946: #endif
947:         MatCholeskyFactorNumeric(F,A,NULL);
948:         S_lower_triangular = PETSC_TRUE;
949:       } else {
950:         MatLUFactorSymbolic(F,A,NULL,NULL,NULL);
951: #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
952:         MatMumpsSetIcntl(F,19,3);
953: #endif
954:         MatLUFactorNumeric(F,A,NULL);
955:       }
956:       MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");

958: #if defined(SUB_SCHURS_DEBUG)
959:       {
960:         PetscViewer viewer;
961:         char        filename[256];
962:         Mat         S;
963:         IS          is;

965:         sprintf(filename,"sub_schurs_Schur_r%d.m",PetscGlobalRank);
966:         PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);
967:         PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
968:         PetscObjectSetName((PetscObject)A,"A");
969:         MatView(A,viewer);
970:         MatFactorCreateSchurComplement(F,&S,NULL);
971:         PetscObjectSetName((PetscObject)S,"S");
972:         MatView(S,viewer);
973:         MatDestroy(&S);
974:         ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);
975:         PetscObjectSetName((PetscObject)is,"I");
976:         ISView(is,viewer);
977:         ISDestroy(&is);
978:         ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);
979:         PetscObjectSetName((PetscObject)is,"B");
980:         ISView(is,viewer);
981:         ISDestroy(&is);
982:         PetscViewerDestroy(&viewer);
983:       }
984: #endif

986:       /* get explicit Schur Complement computed during numeric factorization */
987:       MatFactorGetSchurComplement(F,&S_all,NULL);
988:       MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);
989:       MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);

991:       /* we can reuse the solvers if we are not using the economic version */
992:       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
993:       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
994:       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
995:         reuse_solvers = factor_workaround = PETSC_FALSE;

997:       solver_S = PETSC_TRUE;

999:       /* update the Schur complement with the change of basis on the pressures */
1000:       if (benign_n) {
1001:         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1002:         Vec         v;
1003:         PetscInt    sizeA;

1005:         MatDenseGetArray(S_all,&S_data);
1006:         MatCreateVecs(A,&v,NULL);
1007:         VecGetSize(v,&sizeA);
1008: #if defined(PETSC_HAVE_MUMPS)
1009:         MatMumpsSetIcntl(F,26,0);
1010: #endif
1011: #if defined(PETSC_HAVE_MKL_PARDISO)
1012:         MatMkl_PardisoSetCntl(F,70,1);
1013: #endif
1014:         MatDenseGetArray(cs_AIB_mat,&cs_AIB);
1015:         MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);
1016:         for (i=0;i<benign_n;i++) {
1017:           Vec            benign_AIIm1_ones;
1018:           PetscScalar    *array,sum = 0.,one = 1.;
1019:           const PetscInt *idxs;
1020:           PetscInt       j,nz;
1021:           PetscBLASInt   B_k,B_n;

1023:           VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);
1024:           VecCopy(benign_AIIm1_ones,v);
1025:           MatSolve(F,v,benign_AIIm1_ones);
1026:           ISGetLocalSize(is_p_r[i],&nz);
1027:           ISGetIndices(is_p_r[i],&idxs);
1028:           /* p0 dof (eliminated) is excluded from the sum */
1029:           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1030:           ISRestoreIndices(is_p_r[i],&idxs);
1031:           MatMult(A,benign_AIIm1_ones,v);
1032:           VecGetArray(v,&array);
1033:           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1034:           /* cs_AIB already scaled by 1./nz */
1035:           B_k = 1;
1036:           PetscBLASIntCast(size_schur,&B_n);
1037:           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1038:           sum  = 1.;
1039:           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));
1040:           VecRestoreArray(v,&array);
1041:           /* set p0 entry of AIIm1_ones to zero */
1042:           ISGetIndices(is_p_r[i],&idxs);
1043:           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1044:           ISRestoreIndices(is_p_r[i],&idxs);
1045:           VecDestroy(&benign_AIIm1_ones);
1046:         }
1047:         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1048:           PetscInt k,j;
1049:           for (k=0;k<size_schur;k++) {
1050:             for (j=k;j<size_schur;j++) {
1051:               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1052:             }
1053:           }
1054:         }

1056:         /* restore defaults */
1057: #if defined(PETSC_HAVE_MUMPS)
1058:         MatMumpsSetIcntl(F,26,-1);
1059: #endif
1060: #if defined(PETSC_HAVE_MKL_PARDISO)
1061:         MatMkl_PardisoSetCntl(F,70,0);
1062: #endif
1063:         MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);
1064:         MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);
1065:         VecDestroy(&v);
1066:         MatDenseRestoreArray(S_all,&S_data);
1067:       }
1068:       if (!reuse_solvers) {
1069:         for (i=0;i<benign_n;i++) {
1070:           ISDestroy(&is_p_r[i]);
1071:         }
1072:         PetscFree(is_p_r);
1073:         MatDestroy(&cs_AIB_mat);
1074:         MatDestroy(&benign_AIIm1_ones_mat);
1075:       }
1076:     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1077:       MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);
1078:       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1079:       factor_workaround = PETSC_FALSE;
1080:       solver_S = PETSC_FALSE;
1081:     }

1083:     if (reuse_solvers) {
1084:       Mat                A_II,Afake;
1085:       Vec                vec1_B;
1086:       PCBDDCReuseSolvers msolv_ctx;
1087:       PetscInt           n_R;

1089:       if (sub_schurs->reuse_solver) {
1090:         PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1091:       } else {
1092:         PetscNew(&sub_schurs->reuse_solver);
1093:       }
1094:       msolv_ctx = sub_schurs->reuse_solver;
1095:       MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);
1096:       PetscObjectReference((PetscObject)F);
1097:       msolv_ctx->F = F;
1098:       MatCreateVecs(F,&msolv_ctx->sol,NULL);
1099:       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1100:       {
1101:         PetscScalar *array;
1102:         PetscInt    n;

1104:         VecGetLocalSize(msolv_ctx->sol,&n);
1105:         VecGetArray(msolv_ctx->sol,&array);
1106:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);
1107:         VecRestoreArray(msolv_ctx->sol,&array);
1108:       }
1109:       msolv_ctx->has_vertices = schur_has_vertices;

1111:       /* interior solver */
1112:       PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);
1113:       PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);
1114:       PCSetType(msolv_ctx->interior_solver,PCSHELL);
1115:       PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);
1116:       PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);
1117:       PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);

1119:       /* correction solver */
1120:       PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);
1121:       PCSetType(msolv_ctx->correction_solver,PCSHELL);
1122:       PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);
1123:       PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);
1124:       PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);

1126:       /* scatter and vecs for Schur complement solver */
1127:       MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);
1128:       MatCreateVecs(sub_schurs->S,&vec1_B,NULL);
1129:       if (!schur_has_vertices) {
1130:         ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);
1131:         VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);
1132:         PetscObjectReference((PetscObject)is_A_all);
1133:         msolv_ctx->is_R = is_A_all;
1134:       } else {
1135:         IS              is_B_all;
1136:         const PetscInt* idxs;
1137:         PetscInt        dual,n_v,n;

1139:         ISGetLocalSize(sub_schurs->is_vertices,&n_v);
1140:         dual = size_schur - n_v;
1141:         ISGetLocalSize(is_A_all,&n);
1142:         ISGetIndices(is_A_all,&idxs);
1143:         ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);
1144:         ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);
1145:         ISDestroy(&is_B_all);
1146:         ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);
1147:         VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);
1148:         ISDestroy(&is_B_all);
1149:         ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);
1150:         ISRestoreIndices(is_A_all,&idxs);
1151:       }
1152:       ISGetLocalSize(msolv_ctx->is_R,&n_R);
1153:       MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);
1154:       MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);
1155:       MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);
1156:       PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);
1157:       MatDestroy(&Afake);
1158:       VecDestroy(&vec1_B);

1160:       /* communicate benign info to solver context */
1161:       if (benign_n) {
1162:         PetscScalar *array;

1164:         msolv_ctx->benign_n = benign_n;
1165:         msolv_ctx->benign_zerodiag_subs = is_p_r;
1166:         PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);
1167:         msolv_ctx->benign_csAIB = cs_AIB_mat;
1168:         MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);
1169:         VecGetArray(msolv_ctx->benign_corr_work,&array);
1170:         VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);
1171:         VecRestoreArray(msolv_ctx->benign_corr_work,&array);
1172:         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1173:       }
1174:     } else {
1175:       if (sub_schurs->reuse_solver) {
1176:         PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1177:       }
1178:       PetscFree(sub_schurs->reuse_solver);
1179:     }
1180:     MatDestroy(&A);
1181:     ISDestroy(&is_A_all);

1183:     /* Work arrays */
1184:     PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);

1186:     /* matrices for deluxe scaling and adaptive selection */
1187:     if (compute_Stilda) {
1188:       if (!sub_schurs->sum_S_Ej_tilda_all) {
1189:         MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);
1190:       }
1191:       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1192:         MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);
1193:       }
1194:     }

1196:     /* S_Ej_all */
1197:     cum = cum2 = 0;
1198:     MatDenseGetArray(S_all,&S_data);
1199:     for (i=0;i<sub_schurs->n_subs;i++) {
1200:       PetscInt j;

1202:       /* get S_E */
1203:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1204:       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1205:         PetscInt k;
1206:         for (k=0;k<subset_size;k++) {
1207:           for (j=k;j<subset_size;j++) {
1208:             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1209:             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1210:           }
1211:         }
1212:       } else { /* just copy to workspace */
1213:         PetscInt k;
1214:         for (k=0;k<subset_size;k++) {
1215:           for (j=0;j<subset_size;j++) {
1216:             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1217:           }
1218:         }
1219:       }
1220:       /* insert S_E values */
1221:       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1222:       if (sub_schurs->change) {
1223:         Mat change_sub,SEj,T;

1225:         /* change basis */
1226:         KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1227:         MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1228:         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1229:           Mat T2;
1230:           MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1231:           MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1232:           MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1233:           MatDestroy(&T2);
1234:         } else {
1235:           MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1236:         }
1237:         MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1238:         MatDestroy(&T);
1239:         MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);
1240:         MatDestroy(&SEj);
1241:       }
1242:       if (deluxe) {
1243:         MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1244:         /* if adaptivity is requested, invert S_E blocks */
1245:         if (compute_Stilda) {
1246:           PetscInt k;

1248:           PetscBLASIntCast(subset_size,&B_N);
1249:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1250:           if (use_potr) {
1251:             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1252:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1253:             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1254:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1255:             for (k=0;k<subset_size;k++) {
1256:               for (j=k;j<subset_size;j++) {
1257:                 work[j*subset_size+k] = work[k*subset_size+j];
1258:               }
1259:             }
1260:           } else if (use_sytr) {
1261:             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1262:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1263:             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,work,&B_N,pivots,Bwork,&B_ierr));
1264:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1265:             for (k=0;k<subset_size;k++) {
1266:               for (j=k;j<subset_size;j++) {
1267:                 work[j*subset_size+k] = work[k*subset_size+j];
1268:               }
1269:             }
1270:           } else {
1271:             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1272:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1273:             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1274:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1275:           }
1276:           PetscFPTrapPop();
1277:           MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1278:         }
1279:       } else if (compute_Stilda) { /* not using deluxe */
1280:         Mat         SEj;
1281:         Vec         D;
1282:         PetscScalar *array;

1284:         MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1285:         VecGetArray(Dall,&array);
1286:         VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);
1287:         VecRestoreArray(Dall,&array);
1288:         VecShift(D,-1.);
1289:         MatDiagonalScale(SEj,D,D);
1290:         MatDestroy(&SEj);
1291:         VecDestroy(&D);
1292:         MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1293:       }
1294:       cum += subset_size;
1295:       cum2 += subset_size*(size_schur + 1);
1296:     }
1297:     MatDenseRestoreArray(S_all,&S_data);

1299:     if (solver_S) {
1300:       MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1301:     }

1303:     schur_factor = NULL;
1304:     if (compute_Stilda && size_active_schur) {

1306:       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1307:         PetscInt j;
1308:         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1309:         MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);
1310:       } else {
1311:         Mat S_all_inv=NULL;
1312:         if (solver_S) {
1313:           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1314:              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 */
1315:           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1316:             PetscScalar *data;
1317:             PetscInt     nd = 0;

1319:             if (!use_potr) {
1320:               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1321:             }
1322:             MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1323:             MatDenseGetArray(S_all_inv,&data);
1324:             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1325:               ISGetLocalSize(sub_schurs->is_dir,&nd);
1326:             }

1328:             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1329:             if (schur_has_vertices) {
1330:               Mat          M;
1331:               PetscScalar *tdata;
1332:               PetscInt     nv = 0, news;

1334:               ISGetLocalSize(sub_schurs->is_vertices,&nv);
1335:               news = size_active_schur + nv;
1336:               PetscCalloc1(news*news,&tdata);
1337:               for (i=0;i<size_active_schur;i++) {
1338:                 PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));
1339:                 PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));
1340:               }
1341:               for (i=0;i<nv;i++) {
1342:                 PetscInt k = i+size_active_schur;
1343:                 PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));
1344:               }

1346:               MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);
1347:               MatSetOption(M,MAT_SPD,PETSC_TRUE);
1348:               MatCholeskyFactor(M,NULL,NULL);
1349:               /* save the factors */
1350:               cum = 0;
1351:               PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);
1352:               for (i=0;i<size_active_schur;i++) {
1353:                 PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));
1354:                 cum += size_active_schur - i;
1355:               }
1356:               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1357:               MatSeqDenseInvertFactors_Private(M);
1358:               /* move back just the active dofs to the Schur complement */
1359:               for (i=0;i<size_active_schur;i++) {
1360:                 PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));
1361:               }
1362:               PetscFree(tdata);
1363:               MatDestroy(&M);
1364:             } else { /* we can factorize and invert just the activedofs */
1365:               Mat         M;
1366:               PetscScalar *aux;

1368:               PetscMalloc1(nd,&aux);
1369:               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1370:               MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);
1371:               MatSeqDenseSetLDA(M,size_schur);
1372:               MatSetOption(M,MAT_SPD,PETSC_TRUE);
1373:               MatCholeskyFactor(M,NULL,NULL);
1374:               MatSeqDenseInvertFactors_Private(M);
1375:               MatDestroy(&M);
1376:               MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);
1377:               MatZeroEntries(M);
1378:               MatDestroy(&M);
1379:               MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);
1380:               MatSeqDenseSetLDA(M,size_schur);
1381:               MatZeroEntries(M);
1382:               MatDestroy(&M);
1383:               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
1384:               PetscFree(aux);
1385:             }
1386:             MatDenseRestoreArray(S_all_inv,&data);
1387:           } else { /* use MatFactor calls to invert S */
1388:             MatFactorInvertSchurComplement(F);
1389:             MatFactorGetSchurComplement(F,&S_all_inv,NULL);
1390:           }
1391:         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1392:           PetscObjectReference((PetscObject)S_all);
1393:           S_all_inv = S_all;
1394:           MatDenseGetArray(S_all_inv,&S_data);
1395:           PetscBLASIntCast(size_schur,&B_N);
1396:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1397:           if (use_potr) {
1398:             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1399:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1400:             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1401:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1402:           } else if (use_sytr) {
1403:             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1404:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1405:             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1406:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1407:           } else {
1408:             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1409:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1410:             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1411:             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1412:           }
1413:           PetscFPTrapPop();
1414:           MatDenseRestoreArray(S_all_inv,&S_data);
1415:         }
1416:         /* S_Ej_tilda_all */
1417:         cum = cum2 = 0;
1418:         MatDenseGetArray(S_all_inv,&S_data);
1419:         for (i=0;i<sub_schurs->n_subs;i++) {
1420:           PetscInt j;

1422:           ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1423:           /* get (St^-1)_E */
1424:           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1425:              will be properly accessed later during adaptive selection */
1426:           if (S_lower_triangular) {
1427:             PetscInt k;
1428:             if (sub_schurs->change) {
1429:               for (k=0;k<subset_size;k++) {
1430:                 for (j=k;j<subset_size;j++) {
1431:                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1432:                   work[j*subset_size+k] = work[k*subset_size+j];
1433:                 }
1434:               }
1435:             } else {
1436:               for (k=0;k<subset_size;k++) {
1437:                 for (j=k;j<subset_size;j++) {
1438:                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1439:                 }
1440:               }
1441:             }
1442:           } else {
1443:             PetscInt k;
1444:             for (k=0;k<subset_size;k++) {
1445:               for (j=0;j<subset_size;j++) {
1446:                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1447:               }
1448:             }
1449:           }
1450:           if (sub_schurs->change) {
1451:             Mat change_sub,SEj,T;

1453:             /* change basis */
1454:             KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);
1455:             MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);
1456:             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1457:               Mat T2;
1458:               MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);
1459:               MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1460:               MatDestroy(&T2);
1461:               MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);
1462:             } else {
1463:               MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);
1464:             }
1465:             MatCopy(T,SEj,SAME_NONZERO_PATTERN);
1466:             MatDestroy(&T);
1467:             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1468:             MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);
1469:             MatDestroy(&SEj);
1470:           }
1471:           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1472:           MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);
1473:           cum += subset_size;
1474:           cum2 += subset_size*(size_schur + 1);
1475:         }
1476:         MatDenseRestoreArray(S_all_inv,&S_data);
1477:         if (solver_S) {
1478:           if (schur_has_vertices) {
1479:             MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);
1480:           } else {
1481:             MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);
1482:           }
1483:         }
1484:         MatDestroy(&S_all_inv);
1485:       }

1487:       /* move back factors if needed */
1488:       if (schur_has_vertices) {
1489:         Mat      S_tmp;
1490:         PetscInt nd = 0;

1492:         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1493:         MatFactorGetSchurComplement(F,&S_tmp,NULL);
1494:         if (use_potr) {
1495:           PetscScalar *data;

1497:           MatDenseGetArray(S_tmp,&data);
1498:           PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));

1500:           if (S_lower_triangular) {
1501:             cum = 0;
1502:             for (i=0;i<size_active_schur;i++) {
1503:               PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));
1504:               cum += size_active_schur-i;
1505:             }
1506:           } else {
1507:             PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));
1508:           }
1509:           if (sub_schurs->is_dir) {
1510:             ISGetLocalSize(sub_schurs->is_dir,&nd);
1511:             for (i=0;i<nd;i++) {
1512:               data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1513:             }
1514:           }
1515:           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1516:              set the diagonal entry of the Schur factor to a very large value */
1517:           for (i=size_active_schur+nd;i<size_schur;i++) {
1518:             data[i*(size_schur+1)] = infty;
1519:           }
1520:           MatDenseRestoreArray(S_tmp,&data);
1521:         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1522:         MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);
1523:       }
1524:     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1525:       PetscScalar *data;
1526:       PetscInt    nd = 0;

1528:       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1529:         ISGetLocalSize(sub_schurs->is_dir,&nd);
1530:       }
1531:       MatFactorGetSchurComplement(F,&S_all,NULL);
1532:       MatDenseGetArray(S_all,&data);
1533:       for (i=0;i<size_active_schur;i++) {
1534:         PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));
1535:       }
1536:       for (i=size_active_schur+nd;i<size_schur;i++) {
1537:         PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));
1538:         data[i*(size_schur+1)] = infty;
1539:       }
1540:       MatDenseRestoreArray(S_all,&data);
1541:       MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);
1542:     }
1543:     PetscFree2(dummy_idx,work);
1544:     PetscFree(schur_factor);
1545:     VecDestroy(&Dall);
1546:   }
1547:   PetscFree(nnz);
1548:   MatDestroy(&F);
1549:   ISDestroy(&is_I_layer);
1550:   MatDestroy(&S_all);
1551:   MatDestroy(&A_BB);
1552:   MatDestroy(&A_IB);
1553:   MatDestroy(&A_BI);
1554:   MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1555:   MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);
1556:   if (compute_Stilda) {
1557:     MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1558:     MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);
1559:     if (deluxe) {
1560:       MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1561:       MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);
1562:     }
1563:   }
1564:   /* Global matrix of all assembled Schur on subsets */
1565:   MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);
1566:   MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);
1567:   MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);

1569:   /* Get local part of (\sum_j S_Ej) */
1570:   MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);
1571:   if (!sub_schurs->sum_S_Ej_all) {
1572:     MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);
1573:   } else {
1574:     MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);
1575:   }

1577:   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1578:   if (compute_Stilda) {
1579:     MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);
1580:     MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);
1581:     MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);
1582:     MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);
1583:     if (deluxe) {
1584:       MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);
1585:       MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);
1586:       MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);
1587:       MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);
1588:     } else {
1589:       PetscScalar *array;
1590:       PetscInt    cum;

1592:       MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1593:       cum = 0;
1594:       for (i=0;i<sub_schurs->n_subs;i++) {
1595:         ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1596:         PetscBLASIntCast(subset_size,&B_N);
1597:         PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1598:         if (use_potr) {
1599:           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1600:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1601:           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1602:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1603:         } else if (use_sytr) {
1604:           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1605:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1606:           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1607:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1608:         } else {
1609:           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1610:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1611:           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1612:           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1613:         }
1614:         PetscFPTrapPop();
1615:         cum += subset_size*subset_size;
1616:       }
1617:       MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);
1618:       PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);
1619:       MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1620:       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1621:     }
1622:   }
1623:   MatDestroySubMatrices(1,&submats);
1624: #if defined(SUB_SCHURS_DEBUG)
1625:   {
1626:     PetscViewer viewer;
1627:     char        filename[256];
1628:     PetscInt    cum;

1630:     sprintf(filename,"sub_schurs_mats_r%d.m",PetscGlobalRank);
1631:     PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);
1632:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
1633:     if (sub_schurs->S_Ej_all) {
1634:       PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");
1635:       MatView(sub_schurs->S_Ej_all,viewer);
1636:     }
1637:     if (sub_schurs->sum_S_Ej_all) {
1638:       PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");
1639:       MatView(sub_schurs->sum_S_Ej_all,viewer);
1640:     }
1641:     if (sub_schurs->sum_S_Ej_inv_all) {
1642:       PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");
1643:       MatView(sub_schurs->sum_S_Ej_inv_all,viewer);
1644:     }
1645:     if (sub_schurs->sum_S_Ej_tilda_all) {
1646:       PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");
1647:       MatView(sub_schurs->sum_S_Ej_tilda_all,viewer);
1648:     }
1649:     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
1650:       IS   is;
1651:       char name[16];

1653:       sprintf(name,"IE%d",i);
1654:       ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);
1655:       ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);
1656:       PetscObjectSetName((PetscObject)is,name);
1657:       ISView(is,viewer);
1658:       ISDestroy(&is);
1659:       cum += subset_size;
1660:     }
1661:     PetscViewerDestroy(&viewer);
1662:   }
1663: #endif

1665:   /* free workspace */
1666:   PetscFree(submats);
1667:   PetscFree2(Bwork,pivots);
1668:   MatDestroy(&global_schur_subsets);
1669:   MatDestroy(&work_mat);
1670:   ISDestroy(&all_subsets_n);
1671:   PetscCommDestroy(&comm_n);
1672:   return(0);
1673: }

1675: PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1676: {
1677:   IS              *faces,*edges,*all_cc,vertices;
1678:   PetscInt        i,n_faces,n_edges,n_all_cc;
1679:   PetscBool       is_sorted,ispetsc;
1680:   PetscErrorCode  ierr;

1683:   ISSorted(is_I,&is_sorted);
1684:   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1685:   ISSorted(is_B,&is_sorted);
1686:   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");

1688:   /* reset any previous data */
1689:   PCBDDCSubSchursReset(sub_schurs);

1691:   /* get index sets for faces and edges (already sorted by global ordering) */
1692:   PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1693:   n_all_cc = n_faces+n_edges;
1694:   PetscBTCreate(n_all_cc,&sub_schurs->is_edge);
1695:   PetscMalloc1(n_all_cc,&all_cc);
1696:   for (i=0;i<n_faces;i++) {
1697:     if (copycc) {
1698:       ISDuplicate(faces[i],&all_cc[i]);
1699:     } else {
1700:       PetscObjectReference((PetscObject)faces[i]);
1701:       all_cc[i] = faces[i];
1702:     }
1703:   }
1704:   for (i=0;i<n_edges;i++) {
1705:     if (copycc) {
1706:       ISDuplicate(edges[i],&all_cc[n_faces+i]);
1707:     } else {
1708:       PetscObjectReference((PetscObject)edges[i]);
1709:       all_cc[n_faces+i] = edges[i];
1710:     }
1711:     PetscBTSet(sub_schurs->is_edge,n_faces+i);
1712:   }
1713:   PetscObjectReference((PetscObject)vertices);
1714:   sub_schurs->is_vertices = vertices;
1715:   PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);
1716:   sub_schurs->is_dir = NULL;
1717:   PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);

1719:   /* Determine if MatFactor can be used */
1720:   PetscStrallocpy(prefix,&sub_schurs->prefix);
1721: #if defined(PETSC_HAVE_MUMPS)
1722:   PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);
1723: #elif defined(PETSC_HAVE_MKL_PARDISO)
1724:   PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);
1725: #else
1726:   PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);
1727: #endif
1728: #if defined(PETSC_USE_COMPLEX)
1729:   sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1730: #else
1731:   sub_schurs->is_hermitian = PETSC_TRUE;
1732: #endif
1733:   sub_schurs->is_posdef    = PETSC_TRUE;
1734:   sub_schurs->is_symmetric = PETSC_TRUE;
1735:   PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");
1736:   PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,64,NULL);
1737:   PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);
1738:   PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);
1739:   PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);
1740:   PetscOptionsEnd();
1741:   PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);
1742:   sub_schurs->schur_explicit = (PetscBool)!ispetsc;

1744:   /* for reals, symmetric and hermitian are synonims */
1745: #if !defined(PETSC_USE_COMPLEX)
1746:   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1747:   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1748: #endif

1750:   PetscObjectReference((PetscObject)is_I);
1751:   sub_schurs->is_I = is_I;
1752:   PetscObjectReference((PetscObject)is_B);
1753:   sub_schurs->is_B = is_B;
1754:   PetscObjectReference((PetscObject)graph->l2gmap);
1755:   sub_schurs->l2gmap = graph->l2gmap;
1756:   PetscObjectReference((PetscObject)BtoNmap);
1757:   sub_schurs->BtoNmap = BtoNmap;
1758:   sub_schurs->n_subs = n_all_cc;
1759:   sub_schurs->is_subs = all_cc;
1760:   sub_schurs->S_Ej_all = NULL;
1761:   sub_schurs->sum_S_Ej_all = NULL;
1762:   sub_schurs->sum_S_Ej_inv_all = NULL;
1763:   sub_schurs->sum_S_Ej_tilda_all = NULL;
1764:   sub_schurs->is_Ej_all = NULL;
1765:   return(0);
1766: }

1768: PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1769: {
1770:   PCBDDCSubSchurs schurs_ctx;
1771:   PetscErrorCode  ierr;

1774:   PetscNew(&schurs_ctx);
1775:   schurs_ctx->n_subs = 0;
1776:   *sub_schurs = schurs_ctx;
1777:   return(0);
1778: }

1780: PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1781: {
1782:   PetscInt       i;

1786:   if (!sub_schurs) return(0);
1787:   PetscFree(sub_schurs->prefix);
1788:   MatDestroy(&sub_schurs->A);
1789:   MatDestroy(&sub_schurs->S);
1790:   ISDestroy(&sub_schurs->is_I);
1791:   ISDestroy(&sub_schurs->is_B);
1792:   ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);
1793:   ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);
1794:   MatDestroy(&sub_schurs->S_Ej_all);
1795:   MatDestroy(&sub_schurs->sum_S_Ej_all);
1796:   MatDestroy(&sub_schurs->sum_S_Ej_inv_all);
1797:   MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);
1798:   ISDestroy(&sub_schurs->is_Ej_all);
1799:   ISDestroy(&sub_schurs->is_vertices);
1800:   ISDestroy(&sub_schurs->is_dir);
1801:   PetscBTDestroy(&sub_schurs->is_edge);
1802:   for (i=0;i<sub_schurs->n_subs;i++) {
1803:     ISDestroy(&sub_schurs->is_subs[i]);
1804:   }
1805:   if (sub_schurs->n_subs) {
1806:     PetscFree(sub_schurs->is_subs);
1807:   }
1808:   if (sub_schurs->reuse_solver) {
1809:     PCBDDCReuseSolversReset(sub_schurs->reuse_solver);
1810:   }
1811:   PetscFree(sub_schurs->reuse_solver);
1812:   if (sub_schurs->change) {
1813:     for (i=0;i<sub_schurs->n_subs;i++) {
1814:       KSPDestroy(&sub_schurs->change[i]);
1815:       ISDestroy(&sub_schurs->change_primal_sub[i]);
1816:     }
1817:   }
1818:   PetscFree(sub_schurs->change);
1819:   PetscFree(sub_schurs->change_primal_sub);
1820:   sub_schurs->n_subs = 0;
1821:   return(0);
1822: }

1824: PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1825: {

1829:   PCBDDCSubSchursReset(*sub_schurs);
1830:   PetscFree(*sub_schurs);
1831:   return(0);
1832: }

1834: PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1835: {
1836:   PetscInt       i,j,n;

1840:   n = 0;
1841:   for (i=-n_prev;i<0;i++) {
1842:     PetscInt start_dof = queue_tip[i];
1843:     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1844:       PetscInt dof = adjncy[j];
1845:       if (!PetscBTLookup(touched,dof)) {
1846:         PetscBTSet(touched,dof);
1847:         queue_tip[n] = dof;
1848:         n++;
1849:       }
1850:     }
1851:   }
1852:   *n_added = n;
1853:   return(0);
1854: }