xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 6dba178d53d088f8af27c39b86adaa0496f0d2e2)
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__
11683d3df6SStefano Zampini #define __FUNCT__ "PCBDDCMumpsReuseSolve_Private"
12683d3df6SStefano Zampini static PetscErrorCode PCBDDCMumpsReuseSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
13d62866d3SStefano Zampini {
14d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
15683d3df6SStefano Zampini   PetscScalar      *array,*array_mumps;
16683d3df6SStefano Zampini   PetscInt         n;
17683d3df6SStefano Zampini   PetscBool        copy = PETSC_FALSE;
18d62866d3SStefano Zampini   PetscErrorCode   ierr;
19d62866d3SStefano Zampini 
20d62866d3SStefano Zampini   PetscFunctionBegin;
21d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
22683d3df6SStefano Zampini   if (full) {
23d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
24d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
25d62866d3SStefano Zampini #endif
26683d3df6SStefano Zampini     copy = ctx->has_vertices;
27683d3df6SStefano Zampini   } else {
28*6dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
29683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
30*6dba178dSStefano Zampini #endif
31683d3df6SStefano Zampini     copy = PETSC_TRUE;
32683d3df6SStefano Zampini   }
33683d3df6SStefano Zampini   ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
34683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
35683d3df6SStefano Zampini   if (copy) {
36683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
37683d3df6SStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
38683d3df6SStefano Zampini     ierr = PetscMemcpy(array_mumps,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
39683d3df6SStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
40683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
41683d3df6SStefano Zampini 
42683d3df6SStefano Zampini     if (transpose) {
43683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
44683d3df6SStefano Zampini     } else {
45683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
46683d3df6SStefano Zampini     }
47683d3df6SStefano Zampini 
48683d3df6SStefano Zampini     /* get back data to caller worskpace */
49683d3df6SStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
50683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
51683d3df6SStefano Zampini     ierr = PetscMemcpy(array,array_mumps,n*sizeof(PetscScalar));CHKERRQ(ierr);
52683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
53683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
54683d3df6SStefano Zampini   } else {
55e28d306cSStefano Zampini     if (transpose) {
56e28d306cSStefano Zampini       ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
57e28d306cSStefano Zampini     } else {
586816873aSStefano Zampini       ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
59e28d306cSStefano Zampini     }
60683d3df6SStefano Zampini   }
61d62866d3SStefano Zampini   PetscFunctionReturn(0);
62d62866d3SStefano Zampini }
63d62866d3SStefano Zampini 
64d62866d3SStefano Zampini #undef __FUNCT__
65e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
66e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
67e28d306cSStefano Zampini {
68e28d306cSStefano Zampini   PetscErrorCode   ierr;
69e28d306cSStefano Zampini 
70e28d306cSStefano Zampini   PetscFunctionBegin;
71683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
72e28d306cSStefano Zampini   PetscFunctionReturn(0);
73e28d306cSStefano Zampini }
74e28d306cSStefano Zampini 
75e28d306cSStefano Zampini #undef __FUNCT__
76e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose"
77e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol)
78e28d306cSStefano Zampini {
79e28d306cSStefano Zampini   PetscErrorCode   ierr;
80e28d306cSStefano Zampini 
81e28d306cSStefano Zampini   PetscFunctionBegin;
82683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
83683d3df6SStefano Zampini   PetscFunctionReturn(0);
84683d3df6SStefano Zampini }
85683d3df6SStefano Zampini 
86683d3df6SStefano Zampini #undef __FUNCT__
87683d3df6SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
88683d3df6SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
89683d3df6SStefano Zampini {
90683d3df6SStefano Zampini   PetscErrorCode   ierr;
91683d3df6SStefano Zampini 
92683d3df6SStefano Zampini   PetscFunctionBegin;
93683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
94683d3df6SStefano Zampini   PetscFunctionReturn(0);
95683d3df6SStefano Zampini }
96683d3df6SStefano Zampini 
97683d3df6SStefano Zampini #undef __FUNCT__
98683d3df6SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose"
99683d3df6SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol)
100683d3df6SStefano Zampini {
101683d3df6SStefano Zampini   PetscErrorCode   ierr;
102683d3df6SStefano Zampini 
103683d3df6SStefano Zampini   PetscFunctionBegin;
104683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
105e28d306cSStefano Zampini   PetscFunctionReturn(0);
106e28d306cSStefano Zampini }
107e28d306cSStefano Zampini 
108e28d306cSStefano Zampini #undef __FUNCT__
109d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCReuseMumpsReset"
110d62866d3SStefano Zampini static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
111d62866d3SStefano Zampini {
112d62866d3SStefano Zampini   PetscErrorCode ierr;
113d62866d3SStefano Zampini 
114d62866d3SStefano Zampini   PetscFunctionBegin;
115d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
116d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
117d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
118d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
119d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
12053892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
12153892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
122d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
12353892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
12453892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
125d62866d3SStefano Zampini   PetscFunctionReturn(0);
126d62866d3SStefano Zampini }
127d5574798SStefano Zampini 
128d5574798SStefano Zampini #undef __FUNCT__
1293202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
1305ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
1313202ece2SStefano Zampini {
1323202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
1333202ece2SStefano Zampini   KSP            ksp;
1343202ece2SStefano Zampini   PC             pc;
1353202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
1363202ece2SStefano Zampini   PetscReal      fill = 2.0;
137f11841e3SStefano Zampini   PetscInt       n_I;
1383202ece2SStefano Zampini   PetscMPIInt    size;
1393202ece2SStefano Zampini   PetscErrorCode ierr;
1403202ece2SStefano Zampini 
1413202ece2SStefano Zampini   PetscFunctionBegin;
1423202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
1433202ece2SStefano Zampini   if (size != 1) {
1443202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
1453202ece2SStefano Zampini   }
146f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
147f11841e3SStefano Zampini     PetscBool Sdense;
148f11841e3SStefano Zampini 
149f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
150f11841e3SStefano Zampini     if (!Sdense) {
151683d3df6SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should be dense");
152f11841e3SStefano Zampini     }
153f11841e3SStefano Zampini   }
1543202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
1553202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
1563202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1573202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
1583202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
1593202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
1603202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
1613202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
162f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
163f11841e3SStefano Zampini   if (n_I) {
1643202ece2SStefano Zampini     if (!Bdense) {
1653202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
1663202ece2SStefano Zampini     } else {
1673202ece2SStefano Zampini       Bd = B;
1683202ece2SStefano Zampini     }
1693202ece2SStefano Zampini 
1703202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1713202ece2SStefano Zampini       Mat fact;
1723202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1733202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1743202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1753202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1763202ece2SStefano Zampini     } else {
17707b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
17807b1e237SStefano Zampini 
17907b1e237SStefano Zampini       if (ex) {
1803202ece2SStefano Zampini         Mat Ainvd;
1813202ece2SStefano Zampini 
1823202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1833202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1843202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
18507b1e237SStefano Zampini       } else {
18607b1e237SStefano Zampini         Vec         sol,rhs;
18707b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
18807b1e237SStefano Zampini         PetscInt    i,nrhs,n;
18907b1e237SStefano Zampini 
19007b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
19107b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
19207b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
19307b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
19407b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
19507b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
19607b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
19707b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
19807b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
19907b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
20007b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
20107b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
20207b1e237SStefano Zampini         }
20307b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
20407b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
20507b1e237SStefano Zampini       }
2063202ece2SStefano Zampini     }
2075ec10c6aSStefano Zampini     if (!Bdense & !issym) {
2083202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2093202ece2SStefano Zampini     }
2105ec10c6aSStefano Zampini 
2115ec10c6aSStefano Zampini     if (!issym) {
2123202ece2SStefano Zampini       if (!Cdense) {
2133202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
2143202ece2SStefano Zampini       } else {
2153202ece2SStefano Zampini         Cd = C;
2163202ece2SStefano Zampini       }
2175ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2183202ece2SStefano Zampini       if (!Cdense) {
2193202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
2203202ece2SStefano Zampini       }
2215ec10c6aSStefano Zampini     } else {
2225ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2235ec10c6aSStefano Zampini       if (!Bdense) {
2245ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2255ec10c6aSStefano Zampini       }
2265ec10c6aSStefano Zampini     }
2275ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
228f11841e3SStefano Zampini   }
2293202ece2SStefano Zampini 
2303202ece2SStefano Zampini   if (D) {
2313202ece2SStefano Zampini     Mat       Dd;
2323202ece2SStefano Zampini     PetscBool Ddense;
2333202ece2SStefano Zampini 
2343202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
2353202ece2SStefano Zampini     if (!Ddense) {
2363202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
2373202ece2SStefano Zampini     } else {
2383202ece2SStefano Zampini       Dd = D;
2393202ece2SStefano Zampini     }
240f11841e3SStefano Zampini     if (n_I) {
2413202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
242f11841e3SStefano Zampini     } else {
243f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
244f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
245f11841e3SStefano Zampini       } else {
246f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
247f11841e3SStefano Zampini       }
248f11841e3SStefano Zampini     }
2493202ece2SStefano Zampini     if (!Ddense) {
2503202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
2513202ece2SStefano Zampini     }
2523202ece2SStefano Zampini   } else {
2533202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
2543202ece2SStefano Zampini   }
2553202ece2SStefano Zampini   PetscFunctionReturn(0);
2563202ece2SStefano Zampini }
25734a97f8cSStefano Zampini 
25834a97f8cSStefano Zampini #undef __FUNCT__
2591580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
260683d3df6SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick)
261b1b3d7a2SStefano Zampini {
262be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
263be83ff47SStefano Zampini   Mat                    S_all;
2645a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
2655db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
266dc456d91SStefano Zampini   IS                     is_I,is_I_layer;
267dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
268dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
269dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
2705a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
271683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
27208122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
27306a4b1faSStefano Zampini   PetscScalar            *Bwork;
2745a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2755a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2765a95e1ceSStefano Zampini   MPI_Comm               comm_n;
27706a4e24aSStefano Zampini   PetscBool              use_getr = PETSC_FALSE;
278b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
279b1b3d7a2SStefano Zampini 
280b1b3d7a2SStefano Zampini   PetscFunctionBegin;
281a64f4aa4SStefano Zampini   /* update info in sub_schurs */
282a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
283a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
2843301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2853301b35fSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
286a64f4aa4SStefano Zampini   if (Ain) {
287a64f4aa4SStefano Zampini     PetscBool isseqaij;
2883301b35fSStefano Zampini     /* determine if we are dealing with hermitian positive definite problems */
2893301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
2903301b35fSStefano Zampini     if (Ain->symmetric_set) {
2913301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->symmetric;
2923301b35fSStefano Zampini     }
2933301b35fSStefano Zampini #else
2943301b35fSStefano Zampini     if (Ain->hermitian_set) {
2953301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->hermitian;
2963301b35fSStefano Zampini     }
2973301b35fSStefano Zampini #endif
2983301b35fSStefano Zampini     if (Ain->spd_set) {
2993301b35fSStefano Zampini       sub_schurs->is_posdef = Ain->spd;
3003301b35fSStefano Zampini     }
301a64f4aa4SStefano Zampini 
3023301b35fSStefano Zampini     /* check */
303a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
3041dac5d76SStefano Zampini     if (compute_Stilda && (!Ain->hermitian_set || !Ain->spd_set)) { /* these are lazy checks, maybe I should throw an error if they are not set */
3053301b35fSStefano Zampini       PetscInt lsize;
3063301b35fSStefano Zampini 
3073301b35fSStefano Zampini       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
3083301b35fSStefano Zampini       if (lsize) {
3093301b35fSStefano Zampini         PetscScalar val;
3103301b35fSStefano Zampini         PetscReal   norm;
3113301b35fSStefano Zampini         Vec         vec1,vec2,vec3;
3123301b35fSStefano Zampini 
3133301b35fSStefano Zampini         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
3143301b35fSStefano Zampini         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
3153301b35fSStefano Zampini         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
3163301b35fSStefano Zampini         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
3173301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
3183301b35fSStefano Zampini         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3193301b35fSStefano Zampini #else
3203301b35fSStefano Zampini         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3213301b35fSStefano Zampini #endif
3223301b35fSStefano Zampini         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
3233301b35fSStefano Zampini         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
3241dac5d76SStefano Zampini         if (!Ain->hermitian_set) {
325862806e4SStefano Zampini           if (norm > PetscSqrtReal(PETSC_SMALL)) {
3263301b35fSStefano Zampini             sub_schurs->is_hermitian = PETSC_FALSE;
3273301b35fSStefano Zampini           } else {
3283301b35fSStefano Zampini             sub_schurs->is_hermitian = PETSC_TRUE;
3293301b35fSStefano Zampini           }
3301dac5d76SStefano Zampini         }
3311dac5d76SStefano Zampini         if (!Ain->spd_set && !benign_trick) {
3323301b35fSStefano Zampini           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
3333301b35fSStefano Zampini           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
33406a4e24aSStefano Zampini         }
3353301b35fSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
3363301b35fSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
3373301b35fSStefano Zampini         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
3383301b35fSStefano Zampini       } else {
3393301b35fSStefano Zampini         sub_schurs->is_hermitian = PETSC_TRUE;
3403301b35fSStefano Zampini         sub_schurs->is_posdef = PETSC_TRUE;
3413301b35fSStefano Zampini       }
3423301b35fSStefano Zampini     }
343a64f4aa4SStefano Zampini     if (isseqaij) {
344a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
345a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
3463301b35fSStefano Zampini     } else {
347a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
348a64f4aa4SStefano Zampini     }
349a64f4aa4SStefano Zampini   }
3503301b35fSStefano Zampini 
351a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
352a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
353a64f4aa4SStefano Zampini   if (sub_schurs->use_mumps) {
354a64f4aa4SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
355a64f4aa4SStefano Zampini   }
356a64f4aa4SStefano Zampini 
3575a95e1ceSStefano Zampini   /* preliminary checks */
3585a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
3595a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
3605a95e1ceSStefano Zampini   }
3615a95e1ceSStefano Zampini 
3625a95e1ceSStefano Zampini   /* restrict work on active processes */
3635a95e1ceSStefano Zampini   color = 0;
3645a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
3655a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
3665a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
3675a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
3685a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
3695a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
3705a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
3715a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3725a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
3735a95e1ceSStefano Zampini     PetscFunctionReturn(0);
3745a95e1ceSStefano Zampini   }
37581ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
3765a95e1ceSStefano Zampini 
377b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
378883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
379a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
3803301b35fSStefano Zampini     PetscBool isseqsbaij;
381a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
3823301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
3833301b35fSStefano Zampini     if (isseqsbaij) {
384a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
385a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
386a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
387a64f4aa4SStefano Zampini     } else {
388a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
389a64f4aa4SStefano Zampini       A_BB = tA_BB;
390a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
391a64f4aa4SStefano Zampini       A_IB = tA_IB;
392a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
393a64f4aa4SStefano Zampini       A_BI = tA_BI;
394f11841e3SStefano Zampini     }
395a58a30b4SStefano Zampini   } else {
3965a95e1ceSStefano Zampini     A_II = NULL;
3975a95e1ceSStefano Zampini     A_IB = NULL;
3985a95e1ceSStefano Zampini     A_BI = NULL;
3995a95e1ceSStefano Zampini     A_BB = NULL;
400b1b3d7a2SStefano Zampini   }
4015a95e1ceSStefano Zampini   S_all = NULL;
402b1b3d7a2SStefano Zampini 
403b1b3d7a2SStefano Zampini   /* determine interior problems */
4043dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4053dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
406b1b3d7a2SStefano Zampini     PetscBT                touched;
407b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
408b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
409b1b3d7a2SStefano Zampini 
4103dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
4113dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
4123dc780c3SStefano Zampini     }
413b1b3d7a2SStefano Zampini     /* get sizes */
414b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
415b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
416b1b3d7a2SStefano Zampini 
417b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
418b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
419b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
420b1b3d7a2SStefano Zampini 
421b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
422b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
423b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
424b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
425b1b3d7a2SStefano Zampini     }
426b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
427b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
428b1b3d7a2SStefano Zampini 
429b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
430b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
431b1b3d7a2SStefano Zampini     n_prev_added = n_B;
432b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
433b1b3d7a2SStefano Zampini       PetscInt n_added;
434b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
435b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
436b1b3d7a2SStefano 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);
437b1b3d7a2SStefano Zampini       }
438b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
439b1b3d7a2SStefano Zampini       n_prev_added = n_added;
440b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
441b1b3d7a2SStefano Zampini       if (!n_added) break;
442b1b3d7a2SStefano Zampini     }
443b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
444b1b3d7a2SStefano Zampini 
445883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
446a9b99552SStefano 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);
447b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
448a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
449883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
450883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
451b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
452b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
453a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
454b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
455b1b3d7a2SStefano Zampini 
456b1b3d7a2SStefano Zampini       /* II block */
457b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
458b1b3d7a2SStefano Zampini     }
459b1b3d7a2SStefano Zampini   } else {
460b1b3d7a2SStefano Zampini     PetscInt n_I;
461b1b3d7a2SStefano Zampini 
462b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
463b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
464a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
465b1b3d7a2SStefano Zampini 
466b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
467883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
468b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
469b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
470b1b3d7a2SStefano Zampini 
471b1b3d7a2SStefano Zampini       /* II block is the same */
472b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
473b1b3d7a2SStefano Zampini       AE_II = A_II;
474b1b3d7a2SStefano Zampini     }
475b1b3d7a2SStefano Zampini   }
4765a95e1ceSStefano Zampini 
477883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
4785a95e1ceSStefano Zampini   max_subset_size = 0;
479883469d8SStefano Zampini   local_size = 0;
4805a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4815a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4825a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
483883469d8SStefano Zampini     local_size += subset_size;
484883469d8SStefano Zampini   }
485883469d8SStefano Zampini 
486883469d8SStefano Zampini   /* Work arrays for local indices */
487883469d8SStefano Zampini   extra = 0;
488683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
48981ea8064SStefano Zampini   if (sub_schurs->use_mumps && is_I_layer) {
490a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
491883469d8SStefano Zampini   }
492683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
493883469d8SStefano Zampini   if (extra) {
494883469d8SStefano Zampini     const PetscInt *idxs;
495a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
496883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
497a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
498883469d8SStefano Zampini   }
499883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
500dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
501dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
502883469d8SStefano Zampini 
503883469d8SStefano Zampini   /* Get local indices in local numbering */
504883469d8SStefano Zampini   local_size = 0;
5055a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
506883469d8SStefano Zampini     PetscInt j;
507883469d8SStefano Zampini     const    PetscInt *idxs;
508883469d8SStefano Zampini 
5095a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5105a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
511eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
512eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
513eb595f79SStefano Zampini     auxnum2[i] = subset_size;
514883469d8SStefano Zampini     /* subset indices in local numbering */
515883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5165a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
517883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
518883469d8SStefano Zampini     local_size += subset_size;
519883469d8SStefano Zampini   }
520883469d8SStefano Zampini 
5215a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
522d2627357SStefano Zampini   Bwork = NULL;
523d2627357SStefano Zampini   pivots = NULL;
52406a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
525d2627357SStefano Zampini     PetscScalar lwork;
526d2627357SStefano Zampini 
52706a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
528d2627357SStefano Zampini     B_lwork = -1;
529d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
530d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
531d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
532d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
533d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
534d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
535d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
536d2627357SStefano Zampini   }
537d2627357SStefano Zampini 
538d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
539dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
540dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
541dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
542dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
543dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
544dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
545dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
546dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
547dc456d91SStefano Zampini   if (i != local_size) {
548dc456d91SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size);
549eb595f79SStefano Zampini   }
550dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
551e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
552d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
553d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
554d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
555d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
5562972d61bSStefano Zampini 
5575a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
5585a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
5595a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
5605a95e1ceSStefano Zampini 
5615a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
5625a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
5635a95e1ceSStefano Zampini     if (subset_size != local_size) {
5645a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
5655a95e1ceSStefano Zampini     }
5665a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
567b1b3d7a2SStefano Zampini   }
568b1b3d7a2SStefano Zampini 
5695a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
5705a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
5715a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
5725a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
5735a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
5745a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
575aa83b6aeSStefano Zampini   }
576b1b3d7a2SStefano Zampini 
5775a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
578be83ff47SStefano Zampini   F = NULL;
579683d3df6SStefano Zampini   if (!sub_schurs->use_mumps) { /* this code branch is used when MUMPS is not present; it is not very efficient, unless the economic version of the scaling is required */
5805a95e1ceSStefano Zampini     Mat         S_Ej_expl;
5815a95e1ceSStefano Zampini     PetscScalar *work;
5825a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
5835a95e1ceSStefano Zampini     PetscBool   Sdense;
5845a95e1ceSStefano Zampini 
5855a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
5865a95e1ceSStefano Zampini     local_size = 0;
587b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
5885a95e1ceSStefano Zampini       IS  is_subset_B;
5895a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
5905a95e1ceSStefano Zampini 
5915a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
5925a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
5935a95e1ceSStefano Zampini       /* EE block */
5945a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
5955a95e1ceSStefano Zampini       /* IE block */
5965a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
5975a95e1ceSStefano Zampini       /* EI block */
5985a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
5995a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
6005a95e1ceSStefano Zampini       } else {
6015a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
6025a95e1ceSStefano Zampini       }
603a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
6045a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
6055a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
6065a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
6075a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
608b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
609b1b3d7a2SStefano Zampini         KSP ksp;
610b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
6115a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
612b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
613b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
614b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
615b1b3d7a2SStefano Zampini         KSPType   ksp_type;
616b1b3d7a2SStefano Zampini         PetscInt  n_internal;
6175a95e1ceSStefano Zampini         PetscBool ispcnone;
618b1b3d7a2SStefano Zampini 
619b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
6205a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
621b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
622b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
623b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
624b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
6255a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
6265a95e1ceSStefano Zampini         if (!ispcnone) {
6275a95e1ceSStefano Zampini           PCType pc_type;
628b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
629b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
6305a95e1ceSStefano Zampini         } else {
6315a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
6325a95e1ceSStefano Zampini         }
633b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
634b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
635b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
636b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
637b1b3d7a2SStefano Zampini           if (solver) {
638b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
639b1b3d7a2SStefano Zampini           }
640b1b3d7a2SStefano Zampini         }
641b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
642b1b3d7a2SStefano Zampini       }
6435a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6445a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
6455a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
6465a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
6475a95e1ceSStefano Zampini       if (Sdense) {
6485a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
6495a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
650b1b3d7a2SStefano Zampini         }
6515a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6525a95e1ceSStefano Zampini       } else {
6535a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
6545a95e1ceSStefano Zampini       }
6555a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
656a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
6575a95e1ceSStefano Zampini       local_size += subset_size;
6585a95e1ceSStefano Zampini     }
6595a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
660b1b3d7a2SStefano Zampini     /* free */
661b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
662b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
6635a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
664883469d8SStefano Zampini   } else {
665be83ff47SStefano Zampini     Mat         A;
666883469d8SStefano Zampini     IS          is_A_all;
667683d3df6SStefano Zampini     PetscScalar *work,*S_data,*schur_factor;
668*6dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
669683d3df6SStefano Zampini     PetscBool   economic,mumps_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
670883469d8SStefano Zampini 
671683d3df6SStefano Zampini     /* get sizes */
67281ea8064SStefano Zampini     n_I = 0;
67381ea8064SStefano Zampini     if (is_I_layer) {
674a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
67581ea8064SStefano Zampini     }
676683d3df6SStefano Zampini     economic = PETSC_FALSE;
677683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
678683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
679683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
680d62866d3SStefano Zampini 
681683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
682683d3df6SStefano Zampini     size_active_schur = local_size;
683683d3df6SStefano Zampini     cum = n_I+size_active_schur;
684683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
685683d3df6SStefano Zampini       const PetscInt* idxs;
686683d3df6SStefano Zampini       PetscInt        n_dir;
687683d3df6SStefano Zampini 
688683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
689683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
690683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
691683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
692683d3df6SStefano Zampini       cum += n_dir;
693d62866d3SStefano Zampini     }
694683d3df6SStefano Zampini     factor_workaround = PETSC_FALSE;
695683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
696683d3df6SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && compute_Stilda) {
697683d3df6SStefano Zampini       PetscInt n_v;
698683d3df6SStefano Zampini 
699683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
700683d3df6SStefano Zampini       if (n_v) {
701683d3df6SStefano Zampini         const PetscInt* idxs;
702683d3df6SStefano Zampini 
703683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
704683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
705683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
706683d3df6SStefano Zampini         cum += n_v;
707683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
708683d3df6SStefano Zampini       }
709683d3df6SStefano Zampini     }
710683d3df6SStefano Zampini     size_schur = cum - n_I;
711683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
712683d3df6SStefano Zampini     /* get working mat (TODO: factorize without actually permuting)  */
713683d3df6SStefano Zampini     if (cum == n) {
714683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
715683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
716683d3df6SStefano Zampini     } else {
7176816873aSStefano Zampini       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
718683d3df6SStefano Zampini     }
719a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
7203301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
7213301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
7223301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
723883469d8SStefano Zampini 
724683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
725d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
726d47842beSStefano Zampini 
727683d3df6SStefano Zampini     if (n_I) { /* TODO add ordering from options */
7285a05ddb0SStefano Zampini       IS is_schur;
7295a05ddb0SStefano Zampini 
7309ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
731883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
732883469d8SStefano Zampini       } else {
733883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
734883469d8SStefano Zampini       }
735883469d8SStefano Zampini       /* subsets ordered last */
7365a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
7375a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
7385a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
739883469d8SStefano Zampini 
740883469d8SStefano Zampini       /* factorization step */
7419ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
742883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
743be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
744be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
745a0b0af32SStefano Zampini         if (sub_schurs->is_posdef) {
746a0b0af32SStefano Zampini           ierr = MatMumpsSchurComplementSetSym(F,1);CHKERRQ(ierr);
747a0b0af32SStefano Zampini         } else {
748a0b0af32SStefano Zampini           ierr = MatMumpsSchurComplementSetSym(F,2);CHKERRQ(ierr);
749a0b0af32SStefano Zampini         }
750be83ff47SStefano Zampini #endif
751883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
752a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
753883469d8SStefano Zampini       } else {
754883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
755be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
756be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
757be83ff47SStefano Zampini #endif
758883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
759883469d8SStefano Zampini       }
760883469d8SStefano Zampini 
761883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
7625a05ddb0SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
763d5574798SStefano Zampini 
764d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
765683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
766683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
7674a6c6b0dSStefano Zampini       mumps_S = PETSC_TRUE;
768be83ff47SStefano Zampini     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
769be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
770683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
7714a6c6b0dSStefano Zampini       mumps_S = PETSC_FALSE;
772be83ff47SStefano Zampini     }
773be83ff47SStefano Zampini 
774be83ff47SStefano Zampini     if (reuse_solvers) {
775d62866d3SStefano Zampini       Mat              A_II;
77653892102SStefano Zampini       Vec              vec1_B;
777d62866d3SStefano Zampini       PCBDDCReuseMumps msolv_ctx;
778d5574798SStefano Zampini 
779d62866d3SStefano Zampini       if (sub_schurs->reuse_mumps) {
7806816873aSStefano Zampini         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
781e28d306cSStefano Zampini       } else {
782e28d306cSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
783d62866d3SStefano Zampini       }
784e28d306cSStefano Zampini       msolv_ctx = sub_schurs->reuse_mumps;
785be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
786d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
787d5574798SStefano Zampini       msolv_ctx->F = F;
788683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
789683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
790683d3df6SStefano Zampini       {
791683d3df6SStefano Zampini         PetscScalar *array;
792683d3df6SStefano Zampini         PetscInt    n;
793683d3df6SStefano Zampini 
794683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
795683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
796683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
797683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
798683d3df6SStefano Zampini       }
799683d3df6SStefano Zampini       msolv_ctx->has_vertices = factor_workaround;
800d62866d3SStefano Zampini 
801d62866d3SStefano Zampini       /* interior solver */
802683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
803d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
804d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
805d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
806d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
807e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
808d62866d3SStefano Zampini 
809d62866d3SStefano Zampini       /* correction solver */
810683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A),&msolv_ctx->correction_solver);CHKERRQ(ierr);
811d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
812d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
813d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
814d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
815e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
81653892102SStefano Zampini 
81753892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
81853892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
81953892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
820683d3df6SStefano Zampini       if (!factor_workaround) {
82153892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
82253892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
82353892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
82453892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
825683d3df6SStefano Zampini       } else {
826683d3df6SStefano Zampini         IS              is_B_all;
827683d3df6SStefano Zampini         const PetscInt* idxs;
828683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
829683d3df6SStefano Zampini 
830683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
831683d3df6SStefano Zampini         dual = size_schur - n_v;
832683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
833683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
834683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
835683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
836683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
837683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
838683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
839683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
840683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
841683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
842683d3df6SStefano Zampini       }
843683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
844d5574798SStefano Zampini     }
84508122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
84653892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
8475db18549SStefano Zampini 
848be83ff47SStefano Zampini     /* Work arrays */
849be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
8505a94c880SStefano Zampini 
8515a94c880SStefano Zampini     /* matrices for adaptive selection */
85212d906b1SStefano Zampini     if (compute_Stilda) {
8535a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
8545a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
8555a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
8565a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
8575a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
8585a94c880SStefano Zampini       }
8595a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all) {
8605a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
8615a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
8625a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
8635a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
8645a94c880SStefano Zampini       }
86512d906b1SStefano Zampini     }
866d2627357SStefano Zampini 
867be83ff47SStefano Zampini     /* S_Ej_all */
868be83ff47SStefano Zampini     cum = cum2 = 0;
869be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
87065d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
87165d8bf0aSStefano Zampini       PetscInt j;
87265d8bf0aSStefano Zampini 
8735a95e1ceSStefano Zampini       /* get S_E */
874b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
875683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
876be83ff47SStefano Zampini         PetscInt k;
877be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
878be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
879be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
880683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
881be83ff47SStefano Zampini           }
882be83ff47SStefano Zampini         }
88306a4e24aSStefano Zampini       } else { /* just copy to workspace */
884be83ff47SStefano Zampini         PetscInt k;
885be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
886be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
887be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
888be83ff47SStefano Zampini           }
889be83ff47SStefano Zampini         }
8909087bf02SStefano Zampini       }
8915a95e1ceSStefano Zampini       /* insert S_E values */
892be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
8935a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
894a1337663SStefano Zampini 
895683d3df6SStefano Zampini       /* if adaptivity is requested, invert S_E blocks */
896862806e4SStefano Zampini       if (compute_Stilda) {
89708122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
89808122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
89906a4e24aSStefano Zampini         if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
9005a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
9012972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
9025a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
9032972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
904d6462365SStefano Zampini         } else {
905d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
906d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
907d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
908d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
9092972d61bSStefano Zampini         }
91008122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
9115a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
9129087bf02SStefano Zampini       }
913be83ff47SStefano Zampini       cum += subset_size;
914be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
915be83ff47SStefano Zampini     }
916be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
9174a6c6b0dSStefano Zampini 
9184a6c6b0dSStefano Zampini     if (mumps_S) {
9195a05ddb0SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
9204a6c6b0dSStefano Zampini     }
921683d3df6SStefano Zampini 
922683d3df6SStefano Zampini     schur_factor = NULL;
92345951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
924683d3df6SStefano Zampini 
9254a6c6b0dSStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
9264a6c6b0dSStefano Zampini         PetscInt j;
9274a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
9285a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
9294a6c6b0dSStefano Zampini       } else {
930683d3df6SStefano Zampini         Mat S_all_inv=NULL;
9315a05ddb0SStefano Zampini         if (mumps_S) { /* use MatFactor calls to invert S */
932*6dba178dSStefano Zampini             /* let MatFactor factorize the Schur complement */
933*6dba178dSStefano Zampini           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
934683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
935683d3df6SStefano Zampini              The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */
936683d3df6SStefano Zampini           if (factor_workaround) {
937683d3df6SStefano Zampini             PetscScalar *data;
938683d3df6SStefano Zampini             PetscInt nd = 0;
939*6dba178dSStefano Zampini 
940683d3df6SStefano Zampini             if (!sub_schurs->is_posdef) {
941683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
942683d3df6SStefano Zampini             }
943*6dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
944683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
945683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
946683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
947683d3df6SStefano Zampini             }
948683d3df6SStefano Zampini             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
949683d3df6SStefano Zampini             cum = 0;
950683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
951683d3df6SStefano Zampini               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
952683d3df6SStefano Zampini               cum += size_active_schur-i;
953683d3df6SStefano Zampini             }
954683d3df6SStefano Zampini             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
955*6dba178dSStefano Zampini             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
956683d3df6SStefano Zampini             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
957683d3df6SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
958683d3df6SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
959683d3df6SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
960683d3df6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
961683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
962683d3df6SStefano Zampini           } else {
963*6dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
964*6dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
965683d3df6SStefano Zampini           }
9664a6c6b0dSStefano Zampini         } else { /* we need to invert explicitly since we are not using MUMPS for S */
967683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
968683d3df6SStefano Zampini           S_all_inv = S_all;
969683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
970be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
971be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
97206a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
973be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
974be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
975be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
976be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
977d6462365SStefano Zampini           } else {
978d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
979d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
980d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
981d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
982be83ff47SStefano Zampini           }
983be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
984683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
985be83ff47SStefano Zampini         }
986be83ff47SStefano Zampini         /* S_Ej_tilda_all */
987be83ff47SStefano Zampini         cum = cum2 = 0;
988683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
989be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
990be83ff47SStefano Zampini           PetscInt j;
991862806e4SStefano Zampini 
992862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
993be83ff47SStefano Zampini           /* get (St^-1)_E */
99406a4e24aSStefano Zampini           /* Here I don't need to expand to upper triangular since St^-1
99506a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
996a0b0af32SStefano Zampini           if (S_lower_triangular) {
997be83ff47SStefano Zampini             PetscInt k;
998be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
999be83ff47SStefano Zampini               for (j=k;j<subset_size;j++) {
1000be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1001be83ff47SStefano Zampini               }
1002be83ff47SStefano Zampini             }
1003be83ff47SStefano Zampini           } else { /* copy to workspace */
1004be83ff47SStefano Zampini             PetscInt k;
1005be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1006be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1007be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1008be83ff47SStefano Zampini               }
1009be83ff47SStefano Zampini             }
1010be83ff47SStefano Zampini           }
1011be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
10125a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1013be83ff47SStefano Zampini           cum += subset_size;
1014be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1015883469d8SStefano Zampini         }
1016683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
10174a6c6b0dSStefano Zampini         if (mumps_S) {
1018*6dba178dSStefano Zampini           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
10195db18549SStefano Zampini         }
1020683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1021683d3df6SStefano Zampini       }
1022683d3df6SStefano Zampini 
1023683d3df6SStefano Zampini       /* move back factors to Schur data into F */
1024683d3df6SStefano Zampini       if (factor_workaround) {
1025683d3df6SStefano Zampini         Mat S_tmp;
1026683d3df6SStefano Zampini         PetscInt nv,nd=0;
1027683d3df6SStefano Zampini 
1028683d3df6SStefano Zampini         if (!mumps_S) {
1029683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1030683d3df6SStefano Zampini         }
1031683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1032*6dba178dSStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1033683d3df6SStefano Zampini         if (sub_schurs->is_posdef) {
1034683d3df6SStefano Zampini           PetscScalar *data;
1035683d3df6SStefano Zampini 
1036683d3df6SStefano Zampini           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1037683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1038683d3df6SStefano Zampini 
1039683d3df6SStefano Zampini           if (S_lower_triangular) {
1040683d3df6SStefano Zampini             cum = 0;
1041683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1042683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1043683d3df6SStefano Zampini               cum += size_active_schur-i;
1044683d3df6SStefano Zampini 	    }
1045683d3df6SStefano Zampini           } else {
1046683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1047683d3df6SStefano Zampini           }
1048683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1049683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1050683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1051683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1052683d3df6SStefano Zampini 	    }
1053683d3df6SStefano Zampini           }
1054*6dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1055683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1056683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
1057683d3df6SStefano Zampini             data[i*(size_schur+1)] = PETSC_MAX_REAL;
1058683d3df6SStefano Zampini           }
1059683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1060683d3df6SStefano Zampini         } else {
1061683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1062683d3df6SStefano Zampini         }
1063*6dba178dSStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
10649087bf02SStefano Zampini       }
10654a6c6b0dSStefano Zampini     }
10664a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1067683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
10684a6c6b0dSStefano Zampini   }
10695a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1070be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1071a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1072a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1073a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1074a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1075a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
10765db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10775db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10785a95e1ceSStefano Zampini   if (compute_Stilda) {
10795a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10805a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10815a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10825a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108308122e43SStefano Zampini   }
1084a1337663SStefano Zampini 
10855db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
10865db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
10873927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
10883927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
10895a95e1ceSStefano Zampini 
10905db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
10915a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
1092dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
10935a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
10945a94c880SStefano Zampini   } else {
10955a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
10965a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
1097dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
10985a94c880SStefano Zampini   }
109908122e43SStefano Zampini 
1100ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
1101ac632422SStefano Zampini   if (faster_deluxe) {
11025a95e1ceSStefano Zampini     Mat         tmpmat;
11035a95e1ceSStefano Zampini     PetscScalar *array;
11045a95e1ceSStefano Zampini     PetscInt    cum;
11055a95e1ceSStefano Zampini 
11065a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
11075a95e1ceSStefano Zampini     cum = 0;
11085a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
11095a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
11105a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
11115a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
111206a4e24aSStefano Zampini       if (!use_getr) {
11135a95e1ceSStefano Zampini         PetscInt j,k;
11145a95e1ceSStefano Zampini 
11155a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
11165a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
11175a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
11185a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
11195a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
11205a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
11215a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
11225a95e1ceSStefano Zampini           }
11235a95e1ceSStefano Zampini         }
1124d6462365SStefano Zampini       } else {
1125d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1126d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1127d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1128d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
11295a95e1ceSStefano Zampini       }
11305a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
11315a95e1ceSStefano Zampini       cum += subset_size*subset_size;
11325a95e1ceSStefano Zampini     }
11335a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
11345a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
11355a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1136ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
11375a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
11385a95e1ceSStefano Zampini   }
11395a95e1ceSStefano Zampini 
1140f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
11415a95e1ceSStefano Zampini   if (compute_Stilda) {
11425a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1143a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
11445a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1145dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
11465a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
114708122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
11485a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_inv_all;
1149dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
115008122e43SStefano Zampini   }
11513202ece2SStefano Zampini 
11525a95e1ceSStefano Zampini   /* free workspace */
11535a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
115406a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1155a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
11563202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1157dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
11585a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1159b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1160b1b3d7a2SStefano Zampini }
1161b1b3d7a2SStefano Zampini 
1162b1b3d7a2SStefano Zampini #undef __FUNCT__
1163b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
1164a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1165b1b3d7a2SStefano Zampini {
11669bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
11675a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1168b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1169b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1170b1b3d7a2SStefano Zampini 
1171b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1172b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1173b1b3d7a2SStefano Zampini   if (!is_sorted) {
1174b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1175b1b3d7a2SStefano Zampini   }
1176b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1177b1b3d7a2SStefano Zampini   if (!is_sorted) {
1178b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1179b1b3d7a2SStefano Zampini   }
1180b1b3d7a2SStefano Zampini 
1181b1b3d7a2SStefano Zampini   /* reset any previous data */
1182b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1183b1b3d7a2SStefano Zampini 
11845a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
11859bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1186b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
118708122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1188b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1189b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
1190b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
1191b1b3d7a2SStefano Zampini   }
1192b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
1193b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
119408122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1195b1b3d7a2SStefano Zampini   }
1196b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
1197b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
1198d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1199d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1200b1b3d7a2SStefano Zampini 
12016816873aSStefano Zampini   /* Determine if MUMPS can be used */
1202883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
1203883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1204a64f4aa4SStefano Zampini   sub_schurs->use_mumps = PETSC_TRUE;
1205883469d8SStefano Zampini #endif
1206b1b3d7a2SStefano Zampini 
1207b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1208b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1209b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1210b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
12115db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
12125db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
12135db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
12145db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
12155a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1216b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1217a64f4aa4SStefano Zampini   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1218b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1219b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1220b96c3477SStefano Zampini     }
12219bb4a8caSStefano Zampini   }
1222d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
1223b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1224b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
122508122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1226b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1227b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1228b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1229b1b3d7a2SStefano Zampini }
1230b1b3d7a2SStefano Zampini 
1231b1b3d7a2SStefano Zampini #undef __FUNCT__
123234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
123334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
123434a97f8cSStefano Zampini {
123534a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
123634a97f8cSStefano Zampini   PetscErrorCode  ierr;
123734a97f8cSStefano Zampini 
123834a97f8cSStefano Zampini   PetscFunctionBegin;
123934a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
12405ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
124134a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
124234a97f8cSStefano Zampini   PetscFunctionReturn(0);
124334a97f8cSStefano Zampini }
124434a97f8cSStefano Zampini 
124534a97f8cSStefano Zampini #undef __FUNCT__
124634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
124734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
124834a97f8cSStefano Zampini {
124934a97f8cSStefano Zampini   PetscErrorCode ierr;
125034a97f8cSStefano Zampini 
125134a97f8cSStefano Zampini   PetscFunctionBegin;
125234a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
125334a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
125434a97f8cSStefano Zampini   PetscFunctionReturn(0);
125534a97f8cSStefano Zampini }
125634a97f8cSStefano Zampini 
125734a97f8cSStefano Zampini #undef __FUNCT__
125834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
125934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
126034a97f8cSStefano Zampini {
126134a97f8cSStefano Zampini   PetscInt       i;
126234a97f8cSStefano Zampini   PetscErrorCode ierr;
126334a97f8cSStefano Zampini 
126434a97f8cSStefano Zampini   PetscFunctionBegin;
12651e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1266b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1267b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1268b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
12695db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
12705db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
127141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
127241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
127308122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1274a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
12755db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1276d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1277d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
127808122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
127934a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1280b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
128134a97f8cSStefano Zampini   }
12825ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1283b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
12843dc780c3SStefano Zampini   }
1285d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps) {
1286d62866d3SStefano Zampini     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1287d62866d3SStefano Zampini   }
1288d62866d3SStefano Zampini   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
128934a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
129034a97f8cSStefano Zampini   PetscFunctionReturn(0);
129134a97f8cSStefano Zampini }
129234a97f8cSStefano Zampini 
129334a97f8cSStefano Zampini #undef __FUNCT__
129434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
12952a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
129634a97f8cSStefano Zampini {
129734a97f8cSStefano Zampini   PetscInt       i,j,n;
129834a97f8cSStefano Zampini   PetscErrorCode ierr;
129934a97f8cSStefano Zampini 
130034a97f8cSStefano Zampini   PetscFunctionBegin;
130134a97f8cSStefano Zampini   n = 0;
130234a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
130334a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
130434a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
130534a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
130634a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
130734a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
130834a97f8cSStefano Zampini         queue_tip[n] = dof;
130934a97f8cSStefano Zampini         n++;
131034a97f8cSStefano Zampini       }
131134a97f8cSStefano Zampini     }
131234a97f8cSStefano Zampini   }
131334a97f8cSStefano Zampini   *n_added = n;
131434a97f8cSStefano Zampini   PetscFunctionReturn(0);
131534a97f8cSStefano Zampini }
1316