xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision eb595f79196f12f6b379f375fb4dac46f4df4e83)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
308122e43SStefano Zampini #include <petscblaslapack.h>
434a97f8cSStefano Zampini 
53202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
65ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
7d5574798SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorDestroy(PC);
8d5574798SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC,Vec,Vec);
9d5574798SStefano Zampini 
10d5574798SStefano Zampini #undef __FUNCT__
11d5574798SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
12d5574798SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
13d5574798SStefano Zampini {
14d5574798SStefano Zampini   PCBDDCMumpsInterior ctx;
15d5574798SStefano Zampini   PetscScalar         *array,*array_mumps;
16d5574798SStefano Zampini   PetscInt            ival;
17d5574798SStefano Zampini   PetscErrorCode      ierr;
18d5574798SStefano Zampini 
19d5574798SStefano Zampini   PetscFunctionBegin;
20d5574798SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
21d5574798SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
22d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
23d5574798SStefano Zampini   /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */
24d5574798SStefano Zampini   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
25d5574798SStefano Zampini   ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
26d5574798SStefano Zampini   ierr = PetscMemcpy(array_mumps,array,ctx->n*sizeof(PetscScalar));CHKERRQ(ierr);
27d5574798SStefano Zampini   ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
28d5574798SStefano Zampini   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
29d5574798SStefano Zampini 
30d5574798SStefano Zampini   ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
31d5574798SStefano Zampini 
32d5574798SStefano Zampini   /* get back data to caller worskpace */
33d5574798SStefano Zampini   ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
34d5574798SStefano Zampini   ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
35d5574798SStefano Zampini   ierr = PetscMemcpy(array,array_mumps,ctx->n*sizeof(PetscScalar));CHKERRQ(ierr);
36d5574798SStefano Zampini   ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
37d5574798SStefano Zampini   ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
38d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
39d5574798SStefano Zampini   PetscFunctionReturn(0);
40d5574798SStefano Zampini }
41d5574798SStefano Zampini 
42d5574798SStefano Zampini #undef __FUNCT__
43d5574798SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorDestroy"
44d5574798SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorDestroy(PC pc)
45d5574798SStefano Zampini {
46d5574798SStefano Zampini   PCBDDCMumpsInterior ctx;
47d5574798SStefano Zampini   PetscErrorCode ierr;
48d5574798SStefano Zampini 
49d5574798SStefano Zampini   PetscFunctionBegin;
50d5574798SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
51d5574798SStefano Zampini   ierr = MatDestroy(&ctx->F);CHKERRQ(ierr);
52d5574798SStefano Zampini   ierr = VecDestroy(&ctx->sol);CHKERRQ(ierr);
53d5574798SStefano Zampini   ierr = VecDestroy(&ctx->rhs);CHKERRQ(ierr);
54d5574798SStefano Zampini   ierr = PetscFree(ctx);CHKERRQ(ierr);
55d5574798SStefano Zampini   PetscFunctionReturn(0);
56d5574798SStefano Zampini }
573202ece2SStefano Zampini 
583202ece2SStefano Zampini #undef __FUNCT__
593202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
605ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
613202ece2SStefano Zampini {
623202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
633202ece2SStefano Zampini   KSP            ksp;
643202ece2SStefano Zampini   PC             pc;
653202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
663202ece2SStefano Zampini   PetscReal      fill = 2.0;
67f11841e3SStefano Zampini   PetscInt       n_I;
683202ece2SStefano Zampini   PetscMPIInt    size;
693202ece2SStefano Zampini   PetscErrorCode ierr;
703202ece2SStefano Zampini 
713202ece2SStefano Zampini   PetscFunctionBegin;
723202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
733202ece2SStefano Zampini   if (size != 1) {
743202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
753202ece2SStefano Zampini   }
76f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
77f11841e3SStefano Zampini     PetscBool Sdense;
78f11841e3SStefano Zampini 
79f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
80f11841e3SStefano Zampini     if (!Sdense) {
81f11841e3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
82f11841e3SStefano Zampini     }
83f11841e3SStefano Zampini   }
843202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
853202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
863202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
873202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
883202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
893202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
903202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
913202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
92f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
93f11841e3SStefano Zampini   if (n_I) {
943202ece2SStefano Zampini     if (!Bdense) {
953202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
963202ece2SStefano Zampini     } else {
973202ece2SStefano Zampini       Bd = B;
983202ece2SStefano Zampini     }
993202ece2SStefano Zampini 
1003202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1013202ece2SStefano Zampini       Mat fact;
1023202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1033202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1043202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1053202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1063202ece2SStefano Zampini     } else {
10707b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
10807b1e237SStefano Zampini 
10907b1e237SStefano Zampini       if (ex) {
1103202ece2SStefano Zampini         Mat Ainvd;
1113202ece2SStefano Zampini 
1123202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1133202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1143202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
11507b1e237SStefano Zampini       } else {
11607b1e237SStefano Zampini         Vec         sol,rhs;
11707b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
11807b1e237SStefano Zampini         PetscInt    i,nrhs,n;
11907b1e237SStefano Zampini 
12007b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
12107b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
12207b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
12307b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
12407b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
12507b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
12607b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
12707b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
12807b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
12907b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
13007b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
13107b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
13207b1e237SStefano Zampini         }
13307b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
13407b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
13507b1e237SStefano Zampini       }
1363202ece2SStefano Zampini     }
1375ec10c6aSStefano Zampini     if (!Bdense & !issym) {
1383202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
1393202ece2SStefano Zampini     }
1405ec10c6aSStefano Zampini 
1415ec10c6aSStefano Zampini     if (!issym) {
1423202ece2SStefano Zampini       if (!Cdense) {
1433202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
1443202ece2SStefano Zampini       } else {
1453202ece2SStefano Zampini         Cd = C;
1463202ece2SStefano Zampini       }
1475ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
1483202ece2SStefano Zampini       if (!Cdense) {
1493202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
1503202ece2SStefano Zampini       }
1515ec10c6aSStefano Zampini     } else {
1525ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
1535ec10c6aSStefano Zampini       if (!Bdense) {
1545ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
1555ec10c6aSStefano Zampini       }
1565ec10c6aSStefano Zampini     }
1575ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
158f11841e3SStefano Zampini   }
1593202ece2SStefano Zampini 
1603202ece2SStefano Zampini   if (D) {
1613202ece2SStefano Zampini     Mat       Dd;
1623202ece2SStefano Zampini     PetscBool Ddense;
1633202ece2SStefano Zampini 
1643202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
1653202ece2SStefano Zampini     if (!Ddense) {
1663202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
1673202ece2SStefano Zampini     } else {
1683202ece2SStefano Zampini       Dd = D;
1693202ece2SStefano Zampini     }
170f11841e3SStefano Zampini     if (n_I) {
1713202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
172f11841e3SStefano Zampini     } else {
173f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
174f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
175f11841e3SStefano Zampini       } else {
176f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
177f11841e3SStefano Zampini       }
178f11841e3SStefano Zampini     }
1793202ece2SStefano Zampini     if (!Ddense) {
1803202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
1813202ece2SStefano Zampini     }
1823202ece2SStefano Zampini   } else {
1833202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
1843202ece2SStefano Zampini   }
1853202ece2SStefano Zampini   PetscFunctionReturn(0);
1863202ece2SStefano Zampini }
18734a97f8cSStefano Zampini 
18834a97f8cSStefano Zampini #undef __FUNCT__
1891580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
190c3dc5de7SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool use_edges, PetscBool use_faces)
191b1b3d7a2SStefano Zampini {
1925a95e1ceSStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB,AE_II;
193d2627357SStefano Zampini   Mat                    S_all,S_all_inv;
194d2627357SStefano Zampini   Mat                    global_schur_subsets,work_mat;
19508122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
1965db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
1975a95e1ceSStefano Zampini   IS                     is_I,temp_is;
198d648f858SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
199*eb595f79SStefano Zampini   PetscInt               *auxnum1,*auxnum2,*all_local_idx_G_rep;
2005a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
201883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
20208122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
20306a4b1faSStefano Zampini   PetscScalar            *Bwork;
2045a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2055a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2065a95e1ceSStefano Zampini   MPI_Comm               comm_n;
207b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
208b1b3d7a2SStefano Zampini 
209b1b3d7a2SStefano Zampini   PetscFunctionBegin;
2105a95e1ceSStefano Zampini   /* preliminary checks */
2115a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
2125a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
2135a95e1ceSStefano Zampini   }
2145a95e1ceSStefano Zampini   /* determine if we are dealing with hermitian positive definite problems */
2155a95e1ceSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2165a95e1ceSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
2175a95e1ceSStefano Zampini   if (sub_schurs->A) {
2185a95e1ceSStefano Zampini     PetscInt lsize;
2195a95e1ceSStefano Zampini 
2205a95e1ceSStefano Zampini     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
2215a95e1ceSStefano Zampini     if (lsize) {
2225a95e1ceSStefano Zampini       ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
2235a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
2245a95e1ceSStefano Zampini         PetscScalar val;
2255a95e1ceSStefano Zampini         Vec         vec1,vec2;
2265a95e1ceSStefano Zampini 
2275a95e1ceSStefano Zampini         ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
2285a95e1ceSStefano Zampini         ierr = VecSetRandom(vec1,NULL);
2295a95e1ceSStefano Zampini         ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
2305a95e1ceSStefano Zampini         ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
2315a95e1ceSStefano Zampini         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
2325a95e1ceSStefano Zampini         if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
2335a95e1ceSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
2345a95e1ceSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
2355a95e1ceSStefano Zampini       }
2365a95e1ceSStefano Zampini     } else {
2375a95e1ceSStefano Zampini       sub_schurs->is_hermitian = PETSC_TRUE;
2385a95e1ceSStefano Zampini       sub_schurs->is_posdef = PETSC_TRUE;
2395a95e1ceSStefano Zampini     }
2405a95e1ceSStefano Zampini     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
2415a95e1ceSStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported");
2425a95e1ceSStefano Zampini     }
2435a95e1ceSStefano Zampini   }
2445a95e1ceSStefano Zampini   /* restrict work on active processes */
2455a95e1ceSStefano Zampini   color = 0;
2465a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
2475a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
2485a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
2495a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
2505a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
2515a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
2525a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
2535a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
2545a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
2555a95e1ceSStefano Zampini     PetscFunctionReturn(0);
2565a95e1ceSStefano Zampini   }
2575a95e1ceSStefano Zampini 
258b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
259883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
260f11841e3SStefano Zampini     PetscBool isseqaij;
261b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
262f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
263f11841e3SStefano Zampini     if (!isseqaij) {
264f11841e3SStefano Zampini       ierr = MatConvert(A_BB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BB);CHKERRQ(ierr);
265f11841e3SStefano Zampini       ierr = MatConvert(A_IB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_IB);CHKERRQ(ierr);
266f11841e3SStefano Zampini       ierr = MatConvert(A_BI,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BI);CHKERRQ(ierr);
267f11841e3SStefano Zampini     }
268a58a30b4SStefano Zampini   } else {
2695a95e1ceSStefano Zampini     A_II = NULL;
2705a95e1ceSStefano Zampini     A_IB = NULL;
2715a95e1ceSStefano Zampini     A_BI = NULL;
2725a95e1ceSStefano Zampini     A_BB = NULL;
273b1b3d7a2SStefano Zampini   }
2745a95e1ceSStefano Zampini   S_all = NULL;
2755a95e1ceSStefano Zampini   S_all_inv = NULL;
2765a95e1ceSStefano Zampini   S_Ej_tilda_all = NULL;
2775a95e1ceSStefano Zampini   S_Ej_inv_all = NULL;
278b1b3d7a2SStefano Zampini 
279b1b3d7a2SStefano Zampini   /* determine interior problems */
280b96c3477SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
2813dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
2823dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
283b1b3d7a2SStefano Zampini     PetscBT                touched;
284b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
285b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
286b1b3d7a2SStefano Zampini 
2873dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
2883dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
2893dc780c3SStefano Zampini     }
290b1b3d7a2SStefano Zampini     /* get sizes */
291b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
292b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
293b1b3d7a2SStefano Zampini 
294b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
295b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
296b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
297b1b3d7a2SStefano Zampini 
298b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
299b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
300b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
301b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
302b1b3d7a2SStefano Zampini     }
303b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
304b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
305b1b3d7a2SStefano Zampini 
306b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
307b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
308b1b3d7a2SStefano Zampini     n_prev_added = n_B;
309b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
310b1b3d7a2SStefano Zampini       PetscInt n_added;
311b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
312b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
313b1b3d7a2SStefano Zampini         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);
314b1b3d7a2SStefano Zampini       }
315b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
316b1b3d7a2SStefano Zampini       n_prev_added = n_added;
317b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
318b1b3d7a2SStefano Zampini       if (!n_added) break;
319b1b3d7a2SStefano Zampini     }
320b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
321b1b3d7a2SStefano Zampini 
322883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
32368270318SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_I_layer);CHKERRQ(ierr);
324b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
32568270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
326883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
327883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
328b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
329b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
33068270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
331b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
332b1b3d7a2SStefano Zampini 
333b1b3d7a2SStefano Zampini       /* II block */
334b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
335b1b3d7a2SStefano Zampini     }
336b1b3d7a2SStefano Zampini   } else {
337b1b3d7a2SStefano Zampini     PetscInt n_I;
338b1b3d7a2SStefano Zampini 
339b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
340b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
34168270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
342b1b3d7a2SStefano Zampini 
343b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
344883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
345b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
346b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
347b1b3d7a2SStefano Zampini 
348b1b3d7a2SStefano Zampini       /* II block is the same */
349b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
350b1b3d7a2SStefano Zampini       AE_II = A_II;
351b1b3d7a2SStefano Zampini     }
352b1b3d7a2SStefano Zampini   }
3535a95e1ceSStefano Zampini 
354883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
3555a95e1ceSStefano Zampini   max_subset_size = 0;
356883469d8SStefano Zampini   local_size = 0;
3575a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
3585a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
3595a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
360883469d8SStefano Zampini     local_size += subset_size;
361883469d8SStefano Zampini   }
362883469d8SStefano Zampini 
363883469d8SStefano Zampini   /* Work arrays for local indices */
364883469d8SStefano Zampini   extra = 0;
365883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
366883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
367883469d8SStefano Zampini   }
368883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
369883469d8SStefano Zampini   if (extra) {
370883469d8SStefano Zampini     const PetscInt *idxs;
371883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
372883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
373883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
374883469d8SStefano Zampini   }
375883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
376*eb595f79SStefano Zampini   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
377883469d8SStefano Zampini 
378883469d8SStefano Zampini   /* Get local indices in local numbering */
379883469d8SStefano Zampini   local_size = 0;
3805a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
381883469d8SStefano Zampini     PetscInt j;
382883469d8SStefano Zampini     const    PetscInt *idxs;
383883469d8SStefano Zampini 
3845a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
3855a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
386*eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
387*eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
388*eb595f79SStefano Zampini     auxnum2[i] = subset_size;
389883469d8SStefano Zampini     /* subset indices in local numbering */
390883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
3915a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
392883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
393883469d8SStefano Zampini     local_size += subset_size;
394883469d8SStefano Zampini   }
395883469d8SStefano Zampini 
3965a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
397d2627357SStefano Zampini   Bwork = NULL;
398d2627357SStefano Zampini   pivots = NULL;
3995a95e1ceSStefano Zampini   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
400d2627357SStefano Zampini     PetscScalar lwork;
401d2627357SStefano Zampini 
402d2627357SStefano Zampini     B_lwork = -1;
403d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
404d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
405d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
406d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
407d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
408d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
409d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
410d2627357SStefano Zampini   }
411d2627357SStefano Zampini 
412d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
413*eb595f79SStefano Zampini   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
414*eb595f79SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
415*eb595f79SStefano Zampini   local_size = 0;
416*eb595f79SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
417*eb595f79SStefano Zampini     PetscInt j;
418*eb595f79SStefano Zampini     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
419*eb595f79SStefano Zampini   }
420*eb595f79SStefano Zampini   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
421*eb595f79SStefano Zampini   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
4225a95e1ceSStefano Zampini   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
4235a95e1ceSStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
424d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
425d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
426d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
427d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
4282972d61bSStefano Zampini 
4295a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
4305a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
4315a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
4325a95e1ceSStefano Zampini 
4335a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
4345a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
4355a95e1ceSStefano Zampini     if (subset_size != local_size) {
4365a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
4375a95e1ceSStefano Zampini     }
4385a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
439b1b3d7a2SStefano Zampini   }
440b1b3d7a2SStefano Zampini 
4415a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
4425a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
4435a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
4445a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
4455a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
4465a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
447aa83b6aeSStefano Zampini   } else {
4485a95e1ceSStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
449aa83b6aeSStefano Zampini   }
450b1b3d7a2SStefano Zampini 
4515a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
4525a95e1ceSStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
4535a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps) {
4545a95e1ceSStefano Zampini     Mat         S_Ej_expl;
4555a95e1ceSStefano Zampini     PetscScalar *work;
4565a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
4575a95e1ceSStefano Zampini     PetscBool   Sdense;
4585a95e1ceSStefano Zampini 
4595a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
4605a95e1ceSStefano Zampini     local_size = 0;
461b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
4625a95e1ceSStefano Zampini       IS  is_subset_B;
4635a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
4645a95e1ceSStefano Zampini 
4655a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
4665a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
4675a95e1ceSStefano Zampini       /* EE block */
4685a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
4695a95e1ceSStefano Zampini       /* IE block */
4705a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
4715a95e1ceSStefano Zampini       /* EI block */
4725a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
4735a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
4745a95e1ceSStefano Zampini       } else {
4755a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
4765a95e1ceSStefano Zampini       }
4775a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
4785a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
4795a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
4805a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
481b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
482b1b3d7a2SStefano Zampini         KSP ksp;
483b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
4845a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
485b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
486b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
487b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
488b1b3d7a2SStefano Zampini         KSPType   ksp_type;
489b1b3d7a2SStefano Zampini         PetscInt  n_internal;
4905a95e1ceSStefano Zampini         PetscBool ispcnone;
491b1b3d7a2SStefano Zampini 
492b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
4935a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
494b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
495b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
496b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
497b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
4985a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
4995a95e1ceSStefano Zampini         if (!ispcnone) {
5005a95e1ceSStefano Zampini           PCType pc_type;
501b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
502b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
5035a95e1ceSStefano Zampini         } else {
5045a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
5055a95e1ceSStefano Zampini         }
506b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
507b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
508b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
509b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
510b1b3d7a2SStefano Zampini           if (solver) {
511b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
512b1b3d7a2SStefano Zampini           }
513b1b3d7a2SStefano Zampini         }
514b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
515b1b3d7a2SStefano Zampini       }
5165a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5175a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
5185a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
5195a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
5205a95e1ceSStefano Zampini       if (Sdense) {
5215a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
5225a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
523b1b3d7a2SStefano Zampini         }
5245a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
5255a95e1ceSStefano Zampini       } else {
5265a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
5275a95e1ceSStefano Zampini       }
5285a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
5295a95e1ceSStefano Zampini       local_size += subset_size;
5305a95e1ceSStefano Zampini     }
5315a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
532b1b3d7a2SStefano Zampini     /* free */
533b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
534b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
5355a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
536883469d8SStefano Zampini   } else {
537883469d8SStefano Zampini     Mat         A,F;
538883469d8SStefano Zampini     IS          is_A_all;
5395a95e1ceSStefano Zampini     PetscScalar *work;
540d5574798SStefano Zampini     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx;
541883469d8SStefano Zampini 
542883469d8SStefano Zampini     /* get working mat */
543883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
544883469d8SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
545d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
546883469d8SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
547883469d8SStefano Zampini 
54808122e43SStefano Zampini     if (n_I) {
5499ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
550883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
551883469d8SStefano Zampini       } else {
552883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
553883469d8SStefano Zampini       }
554883469d8SStefano Zampini 
555883469d8SStefano Zampini       /* subsets ordered last */
556883469d8SStefano Zampini       ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
557883469d8SStefano Zampini       for (i=0;i<local_size;i++) {
558883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
559883469d8SStefano Zampini       }
5605a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
561883469d8SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
5625a95e1ceSStefano Zampini #endif
563883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
564883469d8SStefano Zampini 
565883469d8SStefano Zampini       /* factorization step */
5669ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
567883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
568883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
569883469d8SStefano Zampini       } else {
570883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
571883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
572883469d8SStefano Zampini       }
573883469d8SStefano Zampini 
574883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
5755a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
576883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
5775a95e1ceSStefano Zampini #endif
578d5574798SStefano Zampini 
579d5574798SStefano Zampini       /* we can reuse the interior solver if we are not using the economic version */
580d5574798SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
581d5574798SStefano Zampini       if (n_I == n_I_all) {
582d5574798SStefano Zampini         PCBDDCMumpsInterior msolv_ctx;
583d5574798SStefano Zampini 
584d5574798SStefano Zampini         ierr = PetscNew(&msolv_ctx);CHKERRQ(ierr);
585d5574798SStefano Zampini         msolv_ctx->n = n_I;
586d5574798SStefano Zampini         ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
587d5574798SStefano Zampini         msolv_ctx->F = F;
588d5574798SStefano Zampini         ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
589d5574798SStefano Zampini         ierr = PCCreate(PETSC_COMM_SELF,&sub_schurs->interior_solver);CHKERRQ(ierr);
590d5574798SStefano Zampini         ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
591d5574798SStefano Zampini         ierr = PCSetOperators(sub_schurs->interior_solver,A_II,A_II);CHKERRQ(ierr);
592d5574798SStefano Zampini         ierr = PCSetType(sub_schurs->interior_solver,PCSHELL);CHKERRQ(ierr);
593d5574798SStefano Zampini         ierr = PCShellSetContext(sub_schurs->interior_solver,msolv_ctx);CHKERRQ(ierr);
594d5574798SStefano Zampini         ierr = PCShellSetApply(sub_schurs->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
595d5574798SStefano Zampini         ierr = PCShellSetDestroy(sub_schurs->interior_solver,PCBDDCMumpsInteriorDestroy);CHKERRQ(ierr);
596d5574798SStefano Zampini       }
597883469d8SStefano Zampini       ierr = MatDestroy(&F);CHKERRQ(ierr);
59808122e43SStefano Zampini     } else {
59908122e43SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
60008122e43SStefano Zampini     }
60108122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
602d2627357SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
6035db18549SStefano Zampini 
60412d906b1SStefano Zampini     if (compute_Stilda) {
605a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
606a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
607a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
608a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
60908122e43SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
61008122e43SStefano Zampini       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
61108122e43SStefano Zampini       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
61208122e43SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
61306a4b1faSStefano Zampini 
61412d906b1SStefano Zampini       /* compute St^-1 */
61512d906b1SStefano Zampini       if (local_size) { /* multilevel guard */
616d2627357SStefano Zampini         PetscScalar *vals;
617f6f667cfSStefano Zampini 
6185a95e1ceSStefano Zampini         ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
619d2627357SStefano Zampini         ierr = MatDuplicate(S_all,MAT_COPY_VALUES,&S_all_inv);CHKERRQ(ierr);
620d2627357SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&vals);CHKERRQ(ierr);
621d2627357SStefano Zampini         if (!sub_schurs->is_hermitian) {
622d2627357SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
623d2627357SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
624d2627357SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
625d2627357SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
626d2627357SStefano Zampini         } else {
627d2627357SStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
628d2627357SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
629d2627357SStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
630d2627357SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
631d2627357SStefano Zampini         }
632d2627357SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&vals);CHKERRQ(ierr);
633d2627357SStefano Zampini       }
63412d906b1SStefano Zampini     }
635d2627357SStefano Zampini 
6369087bf02SStefano Zampini     /* Work arrays */
6379087bf02SStefano Zampini     if (sub_schurs->n_subs == 1) {
6389087bf02SStefano Zampini       ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
6399087bf02SStefano Zampini     } else {
6409087bf02SStefano Zampini       ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
6419087bf02SStefano Zampini     }
6429087bf02SStefano Zampini 
6435a95e1ceSStefano Zampini     local_size = 0;
64465d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
6455a95e1ceSStefano Zampini       Mat S_Ej;
64665d8bf0aSStefano Zampini       IS  is_E;
64765d8bf0aSStefano Zampini       PetscInt j;
64865d8bf0aSStefano Zampini 
6495a95e1ceSStefano Zampini       /* get S_E */
650b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6519087bf02SStefano Zampini       if (sub_schurs->n_subs == 1) {
6529087bf02SStefano Zampini         ierr = MatDenseGetArray(S_all,&work);CHKERRQ(ierr);
6539087bf02SStefano Zampini         S_Ej = NULL;
6549087bf02SStefano Zampini         is_E = NULL;
6559087bf02SStefano Zampini       } else {
6565a95e1ceSStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,local_size,1,&is_E);CHKERRQ(ierr);
6575a95e1ceSStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej);CHKERRQ(ierr);
6585a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
6599087bf02SStefano Zampini       }
6605a95e1ceSStefano Zampini       /* insert S_E values */
661a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
662a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
663a1337663SStefano Zampini       }
6645a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
665a1337663SStefano Zampini 
6665a95e1ceSStefano Zampini       /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
667d2627357SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
6685a95e1ceSStefano Zampini         /* get S_E^-1 */
66908122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
67008122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6712972d61bSStefano Zampini         if (!sub_schurs->is_hermitian) {
6725a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
67308122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
6745a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
67508122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6762972d61bSStefano Zampini         } else {
6775a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
6782972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
6795a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
6802972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
6812972d61bSStefano Zampini         }
68208122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
6835a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6845a95e1ceSStefano Zampini 
6855a95e1ceSStefano Zampini         /* get St_E^-1 */
6869087bf02SStefano Zampini         if (sub_schurs->n_subs == 1) {
6879087bf02SStefano Zampini           ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
6889087bf02SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&work);CHKERRQ(ierr);
6899087bf02SStefano Zampini         } else {
6905a95e1ceSStefano Zampini           ierr = MatGetSubMatrix(S_all_inv,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
6919087bf02SStefano Zampini         }
6925a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6939087bf02SStefano Zampini         if (sub_schurs->n_subs == 1) {
6949087bf02SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&work);CHKERRQ(ierr);
6959087bf02SStefano Zampini         }
6969087bf02SStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
6979087bf02SStefano Zampini       } else if (sub_schurs->n_subs == 1) {
6989087bf02SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
69908122e43SStefano Zampini       }
7005a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
70165d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
702883469d8SStefano Zampini       local_size += subset_size;
703883469d8SStefano Zampini     }
7049087bf02SStefano Zampini     if (sub_schurs->n_subs == 1) {
7059087bf02SStefano Zampini       ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
7069087bf02SStefano Zampini     } else {
7075ec10c6aSStefano Zampini       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
7085db18549SStefano Zampini     }
7099087bf02SStefano Zampini   }
710a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
711a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
712d2627357SStefano Zampini   ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
7135db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7145db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7155a95e1ceSStefano Zampini   if (compute_Stilda) {
716a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
717a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71808122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71908122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72008122e43SStefano Zampini   }
721a1337663SStefano Zampini 
7225db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
7235db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
7243927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
7253927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
7265a95e1ceSStefano Zampini 
7275db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
7285a95e1ceSStefano Zampini   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
729d648f858SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
730d648f858SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
73108122e43SStefano Zampini 
732ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
733ac632422SStefano Zampini   if (faster_deluxe) {
7345a95e1ceSStefano Zampini     Mat         tmpmat;
7355a95e1ceSStefano Zampini     PetscScalar *array;
7365a95e1ceSStefano Zampini     PetscInt    cum;
7375a95e1ceSStefano Zampini 
7385a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
7395a95e1ceSStefano Zampini     cum = 0;
7405a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7415a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7425a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
7435a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
7445a95e1ceSStefano Zampini       if (!sub_schurs->is_hermitian) {
7455a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
7465a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
7475a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
7485a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
7495a95e1ceSStefano Zampini       } else {
7505a95e1ceSStefano Zampini         PetscInt j,k;
7515a95e1ceSStefano Zampini 
7525a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
7535a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7545a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
7555a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7565a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
7575a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
7585a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
7595a95e1ceSStefano Zampini           }
7605a95e1ceSStefano Zampini         }
7615a95e1ceSStefano Zampini       }
7625a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
7635a95e1ceSStefano Zampini       cum += subset_size*subset_size;
7645a95e1ceSStefano Zampini     }
7655a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
7665a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
7675a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
768ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
7695a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
7705a95e1ceSStefano Zampini   }
7715a95e1ceSStefano Zampini 
772f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
7735a95e1ceSStefano Zampini   if (compute_Stilda) {
774a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
775a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
776d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
777d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
77808122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
77908122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
780d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
781d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
78208122e43SStefano Zampini   }
7833202ece2SStefano Zampini 
7845a95e1ceSStefano Zampini   /* free workspace */
78506a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
786a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
787a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
78808122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
7893202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
7905db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
7915a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
792b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
793b1b3d7a2SStefano Zampini }
794b1b3d7a2SStefano Zampini 
795b1b3d7a2SStefano Zampini #undef __FUNCT__
796b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
7975a95e1ceSStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
798b1b3d7a2SStefano Zampini {
7999bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
8005a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
801b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
802b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
803b1b3d7a2SStefano Zampini 
804b1b3d7a2SStefano Zampini   PetscFunctionBegin;
805b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
806b1b3d7a2SStefano Zampini   if (!is_sorted) {
807b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
808b1b3d7a2SStefano Zampini   }
809b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
810b1b3d7a2SStefano Zampini   if (!is_sorted) {
811b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
812b1b3d7a2SStefano Zampini   }
813b1b3d7a2SStefano Zampini 
814b1b3d7a2SStefano Zampini   /* reset any previous data */
815b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
816b1b3d7a2SStefano Zampini 
8175a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
8189bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
819b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
82008122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
82108122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
822b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
823b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
824b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
825b1b3d7a2SStefano Zampini   }
826b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
827b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
82808122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
829b1b3d7a2SStefano Zampini   }
830b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
831b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
832b1b3d7a2SStefano Zampini 
8335a95e1ceSStefano Zampini   /* Determine if MUMPS cane be used */
834883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
835883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
8364c6709b3SStefano Zampini   sub_schurs->use_mumps = (PetscBool)(!!A);
837883469d8SStefano Zampini #endif
838b1b3d7a2SStefano Zampini 
839b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
840aa83b6aeSStefano Zampini   if (A) {
8419ab7bb16SStefano Zampini     PetscBool isseqaij;
8429ab7bb16SStefano Zampini 
8439ab7bb16SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
8449ab7bb16SStefano Zampini     if (isseqaij) {
8451e9c79c2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
8461e9c79c2SStefano Zampini       sub_schurs->A = A;
8479ab7bb16SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
8489ab7bb16SStefano Zampini       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
8499ab7bb16SStefano Zampini     }
8501e9c79c2SStefano Zampini   }
851b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
852b1b3d7a2SStefano Zampini   sub_schurs->S = S;
853b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
854b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
855b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
856b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
8575db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
8585db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
8595db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
8605db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
8615a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
862b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
8639bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
864b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
865b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
866b96c3477SStefano Zampini     }
8679bb4a8caSStefano Zampini   }
8689bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
869b1b3d7a2SStefano Zampini 
870b96c3477SStefano Zampini 
871b96c3477SStefano Zampini   /* allocate space for schur complements */
872b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
873b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
87408122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
875b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
876b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
877b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
878b1b3d7a2SStefano Zampini }
879b1b3d7a2SStefano Zampini 
880b1b3d7a2SStefano Zampini #undef __FUNCT__
88134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
88234a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
88334a97f8cSStefano Zampini {
88434a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
88534a97f8cSStefano Zampini   PetscErrorCode  ierr;
88634a97f8cSStefano Zampini 
88734a97f8cSStefano Zampini   PetscFunctionBegin;
88834a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
8895ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
89034a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
89134a97f8cSStefano Zampini   PetscFunctionReturn(0);
89234a97f8cSStefano Zampini }
89334a97f8cSStefano Zampini 
89434a97f8cSStefano Zampini #undef __FUNCT__
89534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
89634a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
89734a97f8cSStefano Zampini {
89834a97f8cSStefano Zampini   PetscErrorCode ierr;
89934a97f8cSStefano Zampini 
90034a97f8cSStefano Zampini   PetscFunctionBegin;
90134a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
90234a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
90334a97f8cSStefano Zampini   PetscFunctionReturn(0);
90434a97f8cSStefano Zampini }
90534a97f8cSStefano Zampini 
90634a97f8cSStefano Zampini #undef __FUNCT__
90734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
90834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
90934a97f8cSStefano Zampini {
91034a97f8cSStefano Zampini   PetscInt       i;
91134a97f8cSStefano Zampini   PetscErrorCode ierr;
91234a97f8cSStefano Zampini 
91334a97f8cSStefano Zampini   PetscFunctionBegin;
9141e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
915b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
916b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
917b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
9185db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
9195db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
92041c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
92141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
92208122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
923a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
9245db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
9259bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
92608122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
92708122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
92834a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
929b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
93034a97f8cSStefano Zampini   }
93168270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
9325ff63025SStefano Zampini   if (sub_schurs->n_subs) {
933b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
9343dc780c3SStefano Zampini   }
935d5574798SStefano Zampini   ierr = PCDestroy(&sub_schurs->interior_solver);CHKERRQ(ierr);
93634a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
93734a97f8cSStefano Zampini   PetscFunctionReturn(0);
93834a97f8cSStefano Zampini }
93934a97f8cSStefano Zampini 
94034a97f8cSStefano Zampini #undef __FUNCT__
94134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
9422a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
94334a97f8cSStefano Zampini {
94434a97f8cSStefano Zampini   PetscInt       i,j,n;
94534a97f8cSStefano Zampini   PetscErrorCode ierr;
94634a97f8cSStefano Zampini 
94734a97f8cSStefano Zampini   PetscFunctionBegin;
94834a97f8cSStefano Zampini   n = 0;
94934a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
95034a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
95134a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
95234a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
95334a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
95434a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
95534a97f8cSStefano Zampini         queue_tip[n] = dof;
95634a97f8cSStefano Zampini         n++;
95734a97f8cSStefano Zampini       }
95834a97f8cSStefano Zampini     }
95934a97f8cSStefano Zampini   }
96034a97f8cSStefano Zampini   *n_added = n;
96134a97f8cSStefano Zampini   PetscFunctionReturn(0);
96234a97f8cSStefano Zampini }
963