xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision f6f667cfabdb84fae686043df9fba3b425683839)
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;
1995a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
200883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
20108122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
20206a4b1faSStefano Zampini   PetscScalar            *Bwork;
2035a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2045a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2055a95e1ceSStefano Zampini   MPI_Comm               comm_n;
206b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
207b1b3d7a2SStefano Zampini 
208b1b3d7a2SStefano Zampini   PetscFunctionBegin;
2095a95e1ceSStefano Zampini   /* preliminary checks */
2105a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
2115a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
2125a95e1ceSStefano Zampini   }
2135a95e1ceSStefano Zampini   /* determine if we are dealing with hermitian positive definite problems */
2145a95e1ceSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2155a95e1ceSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
2165a95e1ceSStefano Zampini   if (sub_schurs->A) {
2175a95e1ceSStefano Zampini     PetscInt lsize;
2185a95e1ceSStefano Zampini 
2195a95e1ceSStefano Zampini     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
2205a95e1ceSStefano Zampini     if (lsize) {
2215a95e1ceSStefano Zampini       ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
2225a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
2235a95e1ceSStefano Zampini         PetscScalar val;
2245a95e1ceSStefano Zampini         Vec         vec1,vec2;
2255a95e1ceSStefano Zampini 
2265a95e1ceSStefano Zampini         ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
2275a95e1ceSStefano Zampini         ierr = VecSetRandom(vec1,NULL);
2285a95e1ceSStefano Zampini         ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
2295a95e1ceSStefano Zampini         ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
2305a95e1ceSStefano Zampini         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
2315a95e1ceSStefano Zampini         if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
2325a95e1ceSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
2335a95e1ceSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
2345a95e1ceSStefano Zampini       }
2355a95e1ceSStefano Zampini     } else {
2365a95e1ceSStefano Zampini       sub_schurs->is_hermitian = PETSC_TRUE;
2375a95e1ceSStefano Zampini       sub_schurs->is_posdef = PETSC_TRUE;
2385a95e1ceSStefano Zampini     }
2395a95e1ceSStefano Zampini     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
2405a95e1ceSStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported");
2415a95e1ceSStefano Zampini     }
2425a95e1ceSStefano Zampini   }
2435a95e1ceSStefano Zampini   /* restrict work on active processes */
2445a95e1ceSStefano Zampini   color = 0;
2455a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
2465a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
2475a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
2485a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
2495a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
2505a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
2515a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
2525a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
2535a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
2545a95e1ceSStefano Zampini     PetscFunctionReturn(0);
2555a95e1ceSStefano Zampini   }
2565a95e1ceSStefano Zampini 
257b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
258883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
259f11841e3SStefano Zampini     PetscBool isseqaij;
260b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
261f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
262f11841e3SStefano Zampini     if (!isseqaij) {
263f11841e3SStefano Zampini       ierr = MatConvert(A_BB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BB);CHKERRQ(ierr);
264f11841e3SStefano Zampini       ierr = MatConvert(A_IB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_IB);CHKERRQ(ierr);
265f11841e3SStefano Zampini       ierr = MatConvert(A_BI,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BI);CHKERRQ(ierr);
266f11841e3SStefano Zampini     }
267a58a30b4SStefano Zampini   } else {
2685a95e1ceSStefano Zampini     A_II = NULL;
2695a95e1ceSStefano Zampini     A_IB = NULL;
2705a95e1ceSStefano Zampini     A_BI = NULL;
2715a95e1ceSStefano Zampini     A_BB = NULL;
272b1b3d7a2SStefano Zampini   }
2735a95e1ceSStefano Zampini   S_all = NULL;
2745a95e1ceSStefano Zampini   S_all_inv = NULL;
2755a95e1ceSStefano Zampini   S_Ej_tilda_all = NULL;
2765a95e1ceSStefano Zampini   S_Ej_inv_all = NULL;
277b1b3d7a2SStefano Zampini 
278b1b3d7a2SStefano Zampini   /* determine interior problems */
279b96c3477SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
2803dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
2813dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
282b1b3d7a2SStefano Zampini     PetscBT                touched;
283b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
284b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
285b1b3d7a2SStefano Zampini 
2863dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
2873dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
2883dc780c3SStefano Zampini     }
289b1b3d7a2SStefano Zampini     /* get sizes */
290b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
291b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
292b1b3d7a2SStefano Zampini 
293b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
294b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
295b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
296b1b3d7a2SStefano Zampini 
297b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
298b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
299b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
300b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
301b1b3d7a2SStefano Zampini     }
302b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
303b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
304b1b3d7a2SStefano Zampini 
305b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
306b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
307b1b3d7a2SStefano Zampini     n_prev_added = n_B;
308b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
309b1b3d7a2SStefano Zampini       PetscInt n_added;
310b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
311b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
312b1b3d7a2SStefano 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);
313b1b3d7a2SStefano Zampini       }
314b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
315b1b3d7a2SStefano Zampini       n_prev_added = n_added;
316b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
317b1b3d7a2SStefano Zampini       if (!n_added) break;
318b1b3d7a2SStefano Zampini     }
319b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
320b1b3d7a2SStefano Zampini 
321883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
32268270318SStefano 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);
323b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
32468270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
325883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
326883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
327b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
328b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
32968270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
330b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
331b1b3d7a2SStefano Zampini 
332b1b3d7a2SStefano Zampini       /* II block */
333b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
334b1b3d7a2SStefano Zampini     }
335b1b3d7a2SStefano Zampini   } else {
336b1b3d7a2SStefano Zampini     PetscInt n_I;
337b1b3d7a2SStefano Zampini 
338b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
339b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
34068270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
341b1b3d7a2SStefano Zampini 
342b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
343883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
344b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
345b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
346b1b3d7a2SStefano Zampini 
347b1b3d7a2SStefano Zampini       /* II block is the same */
348b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
349b1b3d7a2SStefano Zampini       AE_II = A_II;
350b1b3d7a2SStefano Zampini     }
351b1b3d7a2SStefano Zampini   }
3525a95e1ceSStefano Zampini 
353883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
3545a95e1ceSStefano Zampini   max_subset_size = 0;
355883469d8SStefano Zampini   local_size = 0;
3565a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
3575a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
3585a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
359883469d8SStefano Zampini     local_size += subset_size;
360883469d8SStefano Zampini   }
361883469d8SStefano Zampini 
362883469d8SStefano Zampini   /* Work arrays for local indices */
363883469d8SStefano Zampini   extra = 0;
364883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
365883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
366883469d8SStefano Zampini   }
367883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
368883469d8SStefano Zampini   if (extra) {
369883469d8SStefano Zampini     const PetscInt *idxs;
370883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
371883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
372883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
373883469d8SStefano Zampini   }
374883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
375883469d8SStefano Zampini 
376883469d8SStefano Zampini   /* Get local indices in local numbering */
377883469d8SStefano Zampini   local_size = 0;
3785a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
379883469d8SStefano Zampini     PetscInt j;
380883469d8SStefano Zampini     const    PetscInt *idxs;
381883469d8SStefano Zampini 
3825a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
3835a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
384883469d8SStefano Zampini     /* subset indices in local numbering */
385883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
3865a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
387883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
388883469d8SStefano Zampini     local_size += subset_size;
389883469d8SStefano Zampini   }
390883469d8SStefano Zampini 
3915a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
392d2627357SStefano Zampini   Bwork = NULL;
393d2627357SStefano Zampini   pivots = NULL;
3945a95e1ceSStefano Zampini   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
395d2627357SStefano Zampini     PetscScalar lwork;
396d2627357SStefano Zampini 
397d2627357SStefano Zampini     B_lwork = -1;
398d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
399d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
400d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
401d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
402d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
403d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
404d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
405d2627357SStefano Zampini   }
406d2627357SStefano Zampini 
407d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
4085a95e1ceSStefano Zampini   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,local_size,all_local_idx_N+extra,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr);
4095a95e1ceSStefano Zampini   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
4105a95e1ceSStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
411d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
412d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
413d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
414d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
4152972d61bSStefano Zampini 
4165a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
4175a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
4185a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
4195a95e1ceSStefano Zampini 
4205a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
4215a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
4225a95e1ceSStefano Zampini     if (subset_size != local_size) {
4235a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
4245a95e1ceSStefano Zampini     }
4255a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
426b1b3d7a2SStefano Zampini   }
427b1b3d7a2SStefano Zampini 
4285a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
4295a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
4305a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
4315a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
4325a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
4335a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
434aa83b6aeSStefano Zampini   } else {
4355a95e1ceSStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
436aa83b6aeSStefano Zampini   }
437b1b3d7a2SStefano Zampini 
4385a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
4395a95e1ceSStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
4405a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps) {
4415a95e1ceSStefano Zampini     Mat         S_Ej_expl;
4425a95e1ceSStefano Zampini     PetscScalar *work;
4435a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
4445a95e1ceSStefano Zampini     PetscBool   Sdense;
4455a95e1ceSStefano Zampini 
4465a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
4475a95e1ceSStefano Zampini     local_size = 0;
448b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
4495a95e1ceSStefano Zampini       IS  is_subset_B;
4505a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
4515a95e1ceSStefano Zampini 
4525a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
4535a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
4545a95e1ceSStefano Zampini       /* EE block */
4555a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
4565a95e1ceSStefano Zampini       /* IE block */
4575a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
4585a95e1ceSStefano Zampini       /* EI block */
4595a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
4605a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
4615a95e1ceSStefano Zampini       } else {
4625a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
4635a95e1ceSStefano Zampini       }
4645a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
4655a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
4665a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
4675a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
468b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
469b1b3d7a2SStefano Zampini         KSP ksp;
470b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
4715a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
472b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
473b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
474b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
475b1b3d7a2SStefano Zampini         KSPType   ksp_type;
476b1b3d7a2SStefano Zampini         PetscInt  n_internal;
4775a95e1ceSStefano Zampini         PetscBool ispcnone;
478b1b3d7a2SStefano Zampini 
479b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
4805a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
481b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
482b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
483b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
484b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
4855a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
4865a95e1ceSStefano Zampini         if (!ispcnone) {
4875a95e1ceSStefano Zampini           PCType pc_type;
488b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
489b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
4905a95e1ceSStefano Zampini         } else {
4915a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
4925a95e1ceSStefano Zampini         }
493b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
494b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
495b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
496b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
497b1b3d7a2SStefano Zampini           if (solver) {
498b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
499b1b3d7a2SStefano Zampini           }
500b1b3d7a2SStefano Zampini         }
501b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
502b1b3d7a2SStefano Zampini       }
5035a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5045a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
5055a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
5065a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
5075a95e1ceSStefano Zampini       if (Sdense) {
5085a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
5095a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
510b1b3d7a2SStefano Zampini         }
5115a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
5125a95e1ceSStefano Zampini       } else {
5135a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
5145a95e1ceSStefano Zampini       }
5155a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
5165a95e1ceSStefano Zampini       local_size += subset_size;
5175a95e1ceSStefano Zampini     }
5185a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
519b1b3d7a2SStefano Zampini     /* free */
520b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
521b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
5225a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
523883469d8SStefano Zampini   } else {
524883469d8SStefano Zampini     Mat         A,F;
525883469d8SStefano Zampini     IS          is_A_all;
5265a95e1ceSStefano Zampini     PetscScalar *work;
527d5574798SStefano Zampini     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx;
528883469d8SStefano Zampini 
529883469d8SStefano Zampini     /* get working mat */
530883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
531883469d8SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
532d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
533883469d8SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
534883469d8SStefano Zampini 
53508122e43SStefano Zampini     if (n_I) {
5369ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
537883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
538883469d8SStefano Zampini       } else {
539883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
540883469d8SStefano Zampini       }
541883469d8SStefano Zampini 
542883469d8SStefano Zampini       /* subsets ordered last */
543883469d8SStefano Zampini       ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
544883469d8SStefano Zampini       for (i=0;i<local_size;i++) {
545883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
546883469d8SStefano Zampini       }
5475a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
548883469d8SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
5495a95e1ceSStefano Zampini #endif
550883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
551883469d8SStefano Zampini 
552883469d8SStefano Zampini       /* factorization step */
5539ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
554883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
555883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
556883469d8SStefano Zampini       } else {
557883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
558883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
559883469d8SStefano Zampini       }
560883469d8SStefano Zampini 
561883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
5625a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
563883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
5645a95e1ceSStefano Zampini #endif
565d5574798SStefano Zampini 
566d5574798SStefano Zampini       /* we can reuse the interior solver if we are not using the economic version */
567d5574798SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
568d5574798SStefano Zampini       if (n_I == n_I_all) {
569d5574798SStefano Zampini         PCBDDCMumpsInterior msolv_ctx;
570d5574798SStefano Zampini 
571d5574798SStefano Zampini         ierr = PetscNew(&msolv_ctx);CHKERRQ(ierr);
572d5574798SStefano Zampini         msolv_ctx->n = n_I;
573d5574798SStefano Zampini         ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
574d5574798SStefano Zampini         msolv_ctx->F = F;
575d5574798SStefano Zampini         ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
576d5574798SStefano Zampini         ierr = PCCreate(PETSC_COMM_SELF,&sub_schurs->interior_solver);CHKERRQ(ierr);
577d5574798SStefano Zampini         ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
578d5574798SStefano Zampini         ierr = PCSetOperators(sub_schurs->interior_solver,A_II,A_II);CHKERRQ(ierr);
579d5574798SStefano Zampini         ierr = PCSetType(sub_schurs->interior_solver,PCSHELL);CHKERRQ(ierr);
580d5574798SStefano Zampini         ierr = PCShellSetContext(sub_schurs->interior_solver,msolv_ctx);CHKERRQ(ierr);
581d5574798SStefano Zampini         ierr = PCShellSetApply(sub_schurs->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
582d5574798SStefano Zampini         ierr = PCShellSetDestroy(sub_schurs->interior_solver,PCBDDCMumpsInteriorDestroy);CHKERRQ(ierr);
583d5574798SStefano Zampini       }
584883469d8SStefano Zampini       ierr = MatDestroy(&F);CHKERRQ(ierr);
58508122e43SStefano Zampini     } else {
58608122e43SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
58708122e43SStefano Zampini     }
58808122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
589d2627357SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
5905db18549SStefano Zampini 
5915a95e1ceSStefano Zampini     if (compute_Stilda) { /* TODO PICKUP BETTER NAMES */
592a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
593a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
594a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
595a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
59608122e43SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
59708122e43SStefano Zampini       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
59808122e43SStefano Zampini       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
59908122e43SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
600a1337663SStefano Zampini     }
60106a4b1faSStefano Zampini 
6025a95e1ceSStefano Zampini     /* Get St^-1 */
6035a95e1ceSStefano Zampini     if (compute_Stilda) {
604d2627357SStefano Zampini       PetscScalar *vals;
605*f6f667cfSStefano Zampini 
6065a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
607d2627357SStefano Zampini       ierr = MatDuplicate(S_all,MAT_COPY_VALUES,&S_all_inv);CHKERRQ(ierr);
608d2627357SStefano Zampini       ierr = MatDenseGetArray(S_all_inv,&vals);CHKERRQ(ierr);
609d2627357SStefano Zampini       if (!sub_schurs->is_hermitian) {
610d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
611d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
612d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
613d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
614d2627357SStefano Zampini       } else {
615d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
616d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
617d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
618d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
619d2627357SStefano Zampini       }
620d2627357SStefano Zampini       ierr = MatDenseRestoreArray(S_all_inv,&vals);CHKERRQ(ierr);
621d2627357SStefano Zampini     }
622d2627357SStefano Zampini 
6239087bf02SStefano Zampini     /* Work arrays */
6249087bf02SStefano Zampini     if (sub_schurs->n_subs == 1) {
6259087bf02SStefano Zampini       ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
6269087bf02SStefano Zampini     } else {
6279087bf02SStefano Zampini       ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
6289087bf02SStefano Zampini     }
6299087bf02SStefano Zampini 
6305a95e1ceSStefano Zampini     local_size = 0;
63165d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
6325a95e1ceSStefano Zampini       Mat S_Ej;
63365d8bf0aSStefano Zampini       IS  is_E;
63465d8bf0aSStefano Zampini       PetscInt j;
63565d8bf0aSStefano Zampini 
6365a95e1ceSStefano Zampini       /* get S_E */
637b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6389087bf02SStefano Zampini       if (sub_schurs->n_subs == 1) {
6399087bf02SStefano Zampini         ierr = MatDenseGetArray(S_all,&work);CHKERRQ(ierr);
6409087bf02SStefano Zampini         S_Ej = NULL;
6419087bf02SStefano Zampini         is_E = NULL;
6429087bf02SStefano Zampini       } else {
6435a95e1ceSStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,local_size,1,&is_E);CHKERRQ(ierr);
6445a95e1ceSStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej);CHKERRQ(ierr);
6455a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
6469087bf02SStefano Zampini       }
6475a95e1ceSStefano Zampini       /* insert S_E values */
648a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
649a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
650a1337663SStefano Zampini       }
6515a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
652a1337663SStefano Zampini 
6535a95e1ceSStefano Zampini       /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
654d2627357SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
6555a95e1ceSStefano Zampini         /* get S_E^-1 */
65608122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
65708122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6582972d61bSStefano Zampini         if (!sub_schurs->is_hermitian) {
6595a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
66008122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
6615a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
66208122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6632972d61bSStefano Zampini         } else {
6645a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
6652972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
6665a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
6672972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
6682972d61bSStefano Zampini         }
66908122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
6705a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6715a95e1ceSStefano Zampini 
6725a95e1ceSStefano Zampini         /* get St_E^-1 */
6739087bf02SStefano Zampini         if (sub_schurs->n_subs == 1) {
6749087bf02SStefano Zampini           ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
6759087bf02SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&work);CHKERRQ(ierr);
6769087bf02SStefano Zampini         } else {
6775a95e1ceSStefano Zampini           ierr = MatGetSubMatrix(S_all_inv,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
6789087bf02SStefano Zampini         }
6795a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6809087bf02SStefano Zampini         if (sub_schurs->n_subs == 1) {
6819087bf02SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&work);CHKERRQ(ierr);
6829087bf02SStefano Zampini         }
6839087bf02SStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
6849087bf02SStefano Zampini       } else if (sub_schurs->n_subs == 1) {
6859087bf02SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
68608122e43SStefano Zampini       }
6875a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
68865d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
689883469d8SStefano Zampini       local_size += subset_size;
690883469d8SStefano Zampini     }
6919087bf02SStefano Zampini     if (sub_schurs->n_subs == 1) {
6929087bf02SStefano Zampini       ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
6939087bf02SStefano Zampini     } else {
6945ec10c6aSStefano Zampini       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
6955db18549SStefano Zampini     }
6969087bf02SStefano Zampini   }
697a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
698a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
699d2627357SStefano Zampini   ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
7005db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7015db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7025a95e1ceSStefano Zampini   if (compute_Stilda) {
703a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
704a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70508122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70608122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70708122e43SStefano Zampini   }
708a1337663SStefano Zampini 
7095db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
7105db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
7113927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
7123927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
7135a95e1ceSStefano Zampini 
7145db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
7155a95e1ceSStefano Zampini   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
716d648f858SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
717d648f858SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
71808122e43SStefano Zampini 
719ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
720ac632422SStefano Zampini   if (faster_deluxe) {
7215a95e1ceSStefano Zampini     Mat         tmpmat;
7225a95e1ceSStefano Zampini     PetscScalar *array;
7235a95e1ceSStefano Zampini     PetscInt    cum;
7245a95e1ceSStefano Zampini 
7255a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
7265a95e1ceSStefano Zampini     cum = 0;
7275a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7285a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7295a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
7305a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
7315a95e1ceSStefano Zampini       if (!sub_schurs->is_hermitian) {
7325a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
7335a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
7345a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
7355a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
7365a95e1ceSStefano Zampini       } else {
7375a95e1ceSStefano Zampini         PetscInt j,k;
7385a95e1ceSStefano Zampini 
7395a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
7405a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7415a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
7425a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7435a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
7445a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
7455a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
7465a95e1ceSStefano Zampini           }
7475a95e1ceSStefano Zampini         }
7485a95e1ceSStefano Zampini       }
7495a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
7505a95e1ceSStefano Zampini       cum += subset_size*subset_size;
7515a95e1ceSStefano Zampini     }
7525a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
7535a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
7545a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
755ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
7565a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
7575a95e1ceSStefano Zampini   }
7585a95e1ceSStefano Zampini 
759*f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
7605a95e1ceSStefano Zampini   if (compute_Stilda) {
761a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
762a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
763d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
764d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
76508122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
76608122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
767d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
768d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
76908122e43SStefano Zampini   }
7703202ece2SStefano Zampini 
7715a95e1ceSStefano Zampini   /* free workspace */
77206a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
773a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
774a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
77508122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
7763202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
7775db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
7785a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
779b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
780b1b3d7a2SStefano Zampini }
781b1b3d7a2SStefano Zampini 
782b1b3d7a2SStefano Zampini #undef __FUNCT__
783b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
7845a95e1ceSStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
785b1b3d7a2SStefano Zampini {
7869bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
7875a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
788b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
789b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
790b1b3d7a2SStefano Zampini 
791b1b3d7a2SStefano Zampini   PetscFunctionBegin;
792b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
793b1b3d7a2SStefano Zampini   if (!is_sorted) {
794b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
795b1b3d7a2SStefano Zampini   }
796b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
797b1b3d7a2SStefano Zampini   if (!is_sorted) {
798b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
799b1b3d7a2SStefano Zampini   }
800b1b3d7a2SStefano Zampini 
801b1b3d7a2SStefano Zampini   /* reset any previous data */
802b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
803b1b3d7a2SStefano Zampini 
8045a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
8059bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
806b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
80708122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
80808122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
809b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
810b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
811b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
812b1b3d7a2SStefano Zampini   }
813b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
814b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
81508122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
816b1b3d7a2SStefano Zampini   }
817b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
818b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
819b1b3d7a2SStefano Zampini 
8205a95e1ceSStefano Zampini   /* Determine if MUMPS cane be used */
821883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
822883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
8234c6709b3SStefano Zampini   sub_schurs->use_mumps = (PetscBool)(!!A);
824883469d8SStefano Zampini #endif
825b1b3d7a2SStefano Zampini 
826b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
827aa83b6aeSStefano Zampini   if (A) {
8289ab7bb16SStefano Zampini     PetscBool isseqaij;
8299ab7bb16SStefano Zampini 
8309ab7bb16SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
8319ab7bb16SStefano Zampini     if (isseqaij) {
8321e9c79c2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
8331e9c79c2SStefano Zampini       sub_schurs->A = A;
8349ab7bb16SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
8359ab7bb16SStefano Zampini       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
8369ab7bb16SStefano Zampini     }
8371e9c79c2SStefano Zampini   }
838b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
839b1b3d7a2SStefano Zampini   sub_schurs->S = S;
840b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
841b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
842b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
843b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
8445db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
8455db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
8465db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
8475db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
8485a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
849b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
8509bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
851b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
852b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
853b96c3477SStefano Zampini     }
8549bb4a8caSStefano Zampini   }
8559bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
856b1b3d7a2SStefano Zampini 
857b96c3477SStefano Zampini 
858b96c3477SStefano Zampini   /* allocate space for schur complements */
859b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
860b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
86108122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
862b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
863b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
864b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
865b1b3d7a2SStefano Zampini }
866b1b3d7a2SStefano Zampini 
867b1b3d7a2SStefano Zampini #undef __FUNCT__
86834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
86934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
87034a97f8cSStefano Zampini {
87134a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
87234a97f8cSStefano Zampini   PetscErrorCode  ierr;
87334a97f8cSStefano Zampini 
87434a97f8cSStefano Zampini   PetscFunctionBegin;
87534a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
8765ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
87734a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
87834a97f8cSStefano Zampini   PetscFunctionReturn(0);
87934a97f8cSStefano Zampini }
88034a97f8cSStefano Zampini 
88134a97f8cSStefano Zampini #undef __FUNCT__
88234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
88334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
88434a97f8cSStefano Zampini {
88534a97f8cSStefano Zampini   PetscErrorCode ierr;
88634a97f8cSStefano Zampini 
88734a97f8cSStefano Zampini   PetscFunctionBegin;
88834a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
88934a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
89034a97f8cSStefano Zampini   PetscFunctionReturn(0);
89134a97f8cSStefano Zampini }
89234a97f8cSStefano Zampini 
89334a97f8cSStefano Zampini #undef __FUNCT__
89434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
89534a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
89634a97f8cSStefano Zampini {
89734a97f8cSStefano Zampini   PetscInt       i;
89834a97f8cSStefano Zampini   PetscErrorCode ierr;
89934a97f8cSStefano Zampini 
90034a97f8cSStefano Zampini   PetscFunctionBegin;
9011e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
902b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
903b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
904b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
9055db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
9065db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
90741c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
90841c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
90908122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
910a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
9115db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
9129bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
91308122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
91408122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
91534a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
916b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
91734a97f8cSStefano Zampini   }
91868270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
9195ff63025SStefano Zampini   if (sub_schurs->n_subs) {
920b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
9213dc780c3SStefano Zampini   }
922d5574798SStefano Zampini   ierr = PCDestroy(&sub_schurs->interior_solver);CHKERRQ(ierr);
92334a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
92434a97f8cSStefano Zampini   PetscFunctionReturn(0);
92534a97f8cSStefano Zampini }
92634a97f8cSStefano Zampini 
92734a97f8cSStefano Zampini #undef __FUNCT__
92834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
9292a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
93034a97f8cSStefano Zampini {
93134a97f8cSStefano Zampini   PetscInt       i,j,n;
93234a97f8cSStefano Zampini   PetscErrorCode ierr;
93334a97f8cSStefano Zampini 
93434a97f8cSStefano Zampini   PetscFunctionBegin;
93534a97f8cSStefano Zampini   n = 0;
93634a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
93734a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
93834a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
93934a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
94034a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
94134a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
94234a97f8cSStefano Zampini         queue_tip[n] = dof;
94334a97f8cSStefano Zampini         n++;
94434a97f8cSStefano Zampini       }
94534a97f8cSStefano Zampini     }
94634a97f8cSStefano Zampini   }
94734a97f8cSStefano Zampini   *n_added = n;
94834a97f8cSStefano Zampini   PetscFunctionReturn(0);
94934a97f8cSStefano Zampini }
950