xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision e28d306cefdfc9732ab6ae67ada7c2d6ef2a32cd)
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 PCBDDCMumpsInteriorSolve(PC,Vec,Vec);
8d62866d3SStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC,Vec,Vec);
9d62866d3SStefano Zampini 
10d62866d3SStefano Zampini #undef __FUNCT__
11*e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve_Private"
12*e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose)
13d62866d3SStefano Zampini {
14d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
15d62866d3SStefano Zampini   PetscInt         ival;
16d62866d3SStefano Zampini   PetscErrorCode   ierr;
17d62866d3SStefano Zampini 
18d62866d3SStefano Zampini   PetscFunctionBegin;
19d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
20d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
21d62866d3SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
22d62866d3SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
23d62866d3SStefano Zampini #endif
24*e28d306cSStefano Zampini   if (transpose) {
25*e28d306cSStefano Zampini     ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
26*e28d306cSStefano Zampini   } else {
276816873aSStefano Zampini     ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
28*e28d306cSStefano Zampini   }
29d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
30d62866d3SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
31d62866d3SStefano Zampini #endif
32d62866d3SStefano Zampini   PetscFunctionReturn(0);
33d62866d3SStefano Zampini }
34d62866d3SStefano Zampini 
35d62866d3SStefano Zampini #undef __FUNCT__
36*e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
37*e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
38*e28d306cSStefano Zampini {
39*e28d306cSStefano Zampini   PetscErrorCode   ierr;
40*e28d306cSStefano Zampini 
41*e28d306cSStefano Zampini   PetscFunctionBegin;
42*e28d306cSStefano Zampini   ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr);
43*e28d306cSStefano Zampini   PetscFunctionReturn(0);
44*e28d306cSStefano Zampini }
45*e28d306cSStefano Zampini 
46*e28d306cSStefano Zampini #undef __FUNCT__
47*e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose"
48*e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol)
49*e28d306cSStefano Zampini {
50*e28d306cSStefano Zampini   PetscErrorCode   ierr;
51*e28d306cSStefano Zampini 
52*e28d306cSStefano Zampini   PetscFunctionBegin;
53*e28d306cSStefano Zampini   ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr);
54*e28d306cSStefano Zampini   PetscFunctionReturn(0);
55*e28d306cSStefano Zampini }
56*e28d306cSStefano Zampini 
57*e28d306cSStefano Zampini #undef __FUNCT__
58d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCReuseMumpsReset"
59d62866d3SStefano Zampini static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
60d62866d3SStefano Zampini {
61d62866d3SStefano Zampini   PetscErrorCode ierr;
62d62866d3SStefano Zampini 
63d62866d3SStefano Zampini   PetscFunctionBegin;
64d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
65d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
66d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
67d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
68d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
69d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
706816873aSStefano Zampini   ierr = VecDestroy(&reuse->solB);CHKERRQ(ierr);
716816873aSStefano Zampini   ierr = VecDestroy(&reuse->rhsB);CHKERRQ(ierr);
72d62866d3SStefano Zampini   PetscFunctionReturn(0);
73d62866d3SStefano Zampini }
74d5574798SStefano Zampini 
75d5574798SStefano Zampini #undef __FUNCT__
76*e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve_Private"
77*e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose)
78d5574798SStefano Zampini {
79d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
80d5574798SStefano Zampini   PetscScalar      *array,*array_mumps;
81d5574798SStefano Zampini   PetscInt         ival;
82d5574798SStefano Zampini   PetscErrorCode   ierr;
83d5574798SStefano Zampini 
84d5574798SStefano Zampini   PetscFunctionBegin;
85d5574798SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
86d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
87d5574798SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
88d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
89d62866d3SStefano Zampini #endif
90d5574798SStefano Zampini   /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */
91d5574798SStefano Zampini   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
92d5574798SStefano Zampini   ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
93d62866d3SStefano Zampini   ierr = PetscMemcpy(array_mumps,array,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
94d5574798SStefano Zampini   ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
95d5574798SStefano Zampini   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
96d5574798SStefano Zampini 
97*e28d306cSStefano Zampini   if (transpose) {
98*e28d306cSStefano Zampini     ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
99*e28d306cSStefano Zampini   } else {
100d5574798SStefano Zampini     ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
101*e28d306cSStefano Zampini   }
102d5574798SStefano Zampini 
103d5574798SStefano Zampini   /* get back data to caller worskpace */
104d5574798SStefano Zampini   ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
105d5574798SStefano Zampini   ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
106d62866d3SStefano Zampini   ierr = PetscMemcpy(array,array_mumps,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
107d5574798SStefano Zampini   ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
108d5574798SStefano Zampini   ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
109d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
110d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
111d62866d3SStefano Zampini #endif
112d5574798SStefano Zampini   PetscFunctionReturn(0);
113d5574798SStefano Zampini }
1143202ece2SStefano Zampini 
1153202ece2SStefano Zampini #undef __FUNCT__
116*e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
117*e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
118*e28d306cSStefano Zampini {
119*e28d306cSStefano Zampini   PetscErrorCode   ierr;
120*e28d306cSStefano Zampini 
121*e28d306cSStefano Zampini   PetscFunctionBegin;
122*e28d306cSStefano Zampini   ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr);
123*e28d306cSStefano Zampini   PetscFunctionReturn(0);
124*e28d306cSStefano Zampini }
125*e28d306cSStefano Zampini 
126*e28d306cSStefano Zampini #undef __FUNCT__
127*e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose"
128*e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol)
129*e28d306cSStefano Zampini {
130*e28d306cSStefano Zampini   PetscErrorCode   ierr;
131*e28d306cSStefano Zampini 
132*e28d306cSStefano Zampini   PetscFunctionBegin;
133*e28d306cSStefano Zampini   ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr);
134*e28d306cSStefano Zampini   PetscFunctionReturn(0);
135*e28d306cSStefano Zampini }
136*e28d306cSStefano Zampini 
137*e28d306cSStefano Zampini #undef __FUNCT__
1383202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
1395ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
1403202ece2SStefano Zampini {
1413202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
1423202ece2SStefano Zampini   KSP            ksp;
1433202ece2SStefano Zampini   PC             pc;
1443202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
1453202ece2SStefano Zampini   PetscReal      fill = 2.0;
146f11841e3SStefano Zampini   PetscInt       n_I;
1473202ece2SStefano Zampini   PetscMPIInt    size;
1483202ece2SStefano Zampini   PetscErrorCode ierr;
1493202ece2SStefano Zampini 
1503202ece2SStefano Zampini   PetscFunctionBegin;
1513202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
1523202ece2SStefano Zampini   if (size != 1) {
1533202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
1543202ece2SStefano Zampini   }
155f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
156f11841e3SStefano Zampini     PetscBool Sdense;
157f11841e3SStefano Zampini 
158f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
159f11841e3SStefano Zampini     if (!Sdense) {
160f11841e3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
161f11841e3SStefano Zampini     }
162f11841e3SStefano Zampini   }
1633202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
1643202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
1653202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1663202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
1673202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
1683202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
1693202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
1703202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
171f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
172f11841e3SStefano Zampini   if (n_I) {
1733202ece2SStefano Zampini     if (!Bdense) {
1743202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
1753202ece2SStefano Zampini     } else {
1763202ece2SStefano Zampini       Bd = B;
1773202ece2SStefano Zampini     }
1783202ece2SStefano Zampini 
1793202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1803202ece2SStefano Zampini       Mat fact;
1813202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1823202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1833202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1843202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1853202ece2SStefano Zampini     } else {
18607b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
18707b1e237SStefano Zampini 
18807b1e237SStefano Zampini       if (ex) {
1893202ece2SStefano Zampini         Mat Ainvd;
1903202ece2SStefano Zampini 
1913202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1923202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1933202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
19407b1e237SStefano Zampini       } else {
19507b1e237SStefano Zampini         Vec         sol,rhs;
19607b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
19707b1e237SStefano Zampini         PetscInt    i,nrhs,n;
19807b1e237SStefano Zampini 
19907b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
20007b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
20107b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
20207b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
20307b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
20407b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
20507b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
20607b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
20707b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
20807b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
20907b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
21007b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
21107b1e237SStefano Zampini         }
21207b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
21307b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
21407b1e237SStefano Zampini       }
2153202ece2SStefano Zampini     }
2165ec10c6aSStefano Zampini     if (!Bdense & !issym) {
2173202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2183202ece2SStefano Zampini     }
2195ec10c6aSStefano Zampini 
2205ec10c6aSStefano Zampini     if (!issym) {
2213202ece2SStefano Zampini       if (!Cdense) {
2223202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
2233202ece2SStefano Zampini       } else {
2243202ece2SStefano Zampini         Cd = C;
2253202ece2SStefano Zampini       }
2265ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2273202ece2SStefano Zampini       if (!Cdense) {
2283202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
2293202ece2SStefano Zampini       }
2305ec10c6aSStefano Zampini     } else {
2315ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2325ec10c6aSStefano Zampini       if (!Bdense) {
2335ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2345ec10c6aSStefano Zampini       }
2355ec10c6aSStefano Zampini     }
2365ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
237f11841e3SStefano Zampini   }
2383202ece2SStefano Zampini 
2393202ece2SStefano Zampini   if (D) {
2403202ece2SStefano Zampini     Mat       Dd;
2413202ece2SStefano Zampini     PetscBool Ddense;
2423202ece2SStefano Zampini 
2433202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
2443202ece2SStefano Zampini     if (!Ddense) {
2453202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
2463202ece2SStefano Zampini     } else {
2473202ece2SStefano Zampini       Dd = D;
2483202ece2SStefano Zampini     }
249f11841e3SStefano Zampini     if (n_I) {
2503202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
251f11841e3SStefano Zampini     } else {
252f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
253f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
254f11841e3SStefano Zampini       } else {
255f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
256f11841e3SStefano Zampini       }
257f11841e3SStefano Zampini     }
2583202ece2SStefano Zampini     if (!Ddense) {
2593202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
2603202ece2SStefano Zampini     }
2613202ece2SStefano Zampini   } else {
2623202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
2633202ece2SStefano Zampini   }
2643202ece2SStefano Zampini   PetscFunctionReturn(0);
2653202ece2SStefano Zampini }
26634a97f8cSStefano Zampini 
26734a97f8cSStefano Zampini #undef __FUNCT__
2681580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
2696816873aSStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers,PetscBool use_edges, PetscBool use_faces)
270b1b3d7a2SStefano Zampini {
271be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
272be83ff47SStefano Zampini   Mat                    S_all;
273d2627357SStefano Zampini   Mat                    global_schur_subsets,work_mat;
27408122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
2755db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
276a9b99552SStefano Zampini   IS                     is_I,is_I_layer,temp_is;
277d648f858SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
278eb595f79SStefano Zampini   PetscInt               *auxnum1,*auxnum2,*all_local_idx_G_rep;
2795a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
280883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
28108122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
28206a4b1faSStefano Zampini   PetscScalar            *Bwork;
2835a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2845a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2855a95e1ceSStefano Zampini   MPI_Comm               comm_n;
286b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
287b1b3d7a2SStefano Zampini 
288b1b3d7a2SStefano Zampini   PetscFunctionBegin;
289a64f4aa4SStefano Zampini   /* update info in sub_schurs */
290a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
291a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
292a64f4aa4SStefano Zampini   if (Ain) {
293a64f4aa4SStefano Zampini     PetscBool isseqaij;
294a64f4aa4SStefano Zampini 
295a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
296a64f4aa4SStefano Zampini     if (isseqaij) {
297a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
298a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
299a64f4aa4SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
300a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
301a64f4aa4SStefano Zampini     }
302a64f4aa4SStefano Zampini   }
303a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
304a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
305a64f4aa4SStefano Zampini   if (sub_schurs->use_mumps) {
306a64f4aa4SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
307a64f4aa4SStefano Zampini   }
308a64f4aa4SStefano Zampini 
3095a95e1ceSStefano Zampini   /* preliminary checks */
3105a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
3115a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
3125a95e1ceSStefano Zampini   }
3135a95e1ceSStefano Zampini   /* determine if we are dealing with hermitian positive definite problems */
3145a95e1ceSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
3155a95e1ceSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
3165a95e1ceSStefano Zampini   if (sub_schurs->A) {
3175a95e1ceSStefano Zampini     PetscInt lsize;
3185a95e1ceSStefano Zampini 
3195a95e1ceSStefano Zampini     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
3205a95e1ceSStefano Zampini     if (lsize) {
3215a95e1ceSStefano Zampini       ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
3225a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
3235a95e1ceSStefano Zampini         PetscScalar val;
3245a95e1ceSStefano Zampini         Vec         vec1,vec2;
3255a95e1ceSStefano Zampini 
3265a95e1ceSStefano Zampini         ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
3275a95e1ceSStefano Zampini         ierr = VecSetRandom(vec1,NULL);
3285a95e1ceSStefano Zampini         ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
3295a95e1ceSStefano Zampini         ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
3305a95e1ceSStefano Zampini         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
3315a95e1ceSStefano Zampini         if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
3325a95e1ceSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
3335a95e1ceSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
3345a95e1ceSStefano Zampini       }
3355a95e1ceSStefano Zampini     } else {
3365a95e1ceSStefano Zampini       sub_schurs->is_hermitian = PETSC_TRUE;
3375a95e1ceSStefano Zampini       sub_schurs->is_posdef = PETSC_TRUE;
3385a95e1ceSStefano Zampini     }
3395a95e1ceSStefano Zampini     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
3405a95e1ceSStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported");
3415a95e1ceSStefano Zampini     }
3425a95e1ceSStefano Zampini   }
3435a95e1ceSStefano Zampini   /* restrict work on active processes */
3445a95e1ceSStefano Zampini   color = 0;
3455a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
3465a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
3475a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
3485a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
3495a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
3505a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
3515a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
3525a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3535a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
3545a95e1ceSStefano Zampini     PetscFunctionReturn(0);
3555a95e1ceSStefano Zampini   }
3565a95e1ceSStefano Zampini 
357b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
358883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
359a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
360f11841e3SStefano Zampini     PetscBool isseqaij;
361a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
362a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
363f11841e3SStefano Zampini     if (!isseqaij) {
364a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
365a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
366a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
367a64f4aa4SStefano Zampini     } else {
368a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
369a64f4aa4SStefano Zampini       A_BB = tA_BB;
370a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
371a64f4aa4SStefano Zampini       A_IB = tA_IB;
372a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
373a64f4aa4SStefano Zampini       A_BI = tA_BI;
374f11841e3SStefano Zampini     }
375a58a30b4SStefano Zampini   } else {
3765a95e1ceSStefano Zampini     A_II = NULL;
3775a95e1ceSStefano Zampini     A_IB = NULL;
3785a95e1ceSStefano Zampini     A_BI = NULL;
3795a95e1ceSStefano Zampini     A_BB = NULL;
380b1b3d7a2SStefano Zampini   }
3815a95e1ceSStefano Zampini   S_all = NULL;
3825a95e1ceSStefano Zampini   S_Ej_tilda_all = NULL;
3835a95e1ceSStefano Zampini   S_Ej_inv_all = NULL;
384b1b3d7a2SStefano Zampini 
385b1b3d7a2SStefano Zampini   /* determine interior problems */
3863dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
3873dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
388b1b3d7a2SStefano Zampini     PetscBT                touched;
389b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
390b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
391b1b3d7a2SStefano Zampini 
3923dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
3933dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
3943dc780c3SStefano Zampini     }
395b1b3d7a2SStefano Zampini     /* get sizes */
396b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
397b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
398b1b3d7a2SStefano Zampini 
399b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
400b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
401b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
402b1b3d7a2SStefano Zampini 
403b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
404b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
405b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
406b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
407b1b3d7a2SStefano Zampini     }
408b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
409b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
410b1b3d7a2SStefano Zampini 
411b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
412b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
413b1b3d7a2SStefano Zampini     n_prev_added = n_B;
414b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
415b1b3d7a2SStefano Zampini       PetscInt n_added;
416b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
417b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
418b1b3d7a2SStefano 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);
419b1b3d7a2SStefano Zampini       }
420b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
421b1b3d7a2SStefano Zampini       n_prev_added = n_added;
422b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
423b1b3d7a2SStefano Zampini       if (!n_added) break;
424b1b3d7a2SStefano Zampini     }
425b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
426b1b3d7a2SStefano Zampini 
427883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
428a9b99552SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);CHKERRQ(ierr);
429b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
430a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
431883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
432883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
433b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
434b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
435a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
436b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
437b1b3d7a2SStefano Zampini 
438b1b3d7a2SStefano Zampini       /* II block */
439b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
440b1b3d7a2SStefano Zampini     }
441b1b3d7a2SStefano Zampini   } else {
442b1b3d7a2SStefano Zampini     PetscInt n_I;
443b1b3d7a2SStefano Zampini 
444b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
445b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
446a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
447b1b3d7a2SStefano Zampini 
448b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
449883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
450b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
451b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
452b1b3d7a2SStefano Zampini 
453b1b3d7a2SStefano Zampini       /* II block is the same */
454b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
455b1b3d7a2SStefano Zampini       AE_II = A_II;
456b1b3d7a2SStefano Zampini     }
457b1b3d7a2SStefano Zampini   }
4585a95e1ceSStefano Zampini 
459883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
4605a95e1ceSStefano Zampini   max_subset_size = 0;
461883469d8SStefano Zampini   local_size = 0;
4625a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4635a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4645a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
465883469d8SStefano Zampini     local_size += subset_size;
466883469d8SStefano Zampini   }
467883469d8SStefano Zampini 
468883469d8SStefano Zampini   /* Work arrays for local indices */
469883469d8SStefano Zampini   extra = 0;
470883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
471a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
472883469d8SStefano Zampini   }
473883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
474883469d8SStefano Zampini   if (extra) {
475883469d8SStefano Zampini     const PetscInt *idxs;
476a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
477883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
478a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
479883469d8SStefano Zampini   }
480883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
481eb595f79SStefano Zampini   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
482883469d8SStefano Zampini 
483883469d8SStefano Zampini   /* Get local indices in local numbering */
484883469d8SStefano Zampini   local_size = 0;
4855a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
486883469d8SStefano Zampini     PetscInt j;
487883469d8SStefano Zampini     const    PetscInt *idxs;
488883469d8SStefano Zampini 
4895a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4905a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
491eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
492eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
493eb595f79SStefano Zampini     auxnum2[i] = subset_size;
494883469d8SStefano Zampini     /* subset indices in local numbering */
495883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
4965a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
497883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
498883469d8SStefano Zampini     local_size += subset_size;
499883469d8SStefano Zampini   }
500883469d8SStefano Zampini 
5015a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
502d2627357SStefano Zampini   Bwork = NULL;
503d2627357SStefano Zampini   pivots = NULL;
5045a95e1ceSStefano Zampini   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
505d2627357SStefano Zampini     PetscScalar lwork;
506d2627357SStefano Zampini 
507d2627357SStefano Zampini     B_lwork = -1;
508d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
509d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
510d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
511d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
512d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
513d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
514d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
515d2627357SStefano Zampini   }
516d2627357SStefano Zampini 
517d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
518eb595f79SStefano Zampini   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
519eb595f79SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
520eb595f79SStefano Zampini   local_size = 0;
521eb595f79SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
522eb595f79SStefano Zampini     PetscInt j;
523eb595f79SStefano Zampini     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
524eb595f79SStefano Zampini   }
525eb595f79SStefano Zampini   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
526eb595f79SStefano Zampini   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
5275a95e1ceSStefano Zampini   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
5285a95e1ceSStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
529d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
530d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
531d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
532d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
5332972d61bSStefano Zampini 
5345a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
5355a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
5365a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
5375a95e1ceSStefano Zampini 
5385a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
5395a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
5405a95e1ceSStefano Zampini     if (subset_size != local_size) {
5415a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
5425a95e1ceSStefano Zampini     }
5435a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
544b1b3d7a2SStefano Zampini   }
545b1b3d7a2SStefano Zampini 
5465a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
5475a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
5485a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
5495a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
5505a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
5515a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
552aa83b6aeSStefano Zampini   } else {
5535a95e1ceSStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
554aa83b6aeSStefano Zampini   }
555b1b3d7a2SStefano Zampini 
5565a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
557be83ff47SStefano Zampini   F = NULL;
5585a95e1ceSStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
5595a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps) {
5605a95e1ceSStefano Zampini     Mat         S_Ej_expl;
5615a95e1ceSStefano Zampini     PetscScalar *work;
5625a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
5635a95e1ceSStefano Zampini     PetscBool   Sdense;
5645a95e1ceSStefano Zampini 
5655a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
5665a95e1ceSStefano Zampini     local_size = 0;
567b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
5685a95e1ceSStefano Zampini       IS  is_subset_B;
5695a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
5705a95e1ceSStefano Zampini 
5715a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
5725a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
5735a95e1ceSStefano Zampini       /* EE block */
5745a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
5755a95e1ceSStefano Zampini       /* IE block */
5765a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
5775a95e1ceSStefano Zampini       /* EI block */
5785a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
5795a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
5805a95e1ceSStefano Zampini       } else {
5815a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
5825a95e1ceSStefano Zampini       }
583a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
5845a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
5855a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
5865a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
5875a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
588b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
589b1b3d7a2SStefano Zampini         KSP ksp;
590b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
5915a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
592b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
593b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
594b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
595b1b3d7a2SStefano Zampini         KSPType   ksp_type;
596b1b3d7a2SStefano Zampini         PetscInt  n_internal;
5975a95e1ceSStefano Zampini         PetscBool ispcnone;
598b1b3d7a2SStefano Zampini 
599b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
6005a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
601b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
602b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
603b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
604b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
6055a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
6065a95e1ceSStefano Zampini         if (!ispcnone) {
6075a95e1ceSStefano Zampini           PCType pc_type;
608b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
609b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
6105a95e1ceSStefano Zampini         } else {
6115a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
6125a95e1ceSStefano Zampini         }
613b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
614b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
615b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
616b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
617b1b3d7a2SStefano Zampini           if (solver) {
618b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
619b1b3d7a2SStefano Zampini           }
620b1b3d7a2SStefano Zampini         }
621b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
622b1b3d7a2SStefano Zampini       }
6235a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6245a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
6255a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
6265a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
6275a95e1ceSStefano Zampini       if (Sdense) {
6285a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
6295a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
630b1b3d7a2SStefano Zampini         }
6315a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6325a95e1ceSStefano Zampini       } else {
6335a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
6345a95e1ceSStefano Zampini       }
6355a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
636a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
6375a95e1ceSStefano Zampini       local_size += subset_size;
6385a95e1ceSStefano Zampini     }
6395a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
640b1b3d7a2SStefano Zampini     /* free */
641b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
642b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
6435a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
644883469d8SStefano Zampini   } else {
645be83ff47SStefano Zampini     Mat         A;
646883469d8SStefano Zampini     IS          is_A_all;
647be83ff47SStefano Zampini     PetscScalar *work,*S_data;
648be83ff47SStefano Zampini     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2;
649be83ff47SStefano Zampini     PetscBool   restore_S;
650883469d8SStefano Zampini 
651883469d8SStefano Zampini     /* get working mat */
652a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
653d62866d3SStefano Zampini     if (!sub_schurs->is_dir) {
654883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
655d62866d3SStefano Zampini       size_schur = local_size;
656d62866d3SStefano Zampini     } else {
657d62866d3SStefano Zampini       IS list[2];
658d62866d3SStefano Zampini 
659d62866d3SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr);
660d62866d3SStefano Zampini       list[1] = sub_schurs->is_dir;
661d62866d3SStefano Zampini       ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr);
662d62866d3SStefano Zampini       ierr = ISDestroy(&list[0]);CHKERRQ(ierr);
663d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr);
664d62866d3SStefano Zampini       size_schur += local_size;
665d62866d3SStefano Zampini     }
666d62866d3SStefano Zampini     size_active_schur = local_size;
6676816873aSStefano Zampini     ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
668a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
669883469d8SStefano Zampini 
67008122e43SStefano Zampini     if (n_I) {
6719ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
672d62866d3SStefano Zampini         ierr = MatSetOption(A,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
673883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
674883469d8SStefano Zampini       } else {
675d62866d3SStefano Zampini         ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
676883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
677883469d8SStefano Zampini       }
678883469d8SStefano Zampini 
679883469d8SStefano Zampini       /* subsets ordered last */
680d62866d3SStefano Zampini       ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr);
681d62866d3SStefano Zampini       for (i=0;i<size_schur;i++) {
682883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
683883469d8SStefano Zampini       }
6845a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
685d62866d3SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr);
6865a95e1ceSStefano Zampini #endif
687883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
688883469d8SStefano Zampini 
689883469d8SStefano Zampini       /* factorization step */
6909ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
691883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
692be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
693be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
694be83ff47SStefano Zampini #endif
695883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
696883469d8SStefano Zampini       } else {
697883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
698be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
699be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
700be83ff47SStefano Zampini #endif
701883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
702883469d8SStefano Zampini       }
703883469d8SStefano Zampini 
704883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
7055a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
706883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
7075a95e1ceSStefano Zampini #endif
708d5574798SStefano Zampini 
709d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
710d5574798SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
711be83ff47SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all));
712be83ff47SStefano Zampini       restore_S = PETSC_TRUE;
713be83ff47SStefano Zampini     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
714be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
715be83ff47SStefano Zampini       reuse_solvers = PETSC_FALSE;
716be83ff47SStefano Zampini       restore_S = PETSC_FALSE;
717be83ff47SStefano Zampini     }
718be83ff47SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
719be83ff47SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
720be83ff47SStefano Zampini 
721be83ff47SStefano Zampini     if (reuse_solvers) {
722d62866d3SStefano Zampini       Mat              A_II;
723d62866d3SStefano Zampini       PCBDDCReuseMumps msolv_ctx;
724d5574798SStefano Zampini 
725d62866d3SStefano Zampini       if (sub_schurs->reuse_mumps) {
7266816873aSStefano Zampini         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
727*e28d306cSStefano Zampini       } else {
728*e28d306cSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
729d62866d3SStefano Zampini       }
730*e28d306cSStefano Zampini       msolv_ctx = sub_schurs->reuse_mumps;
731be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
732be83ff47SStefano Zampini       ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr);
733d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
734d5574798SStefano Zampini       msolv_ctx->F = F;
735d5574798SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
736d62866d3SStefano Zampini 
737d62866d3SStefano Zampini       /* interior solver */
738d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr);
739d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
740d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
741d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
742d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
743*e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
744d62866d3SStefano Zampini 
745d62866d3SStefano Zampini       /* correction solver */
746d62866d3SStefano Zampini       /* auxiliary scatters are needed and are created in PCBDDCSetUpLocalScatters */
747d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr);
748d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
749d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
750d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
751d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
752*e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
753be83ff47SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->solB,&msolv_ctx->rhsB);CHKERRQ(ierr);
754d5574798SStefano Zampini     }
75508122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
7565db18549SStefano Zampini 
757be83ff47SStefano Zampini     /* Work arrays */
758be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
75912d906b1SStefano Zampini     if (compute_Stilda) {
760a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
761d62866d3SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
762a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
763a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
76408122e43SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
765d62866d3SStefano Zampini       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
76608122e43SStefano Zampini       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
76708122e43SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
76812d906b1SStefano Zampini     }
769d2627357SStefano Zampini 
770be83ff47SStefano Zampini     /* S_Ej_all */
771be83ff47SStefano Zampini     cum = cum2 = 0;
772be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
77365d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
77465d8bf0aSStefano Zampini       PetscInt j;
77565d8bf0aSStefano Zampini 
7765a95e1ceSStefano Zampini       /* get S_E */
777b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
778be83ff47SStefano Zampini       if (restore_S && sub_schurs->is_hermitian) { /* When using MUMPS data I need to expand to upper triangular (column oriented) */
779be83ff47SStefano Zampini         PetscInt k;
780be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
781be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
782be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
783be83ff47SStefano Zampini             work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
784be83ff47SStefano Zampini           }
785be83ff47SStefano Zampini         }
786be83ff47SStefano Zampini       } else { /* copy to workspace */
787be83ff47SStefano Zampini         PetscInt k;
788be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
789be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
790be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
791be83ff47SStefano Zampini           }
792be83ff47SStefano Zampini         }
7939087bf02SStefano Zampini       }
7945a95e1ceSStefano Zampini       /* insert S_E values */
795be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
7965a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
797a1337663SStefano Zampini 
798be83ff47SStefano Zampini       /* if adaptivity is requested, invert S_E block */
799d2627357SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
80008122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
80108122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
802be83ff47SStefano Zampini         if (!sub_schurs->is_hermitian) { /* TODO add sytrf/i for symmetric non hermitian */
8035a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
80408122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
8055a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
80608122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
8072972d61bSStefano Zampini         } else {
8085a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
8092972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
8105a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
8112972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
8122972d61bSStefano Zampini         }
81308122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
8145a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8159087bf02SStefano Zampini       }
816be83ff47SStefano Zampini       cum += subset_size;
817be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
818be83ff47SStefano Zampini     }
819be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
820be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
821be83ff47SStefano Zampini     if (restore_S) {
822be83ff47SStefano Zampini       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
823be83ff47SStefano Zampini       ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr);
824be83ff47SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
825be83ff47SStefano Zampini     } else if (compute_Stilda) { /* we need to invert explicitly since we are not using MUMPS */
826be83ff47SStefano Zampini       ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
827be83ff47SStefano Zampini       ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
828be83ff47SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
829be83ff47SStefano Zampini       if (!sub_schurs->is_hermitian) { /* TODO add sytrf/i for symmetric non hermitian */
830be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
831be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
832be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
833be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
834be83ff47SStefano Zampini       } else {
835be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
836be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
837be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
838be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
839be83ff47SStefano Zampini       }
840be83ff47SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
841be83ff47SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
842be83ff47SStefano Zampini     }
843be83ff47SStefano Zampini #endif
844be83ff47SStefano Zampini     /* S_Ej_tilda_all */
845be83ff47SStefano Zampini     cum = cum2 = 0;
846be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
847be83ff47SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
848be83ff47SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
849be83ff47SStefano Zampini       /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
850be83ff47SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
851be83ff47SStefano Zampini         PetscInt j;
852be83ff47SStefano Zampini         /* get (St^-1)_E */
853be83ff47SStefano Zampini         if (sub_schurs->is_hermitian) { /* I need to expand to upper triangular (column oriented) */
854be83ff47SStefano Zampini           PetscInt k;
855be83ff47SStefano Zampini           for (k=0;k<subset_size;k++) {
856be83ff47SStefano Zampini             for (j=k;j<subset_size;j++) {
857be83ff47SStefano Zampini               work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
858be83ff47SStefano Zampini               work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
859be83ff47SStefano Zampini             }
860be83ff47SStefano Zampini           }
861be83ff47SStefano Zampini         } else { /* copy to workspace */
862be83ff47SStefano Zampini           PetscInt k;
863be83ff47SStefano Zampini           for (k=0;k<subset_size;k++) {
864be83ff47SStefano Zampini             for (j=0;j<subset_size;j++) {
865be83ff47SStefano Zampini               work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
866be83ff47SStefano Zampini             }
867be83ff47SStefano Zampini           }
868be83ff47SStefano Zampini         }
869be83ff47SStefano Zampini         for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
8705a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8719087bf02SStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
87208122e43SStefano Zampini       }
873be83ff47SStefano Zampini       cum += subset_size;
874be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
875883469d8SStefano Zampini     }
876be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
8775ec10c6aSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
878be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
879be83ff47SStefano Zampini     if (restore_S) {
880be83ff47SStefano Zampini       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
8815db18549SStefano Zampini     }
882be83ff47SStefano Zampini #endif
8839087bf02SStefano Zampini   }
884be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
885a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
886a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
887a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
888a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
889a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
890a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
8915db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8925db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8935a95e1ceSStefano Zampini   if (compute_Stilda) {
894a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
895a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89608122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89708122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89808122e43SStefano Zampini   }
899a1337663SStefano Zampini 
9005db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
9015db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
9023927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
9033927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
9045a95e1ceSStefano Zampini 
9055db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
9065a95e1ceSStefano Zampini   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
907d648f858SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
908d648f858SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
90908122e43SStefano Zampini 
910ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
911ac632422SStefano Zampini   if (faster_deluxe) {
9125a95e1ceSStefano Zampini     Mat         tmpmat;
9135a95e1ceSStefano Zampini     PetscScalar *array;
9145a95e1ceSStefano Zampini     PetscInt    cum;
9155a95e1ceSStefano Zampini 
9165a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
9175a95e1ceSStefano Zampini     cum = 0;
9185a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
9195a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
9205a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
9215a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
9225a95e1ceSStefano Zampini       if (!sub_schurs->is_hermitian) {
9235a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
9245a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
9255a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
9265a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
9275a95e1ceSStefano Zampini       } else {
9285a95e1ceSStefano Zampini         PetscInt j,k;
9295a95e1ceSStefano Zampini 
9305a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
9315a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
9325a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
9335a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
9345a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
9355a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
9365a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
9375a95e1ceSStefano Zampini           }
9385a95e1ceSStefano Zampini         }
9395a95e1ceSStefano Zampini       }
9405a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
9415a95e1ceSStefano Zampini       cum += subset_size*subset_size;
9425a95e1ceSStefano Zampini     }
9435a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
9445a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
9455a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
946ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
9475a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
9485a95e1ceSStefano Zampini   }
9495a95e1ceSStefano Zampini 
950f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
9515a95e1ceSStefano Zampini   if (compute_Stilda) {
952a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
953a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
954d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
955d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
95608122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
95708122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
958d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
959d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
96008122e43SStefano Zampini   }
9613202ece2SStefano Zampini 
9625a95e1ceSStefano Zampini   /* free workspace */
96306a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
964a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
965a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
96608122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
9673202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
9685db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
9695a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
970b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
971b1b3d7a2SStefano Zampini }
972b1b3d7a2SStefano Zampini 
973b1b3d7a2SStefano Zampini #undef __FUNCT__
974b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
975a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
976b1b3d7a2SStefano Zampini {
9779bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
9785a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
979b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
980b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
981b1b3d7a2SStefano Zampini 
982b1b3d7a2SStefano Zampini   PetscFunctionBegin;
983b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
984b1b3d7a2SStefano Zampini   if (!is_sorted) {
985b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
986b1b3d7a2SStefano Zampini   }
987b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
988b1b3d7a2SStefano Zampini   if (!is_sorted) {
989b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
990b1b3d7a2SStefano Zampini   }
991b1b3d7a2SStefano Zampini 
992b1b3d7a2SStefano Zampini   /* reset any previous data */
993b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
994b1b3d7a2SStefano Zampini 
9955a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
9969bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
997b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
99808122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
99908122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
1000b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1001b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
1002b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
1003b1b3d7a2SStefano Zampini   }
1004b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
1005b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
100608122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1007b1b3d7a2SStefano Zampini   }
1008b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
1009b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
1010d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1011d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1012b1b3d7a2SStefano Zampini 
10136816873aSStefano Zampini   /* Determine if MUMPS can be used */
1014883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
1015883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1016a64f4aa4SStefano Zampini   sub_schurs->use_mumps = PETSC_TRUE;
1017883469d8SStefano Zampini #endif
1018b1b3d7a2SStefano Zampini 
1019b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1020b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1021b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1022b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
10235db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
10245db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
10255db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
10265db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
10275a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1028b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1029a64f4aa4SStefano Zampini   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1030b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1031b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1032b96c3477SStefano Zampini     }
10339bb4a8caSStefano Zampini   }
1034d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
1035b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1036b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
103708122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1038b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1039b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1040b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1041b1b3d7a2SStefano Zampini }
1042b1b3d7a2SStefano Zampini 
1043b1b3d7a2SStefano Zampini #undef __FUNCT__
104434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
104534a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
104634a97f8cSStefano Zampini {
104734a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
104834a97f8cSStefano Zampini   PetscErrorCode  ierr;
104934a97f8cSStefano Zampini 
105034a97f8cSStefano Zampini   PetscFunctionBegin;
105134a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
10525ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
105334a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
105434a97f8cSStefano Zampini   PetscFunctionReturn(0);
105534a97f8cSStefano Zampini }
105634a97f8cSStefano Zampini 
105734a97f8cSStefano Zampini #undef __FUNCT__
105834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
105934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
106034a97f8cSStefano Zampini {
106134a97f8cSStefano Zampini   PetscErrorCode ierr;
106234a97f8cSStefano Zampini 
106334a97f8cSStefano Zampini   PetscFunctionBegin;
106434a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
106534a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
106634a97f8cSStefano Zampini   PetscFunctionReturn(0);
106734a97f8cSStefano Zampini }
106834a97f8cSStefano Zampini 
106934a97f8cSStefano Zampini #undef __FUNCT__
107034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
107134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
107234a97f8cSStefano Zampini {
107334a97f8cSStefano Zampini   PetscInt       i;
107434a97f8cSStefano Zampini   PetscErrorCode ierr;
107534a97f8cSStefano Zampini 
107634a97f8cSStefano Zampini   PetscFunctionBegin;
10771e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1078b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1079b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1080b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
10815db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
10825db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
108341c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
108441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
108508122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1086a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
10875db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1088d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1089d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
109008122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
109108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
109234a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1093b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
109434a97f8cSStefano Zampini   }
10955ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1096b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
10973dc780c3SStefano Zampini   }
1098d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps) {
1099d62866d3SStefano Zampini     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1100d62866d3SStefano Zampini   }
1101d62866d3SStefano Zampini   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
110234a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
110334a97f8cSStefano Zampini   PetscFunctionReturn(0);
110434a97f8cSStefano Zampini }
110534a97f8cSStefano Zampini 
110634a97f8cSStefano Zampini #undef __FUNCT__
110734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
11082a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
110934a97f8cSStefano Zampini {
111034a97f8cSStefano Zampini   PetscInt       i,j,n;
111134a97f8cSStefano Zampini   PetscErrorCode ierr;
111234a97f8cSStefano Zampini 
111334a97f8cSStefano Zampini   PetscFunctionBegin;
111434a97f8cSStefano Zampini   n = 0;
111534a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
111634a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
111734a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
111834a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
111934a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
112034a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
112134a97f8cSStefano Zampini         queue_tip[n] = dof;
112234a97f8cSStefano Zampini         n++;
112334a97f8cSStefano Zampini       }
112434a97f8cSStefano Zampini     }
112534a97f8cSStefano Zampini   }
112634a97f8cSStefano Zampini   *n_added = n;
112734a97f8cSStefano Zampini   PetscFunctionReturn(0);
112834a97f8cSStefano Zampini }
1129