xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision e8ade678444ba846a4c00911e020f3c490b8868c)
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;
27d4933d67SStefano Zampini   } else { /* interior solver */
286dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
29683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
306dba178dSStefano Zampini #endif
31d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
32d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
33d4933d67SStefano Zampini #endif
34683d3df6SStefano Zampini     copy = PETSC_TRUE;
35683d3df6SStefano Zampini   }
36683d3df6SStefano Zampini   ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
37683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
38683d3df6SStefano Zampini   if (copy) {
39683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
40683d3df6SStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
41683d3df6SStefano Zampini     ierr = PetscMemcpy(array_mumps,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
42683d3df6SStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
43683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
44683d3df6SStefano Zampini 
45683d3df6SStefano Zampini     if (transpose) {
46683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
47683d3df6SStefano Zampini     } else {
48683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
49683d3df6SStefano Zampini     }
50683d3df6SStefano Zampini 
51683d3df6SStefano Zampini     /* get back data to caller worskpace */
52683d3df6SStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
53683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
54683d3df6SStefano Zampini     ierr = PetscMemcpy(array,array_mumps,n*sizeof(PetscScalar));CHKERRQ(ierr);
55683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
56683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
57683d3df6SStefano Zampini   } else {
58e28d306cSStefano Zampini     if (transpose) {
59e28d306cSStefano Zampini       ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
60e28d306cSStefano Zampini     } else {
616816873aSStefano Zampini       ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
62e28d306cSStefano Zampini     }
63683d3df6SStefano Zampini   }
64d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
65d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
66d4933d67SStefano Zampini #endif
67d62866d3SStefano Zampini   PetscFunctionReturn(0);
68d62866d3SStefano Zampini }
69d62866d3SStefano Zampini 
70d62866d3SStefano Zampini #undef __FUNCT__
71e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
72e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
73e28d306cSStefano Zampini {
74e28d306cSStefano Zampini   PetscErrorCode   ierr;
75e28d306cSStefano Zampini 
76e28d306cSStefano Zampini   PetscFunctionBegin;
77683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
78e28d306cSStefano Zampini   PetscFunctionReturn(0);
79e28d306cSStefano Zampini }
80e28d306cSStefano Zampini 
81e28d306cSStefano Zampini #undef __FUNCT__
82e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose"
83e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol)
84e28d306cSStefano Zampini {
85e28d306cSStefano Zampini   PetscErrorCode   ierr;
86e28d306cSStefano Zampini 
87e28d306cSStefano Zampini   PetscFunctionBegin;
88683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
89683d3df6SStefano Zampini   PetscFunctionReturn(0);
90683d3df6SStefano Zampini }
91683d3df6SStefano Zampini 
92683d3df6SStefano Zampini #undef __FUNCT__
93683d3df6SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
94683d3df6SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
95683d3df6SStefano Zampini {
96683d3df6SStefano Zampini   PetscErrorCode   ierr;
97683d3df6SStefano Zampini 
98683d3df6SStefano Zampini   PetscFunctionBegin;
99683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
100683d3df6SStefano Zampini   PetscFunctionReturn(0);
101683d3df6SStefano Zampini }
102683d3df6SStefano Zampini 
103683d3df6SStefano Zampini #undef __FUNCT__
104683d3df6SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose"
105683d3df6SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol)
106683d3df6SStefano Zampini {
107683d3df6SStefano Zampini   PetscErrorCode   ierr;
108683d3df6SStefano Zampini 
109683d3df6SStefano Zampini   PetscFunctionBegin;
110683d3df6SStefano Zampini   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
111e28d306cSStefano Zampini   PetscFunctionReturn(0);
112e28d306cSStefano Zampini }
113e28d306cSStefano Zampini 
114e28d306cSStefano Zampini #undef __FUNCT__
115d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCReuseMumpsReset"
116d62866d3SStefano Zampini static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
117d62866d3SStefano Zampini {
118d62866d3SStefano Zampini   PetscErrorCode ierr;
119d62866d3SStefano Zampini 
120d62866d3SStefano Zampini   PetscFunctionBegin;
121d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
122d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
123d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
124d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
125d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
12653892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
12753892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
128d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
12953892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
13053892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
131d62866d3SStefano Zampini   PetscFunctionReturn(0);
132d62866d3SStefano Zampini }
133d5574798SStefano Zampini 
134d5574798SStefano Zampini #undef __FUNCT__
1353202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
1365ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
1373202ece2SStefano Zampini {
1383202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
1393202ece2SStefano Zampini   KSP            ksp;
1403202ece2SStefano Zampini   PC             pc;
1413202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
1423202ece2SStefano Zampini   PetscReal      fill = 2.0;
143f11841e3SStefano Zampini   PetscInt       n_I;
1443202ece2SStefano Zampini   PetscMPIInt    size;
1453202ece2SStefano Zampini   PetscErrorCode ierr;
1463202ece2SStefano Zampini 
1473202ece2SStefano Zampini   PetscFunctionBegin;
1483202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
1493202ece2SStefano Zampini   if (size != 1) {
1503202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
1513202ece2SStefano Zampini   }
152f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
153f11841e3SStefano Zampini     PetscBool Sdense;
154f11841e3SStefano Zampini 
155f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
156f11841e3SStefano Zampini     if (!Sdense) {
157683d3df6SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should be dense");
158f11841e3SStefano Zampini     }
159f11841e3SStefano Zampini   }
1603202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
1613202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
1623202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1633202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
1643202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
1653202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
1663202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
1673202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
168f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
169f11841e3SStefano Zampini   if (n_I) {
1703202ece2SStefano Zampini     if (!Bdense) {
1713202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
1723202ece2SStefano Zampini     } else {
1733202ece2SStefano Zampini       Bd = B;
1743202ece2SStefano Zampini     }
1753202ece2SStefano Zampini 
1763202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1773202ece2SStefano Zampini       Mat fact;
1783202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1793202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1803202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1813202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1823202ece2SStefano Zampini     } else {
18307b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
18407b1e237SStefano Zampini 
18507b1e237SStefano Zampini       if (ex) {
1863202ece2SStefano Zampini         Mat Ainvd;
1873202ece2SStefano Zampini 
1883202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1893202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1903202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
19107b1e237SStefano Zampini       } else {
19207b1e237SStefano Zampini         Vec         sol,rhs;
19307b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
19407b1e237SStefano Zampini         PetscInt    i,nrhs,n;
19507b1e237SStefano Zampini 
19607b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
19707b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
19807b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
19907b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
20007b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
20107b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
20207b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
20307b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
20407b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
20507b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
20607b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
20707b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
20807b1e237SStefano Zampini         }
20907b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
21007b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
21107b1e237SStefano Zampini       }
2123202ece2SStefano Zampini     }
2135ec10c6aSStefano Zampini     if (!Bdense & !issym) {
2143202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2153202ece2SStefano Zampini     }
2165ec10c6aSStefano Zampini 
2175ec10c6aSStefano Zampini     if (!issym) {
2183202ece2SStefano Zampini       if (!Cdense) {
2193202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
2203202ece2SStefano Zampini       } else {
2213202ece2SStefano Zampini         Cd = C;
2223202ece2SStefano Zampini       }
2235ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2243202ece2SStefano Zampini       if (!Cdense) {
2253202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
2263202ece2SStefano Zampini       }
2275ec10c6aSStefano Zampini     } else {
2285ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2295ec10c6aSStefano Zampini       if (!Bdense) {
2305ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2315ec10c6aSStefano Zampini       }
2325ec10c6aSStefano Zampini     }
2335ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
234f11841e3SStefano Zampini   }
2353202ece2SStefano Zampini 
2363202ece2SStefano Zampini   if (D) {
2373202ece2SStefano Zampini     Mat       Dd;
2383202ece2SStefano Zampini     PetscBool Ddense;
2393202ece2SStefano Zampini 
2403202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
2413202ece2SStefano Zampini     if (!Ddense) {
2423202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
2433202ece2SStefano Zampini     } else {
2443202ece2SStefano Zampini       Dd = D;
2453202ece2SStefano Zampini     }
246f11841e3SStefano Zampini     if (n_I) {
2473202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
248f11841e3SStefano Zampini     } else {
249f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
250f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
251f11841e3SStefano Zampini       } else {
252f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
253f11841e3SStefano Zampini       }
254f11841e3SStefano Zampini     }
2553202ece2SStefano Zampini     if (!Ddense) {
2563202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
2573202ece2SStefano Zampini     }
2583202ece2SStefano Zampini   } else {
2593202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
2603202ece2SStefano Zampini   }
2613202ece2SStefano Zampini   PetscFunctionReturn(0);
2623202ece2SStefano Zampini }
26334a97f8cSStefano Zampini 
26434a97f8cSStefano Zampini #undef __FUNCT__
2651580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
266683d3df6SStefano 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)
267b1b3d7a2SStefano Zampini {
268be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
269be83ff47SStefano Zampini   Mat                    S_all;
2705a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
2715db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
272dc456d91SStefano Zampini   IS                     is_I,is_I_layer;
273dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
274dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
275dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
2765a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
277683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
27808122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
27906a4b1faSStefano Zampini   PetscScalar            *Bwork;
2805a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2815a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2825a95e1ceSStefano Zampini   MPI_Comm               comm_n;
28306a4e24aSStefano Zampini   PetscBool              use_getr = PETSC_FALSE;
284b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
285b1b3d7a2SStefano Zampini 
286b1b3d7a2SStefano Zampini   PetscFunctionBegin;
287a64f4aa4SStefano Zampini   /* update info in sub_schurs */
288a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
289a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
2903301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2913301b35fSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
292a64f4aa4SStefano Zampini   if (Ain) {
293a64f4aa4SStefano Zampini     PetscBool isseqaij;
2943301b35fSStefano Zampini     /* determine if we are dealing with hermitian positive definite problems */
2953301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
2963301b35fSStefano Zampini     if (Ain->symmetric_set) {
2973301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->symmetric;
2983301b35fSStefano Zampini     }
2993301b35fSStefano Zampini #else
3003301b35fSStefano Zampini     if (Ain->hermitian_set) {
3013301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->hermitian;
3023301b35fSStefano Zampini     }
3033301b35fSStefano Zampini #endif
3043301b35fSStefano Zampini     if (Ain->spd_set) {
3053301b35fSStefano Zampini       sub_schurs->is_posdef = Ain->spd;
3063301b35fSStefano Zampini     }
307a64f4aa4SStefano Zampini 
3083301b35fSStefano Zampini     /* check */
309a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
3101dac5d76SStefano 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 */
3113301b35fSStefano Zampini       PetscInt lsize;
3123301b35fSStefano Zampini 
3133301b35fSStefano Zampini       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
3143301b35fSStefano Zampini       if (lsize) {
3153301b35fSStefano Zampini         PetscScalar val;
3163301b35fSStefano Zampini         PetscReal   norm;
3173301b35fSStefano Zampini         Vec         vec1,vec2,vec3;
3183301b35fSStefano Zampini 
3193301b35fSStefano Zampini         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
3203301b35fSStefano Zampini         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
3213301b35fSStefano Zampini         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
3223301b35fSStefano Zampini         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
3233301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
3243301b35fSStefano Zampini         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3253301b35fSStefano Zampini #else
3263301b35fSStefano Zampini         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3273301b35fSStefano Zampini #endif
3283301b35fSStefano Zampini         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
3293301b35fSStefano Zampini         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
3301dac5d76SStefano Zampini         if (!Ain->hermitian_set) {
331862806e4SStefano Zampini           if (norm > PetscSqrtReal(PETSC_SMALL)) {
3323301b35fSStefano Zampini             sub_schurs->is_hermitian = PETSC_FALSE;
3333301b35fSStefano Zampini           } else {
3343301b35fSStefano Zampini             sub_schurs->is_hermitian = PETSC_TRUE;
3353301b35fSStefano Zampini           }
3361dac5d76SStefano Zampini         }
3371dac5d76SStefano Zampini         if (!Ain->spd_set && !benign_trick) {
3383301b35fSStefano Zampini           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
3393301b35fSStefano Zampini           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
34006a4e24aSStefano Zampini         }
3413301b35fSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
3423301b35fSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
3433301b35fSStefano Zampini         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
3443301b35fSStefano Zampini       } else {
3453301b35fSStefano Zampini         sub_schurs->is_hermitian = PETSC_TRUE;
3463301b35fSStefano Zampini         sub_schurs->is_posdef = PETSC_TRUE;
3473301b35fSStefano Zampini       }
3483301b35fSStefano Zampini     }
349a64f4aa4SStefano Zampini     if (isseqaij) {
350a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
351a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
3523301b35fSStefano Zampini     } else {
353a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
354a64f4aa4SStefano Zampini     }
355a64f4aa4SStefano Zampini   }
3563301b35fSStefano Zampini 
357a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
358a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
359a64f4aa4SStefano Zampini   if (sub_schurs->use_mumps) {
360a64f4aa4SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
361a64f4aa4SStefano Zampini   }
362a64f4aa4SStefano Zampini 
3635a95e1ceSStefano Zampini   /* preliminary checks */
3645a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
3655a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
3665a95e1ceSStefano Zampini   }
3675a95e1ceSStefano Zampini 
3685a95e1ceSStefano Zampini   /* restrict work on active processes */
3695a95e1ceSStefano Zampini   color = 0;
3705a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
3715a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
3725a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
3735a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
3745a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
3755a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
3765a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
3775a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3785a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
3795a95e1ceSStefano Zampini     PetscFunctionReturn(0);
3805a95e1ceSStefano Zampini   }
38181ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
3825a95e1ceSStefano Zampini 
383b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
384883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
385a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
3863301b35fSStefano Zampini     PetscBool isseqsbaij;
387a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
3883301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
3893301b35fSStefano Zampini     if (isseqsbaij) {
390a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
391a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
392a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
393a64f4aa4SStefano Zampini     } else {
394a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
395a64f4aa4SStefano Zampini       A_BB = tA_BB;
396a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
397a64f4aa4SStefano Zampini       A_IB = tA_IB;
398a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
399a64f4aa4SStefano Zampini       A_BI = tA_BI;
400f11841e3SStefano Zampini     }
401a58a30b4SStefano Zampini   } else {
4025a95e1ceSStefano Zampini     A_II = NULL;
4035a95e1ceSStefano Zampini     A_IB = NULL;
4045a95e1ceSStefano Zampini     A_BI = NULL;
4055a95e1ceSStefano Zampini     A_BB = NULL;
406b1b3d7a2SStefano Zampini   }
4075a95e1ceSStefano Zampini   S_all = NULL;
408b1b3d7a2SStefano Zampini 
409b1b3d7a2SStefano Zampini   /* determine interior problems */
4103dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4113dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
412b1b3d7a2SStefano Zampini     PetscBT                touched;
413b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
414b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
415b1b3d7a2SStefano Zampini 
4163dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
4173dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
4183dc780c3SStefano Zampini     }
419b1b3d7a2SStefano Zampini     /* get sizes */
420b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
421b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
422b1b3d7a2SStefano Zampini 
423b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
424b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
425b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
426b1b3d7a2SStefano Zampini 
427b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
428b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
429b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
430b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
431b1b3d7a2SStefano Zampini     }
432b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
433b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
434b1b3d7a2SStefano Zampini 
435b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
436b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
437b1b3d7a2SStefano Zampini     n_prev_added = n_B;
438b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
439b1b3d7a2SStefano Zampini       PetscInt n_added;
440b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
441b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
442b1b3d7a2SStefano 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);
443b1b3d7a2SStefano Zampini       }
444b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
445b1b3d7a2SStefano Zampini       n_prev_added = n_added;
446b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
447b1b3d7a2SStefano Zampini       if (!n_added) break;
448b1b3d7a2SStefano Zampini     }
449b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
450b1b3d7a2SStefano Zampini 
451883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
452a9b99552SStefano 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);
453b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
454a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
455883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
456883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
457b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
458b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
459a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
460b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
461b1b3d7a2SStefano Zampini 
462b1b3d7a2SStefano Zampini       /* II block */
463b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
464b1b3d7a2SStefano Zampini     }
465b1b3d7a2SStefano Zampini   } else {
466b1b3d7a2SStefano Zampini     PetscInt n_I;
467b1b3d7a2SStefano Zampini 
468b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
469b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
470a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
471b1b3d7a2SStefano Zampini 
472b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
473883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
474b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
475b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
476b1b3d7a2SStefano Zampini 
477b1b3d7a2SStefano Zampini       /* II block is the same */
478b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
479b1b3d7a2SStefano Zampini       AE_II = A_II;
480b1b3d7a2SStefano Zampini     }
481b1b3d7a2SStefano Zampini   }
4825a95e1ceSStefano Zampini 
483883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
4845a95e1ceSStefano Zampini   max_subset_size = 0;
485883469d8SStefano Zampini   local_size = 0;
4865a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4875a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4885a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
489883469d8SStefano Zampini     local_size += subset_size;
490883469d8SStefano Zampini   }
491883469d8SStefano Zampini 
492883469d8SStefano Zampini   /* Work arrays for local indices */
493883469d8SStefano Zampini   extra = 0;
494683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
49581ea8064SStefano Zampini   if (sub_schurs->use_mumps && is_I_layer) {
496a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
497883469d8SStefano Zampini   }
498683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
499883469d8SStefano Zampini   if (extra) {
500883469d8SStefano Zampini     const PetscInt *idxs;
501a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
502883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
503a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
504883469d8SStefano Zampini   }
505883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
506dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
507dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
508883469d8SStefano Zampini 
509883469d8SStefano Zampini   /* Get local indices in local numbering */
510883469d8SStefano Zampini   local_size = 0;
5115a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
512883469d8SStefano Zampini     PetscInt j;
513883469d8SStefano Zampini     const    PetscInt *idxs;
514883469d8SStefano Zampini 
5155a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5165a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
517eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
518eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
519eb595f79SStefano Zampini     auxnum2[i] = subset_size;
520883469d8SStefano Zampini     /* subset indices in local numbering */
521883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5225a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
523883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
524883469d8SStefano Zampini     local_size += subset_size;
525883469d8SStefano Zampini   }
526883469d8SStefano Zampini 
5275a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
528d2627357SStefano Zampini   Bwork = NULL;
529d2627357SStefano Zampini   pivots = NULL;
53006a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
531d2627357SStefano Zampini     PetscScalar lwork;
532d2627357SStefano Zampini 
53306a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
534d2627357SStefano Zampini     B_lwork = -1;
535d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
536d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
537d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
538d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
539d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
540d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
541d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
542d2627357SStefano Zampini   }
543d2627357SStefano Zampini 
544d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
545dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
546dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
547dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
548dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
549dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
550dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
551dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
552dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
553dc456d91SStefano Zampini   if (i != local_size) {
554dc456d91SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size);
555eb595f79SStefano Zampini   }
556dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
557e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
558d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
559d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
560d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
561d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
5622972d61bSStefano Zampini 
5635a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
5645a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
5655a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
5665a95e1ceSStefano Zampini 
5675a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
5685a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
5695a95e1ceSStefano Zampini     if (subset_size != local_size) {
5705a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
5715a95e1ceSStefano Zampini     }
5725a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
573b1b3d7a2SStefano Zampini   }
574b1b3d7a2SStefano Zampini 
5755a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
5765a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
5775a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
5785a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
5795a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
5805a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
581aa83b6aeSStefano Zampini   }
582b1b3d7a2SStefano Zampini 
5835a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
584be83ff47SStefano Zampini   F = NULL;
585683d3df6SStefano 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 */
5865a95e1ceSStefano Zampini     Mat         S_Ej_expl;
5875a95e1ceSStefano Zampini     PetscScalar *work;
5885a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
5895a95e1ceSStefano Zampini     PetscBool   Sdense;
5905a95e1ceSStefano Zampini 
5915a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
5925a95e1ceSStefano Zampini     local_size = 0;
593b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
5945a95e1ceSStefano Zampini       IS  is_subset_B;
5955a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
5965a95e1ceSStefano Zampini 
5975a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
5985a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
5995a95e1ceSStefano Zampini       /* EE block */
6005a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
6015a95e1ceSStefano Zampini       /* IE block */
6025a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
6035a95e1ceSStefano Zampini       /* EI block */
6045a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
6055a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
6065a95e1ceSStefano Zampini       } else {
6075a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
6085a95e1ceSStefano Zampini       }
609a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
6105a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
6115a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
6125a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
6135a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
614b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
615b1b3d7a2SStefano Zampini         KSP ksp;
616b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
6175a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
618b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
619b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
620b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
621b1b3d7a2SStefano Zampini         KSPType   ksp_type;
622b1b3d7a2SStefano Zampini         PetscInt  n_internal;
6235a95e1ceSStefano Zampini         PetscBool ispcnone;
624b1b3d7a2SStefano Zampini 
625b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
6265a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
627b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
628b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
629b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
630b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
6315a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
6325a95e1ceSStefano Zampini         if (!ispcnone) {
6335a95e1ceSStefano Zampini           PCType pc_type;
634b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
635b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
6365a95e1ceSStefano Zampini         } else {
6375a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
6385a95e1ceSStefano Zampini         }
639b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
640b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
641b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
642b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
643b1b3d7a2SStefano Zampini           if (solver) {
644b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
645b1b3d7a2SStefano Zampini           }
646b1b3d7a2SStefano Zampini         }
647b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
648b1b3d7a2SStefano Zampini       }
6495a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6505a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
6515a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
6525a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
6535a95e1ceSStefano Zampini       if (Sdense) {
6545a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
6555a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
656b1b3d7a2SStefano Zampini         }
6575a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6585a95e1ceSStefano Zampini       } else {
6595a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
6605a95e1ceSStefano Zampini       }
6615a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
662a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
6635a95e1ceSStefano Zampini       local_size += subset_size;
6645a95e1ceSStefano Zampini     }
6655a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
666b1b3d7a2SStefano Zampini     /* free */
667b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
668b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
6695a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
670883469d8SStefano Zampini   } else {
671be83ff47SStefano Zampini     Mat         A;
672883469d8SStefano Zampini     IS          is_A_all;
673683d3df6SStefano Zampini     PetscScalar *work,*S_data,*schur_factor;
6746dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
675683d3df6SStefano Zampini     PetscBool   economic,mumps_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
676d4933d67SStefano Zampini     char        solver[256];
677883469d8SStefano Zampini 
678683d3df6SStefano Zampini     /* get sizes */
67981ea8064SStefano Zampini     n_I = 0;
68081ea8064SStefano Zampini     if (is_I_layer) {
681a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
68281ea8064SStefano Zampini     }
683683d3df6SStefano Zampini     economic = PETSC_FALSE;
684683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
685683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
686683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
687d62866d3SStefano Zampini 
688683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
689683d3df6SStefano Zampini     size_active_schur = local_size;
690683d3df6SStefano Zampini     cum = n_I+size_active_schur;
691683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
692683d3df6SStefano Zampini       const PetscInt* idxs;
693683d3df6SStefano Zampini       PetscInt        n_dir;
694683d3df6SStefano Zampini 
695683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
696683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
697683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
698683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
699683d3df6SStefano Zampini       cum += n_dir;
700d62866d3SStefano Zampini     }
701683d3df6SStefano Zampini     factor_workaround = PETSC_FALSE;
702683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
703683d3df6SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && compute_Stilda) {
704683d3df6SStefano Zampini       PetscInt n_v;
705683d3df6SStefano Zampini 
706683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
707683d3df6SStefano Zampini       if (n_v) {
708683d3df6SStefano Zampini         const PetscInt* idxs;
709683d3df6SStefano Zampini 
710683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
711683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
712683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
713683d3df6SStefano Zampini         cum += n_v;
714683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
715683d3df6SStefano Zampini       }
716683d3df6SStefano Zampini     }
717683d3df6SStefano Zampini     size_schur = cum - n_I;
718683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
719683d3df6SStefano Zampini     /* get working mat (TODO: factorize without actually permuting)  */
720683d3df6SStefano Zampini     if (cum == n) {
721683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
722683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
723683d3df6SStefano Zampini     } else {
7246816873aSStefano Zampini       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
725683d3df6SStefano Zampini     }
726a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
7273301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
7283301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
7293301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
730d4933d67SStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
731d4933d67SStefano Zampini     ierr = PetscOptionsGetString("sub_schurs_","-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
732883469d8SStefano Zampini 
733683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
734d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
735d47842beSStefano Zampini 
736683d3df6SStefano Zampini     if (n_I) { /* TODO add ordering from options */
7375a05ddb0SStefano Zampini       IS is_schur;
7385a05ddb0SStefano Zampini 
7399ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
740d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
741883469d8SStefano Zampini       } else {
742d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
743883469d8SStefano Zampini       }
744883469d8SStefano Zampini       /* subsets ordered last */
7455a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
7465a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
7475a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
748883469d8SStefano Zampini 
749883469d8SStefano Zampini       /* factorization step */
7509ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
751883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
752be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
753be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
754be83ff47SStefano Zampini #endif
755*e8ade678SStefano Zampini         if (sub_schurs->is_posdef) {
756*e8ade678SStefano Zampini           ierr = MatFactorSetSchurComplementSolverType(F,1);CHKERRQ(ierr);
757*e8ade678SStefano Zampini         } else {
758*e8ade678SStefano Zampini           ierr = MatFactorSetSchurComplementSolverType(F,2);CHKERRQ(ierr);
759*e8ade678SStefano Zampini         }
760883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
761a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
762883469d8SStefano Zampini       } else {
763883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
764be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
765be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
766be83ff47SStefano Zampini #endif
767883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
768883469d8SStefano Zampini       }
769883469d8SStefano Zampini 
770883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
7715a05ddb0SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
772d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
773683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
774683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
7754a6c6b0dSStefano Zampini       mumps_S = PETSC_TRUE;
776be83ff47SStefano Zampini     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
777be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
778683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
7794a6c6b0dSStefano Zampini       mumps_S = PETSC_FALSE;
780be83ff47SStefano Zampini     }
781be83ff47SStefano Zampini 
782be83ff47SStefano Zampini     if (reuse_solvers) {
783d62866d3SStefano Zampini       Mat              A_II;
78453892102SStefano Zampini       Vec              vec1_B;
785d62866d3SStefano Zampini       PCBDDCReuseMumps msolv_ctx;
786d5574798SStefano Zampini 
787d62866d3SStefano Zampini       if (sub_schurs->reuse_mumps) {
7886816873aSStefano Zampini         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
789e28d306cSStefano Zampini       } else {
790e28d306cSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
791d62866d3SStefano Zampini       }
792e28d306cSStefano Zampini       msolv_ctx = sub_schurs->reuse_mumps;
793be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
794d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
795d5574798SStefano Zampini       msolv_ctx->F = F;
796683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
797683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
798683d3df6SStefano Zampini       {
799683d3df6SStefano Zampini         PetscScalar *array;
800683d3df6SStefano Zampini         PetscInt    n;
801683d3df6SStefano Zampini 
802683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
803683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
804683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
805683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
806683d3df6SStefano Zampini       }
807683d3df6SStefano Zampini       msolv_ctx->has_vertices = factor_workaround;
808d62866d3SStefano Zampini 
809d62866d3SStefano Zampini       /* interior solver */
810683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
811d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
812d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
813d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
814d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
815e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
816d62866d3SStefano Zampini 
817d62866d3SStefano Zampini       /* correction solver */
818683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A),&msolv_ctx->correction_solver);CHKERRQ(ierr);
819d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
820d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
821d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
822d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
823e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
82453892102SStefano Zampini 
82553892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
82653892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
82753892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
828683d3df6SStefano Zampini       if (!factor_workaround) {
82953892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
83053892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
83153892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
83253892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
833683d3df6SStefano Zampini       } else {
834683d3df6SStefano Zampini         IS              is_B_all;
835683d3df6SStefano Zampini         const PetscInt* idxs;
836683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
837683d3df6SStefano Zampini 
838683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
839683d3df6SStefano Zampini         dual = size_schur - n_v;
840683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
841683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
842683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
843683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
844683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
845683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
846683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
847683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
848683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
849683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
850683d3df6SStefano Zampini       }
851683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
852d5574798SStefano Zampini     }
85308122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
85453892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
8555db18549SStefano Zampini 
856be83ff47SStefano Zampini     /* Work arrays */
857be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
8585a94c880SStefano Zampini 
8595a94c880SStefano Zampini     /* matrices for adaptive selection */
86012d906b1SStefano Zampini     if (compute_Stilda) {
8615a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
8625a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
8635a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
8645a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
8655a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
8665a94c880SStefano Zampini       }
8675a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all) {
8685a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
8695a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
8705a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
8715a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
8725a94c880SStefano Zampini       }
87312d906b1SStefano Zampini     }
874d2627357SStefano Zampini 
875be83ff47SStefano Zampini     /* S_Ej_all */
876be83ff47SStefano Zampini     cum = cum2 = 0;
877be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
87865d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
87965d8bf0aSStefano Zampini       PetscInt j;
88065d8bf0aSStefano Zampini 
8815a95e1ceSStefano Zampini       /* get S_E */
882b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
883683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
884be83ff47SStefano Zampini         PetscInt k;
885be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
886be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
887be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
888683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
889be83ff47SStefano Zampini           }
890be83ff47SStefano Zampini         }
89106a4e24aSStefano Zampini       } else { /* just copy to workspace */
892be83ff47SStefano Zampini         PetscInt k;
893be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
894be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
895be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
896be83ff47SStefano Zampini           }
897be83ff47SStefano Zampini         }
8989087bf02SStefano Zampini       }
8995a95e1ceSStefano Zampini       /* insert S_E values */
900be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
9015a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
902a1337663SStefano Zampini 
903683d3df6SStefano Zampini       /* if adaptivity is requested, invert S_E blocks */
904862806e4SStefano Zampini       if (compute_Stilda) {
90508122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
90608122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
90706a4e24aSStefano Zampini         if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
9085a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
9092972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
9105a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
9112972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
912d6462365SStefano Zampini         } else {
913d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
914d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
915d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
916d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
9172972d61bSStefano Zampini         }
91808122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
9195a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
9209087bf02SStefano Zampini       }
921be83ff47SStefano Zampini       cum += subset_size;
922be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
923be83ff47SStefano Zampini     }
924be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
9254a6c6b0dSStefano Zampini 
9264a6c6b0dSStefano Zampini     if (mumps_S) {
9275a05ddb0SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
9284a6c6b0dSStefano Zampini     }
929683d3df6SStefano Zampini 
930683d3df6SStefano Zampini     schur_factor = NULL;
93145951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
932683d3df6SStefano Zampini 
9334a6c6b0dSStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
9344a6c6b0dSStefano Zampini         PetscInt j;
9354a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
9365a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
9374a6c6b0dSStefano Zampini       } else {
938683d3df6SStefano Zampini         Mat S_all_inv=NULL;
9395a05ddb0SStefano Zampini         if (mumps_S) { /* use MatFactor calls to invert S */
9406dba178dSStefano Zampini             /* let MatFactor factorize the Schur complement */
9416dba178dSStefano Zampini           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
942683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
943683d3df6SStefano 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 */
944683d3df6SStefano Zampini           if (factor_workaround) {
945683d3df6SStefano Zampini             PetscScalar *data;
946683d3df6SStefano Zampini             PetscInt nd = 0;
9476dba178dSStefano Zampini 
948683d3df6SStefano Zampini             if (!sub_schurs->is_posdef) {
949683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
950683d3df6SStefano Zampini             }
9516dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
952683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
953683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
954683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
955683d3df6SStefano Zampini             }
956683d3df6SStefano Zampini             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
957683d3df6SStefano Zampini             cum = 0;
958683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
959683d3df6SStefano Zampini               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
960683d3df6SStefano Zampini               cum += size_active_schur-i;
961683d3df6SStefano Zampini             }
962683d3df6SStefano Zampini             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
9636dba178dSStefano Zampini             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
964683d3df6SStefano Zampini             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
965683d3df6SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
966683d3df6SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
967683d3df6SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
968683d3df6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
969683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
970683d3df6SStefano Zampini           } else {
9716dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
9726dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
973683d3df6SStefano Zampini           }
9744a6c6b0dSStefano Zampini         } else { /* we need to invert explicitly since we are not using MUMPS for S */
975683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
976683d3df6SStefano Zampini           S_all_inv = S_all;
977683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
978be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
979be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
98006a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
981be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
982be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
983be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
984be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
985d6462365SStefano Zampini           } else {
986d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
987d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
988d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
989d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
990be83ff47SStefano Zampini           }
991be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
992683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
993be83ff47SStefano Zampini         }
994be83ff47SStefano Zampini         /* S_Ej_tilda_all */
995be83ff47SStefano Zampini         cum = cum2 = 0;
996683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
997be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
998be83ff47SStefano Zampini           PetscInt j;
999862806e4SStefano Zampini 
1000862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1001be83ff47SStefano Zampini           /* get (St^-1)_E */
100206a4e24aSStefano Zampini           /* Here I don't need to expand to upper triangular since St^-1
100306a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1004a0b0af32SStefano Zampini           if (S_lower_triangular) {
1005be83ff47SStefano Zampini             PetscInt k;
1006be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1007be83ff47SStefano Zampini               for (j=k;j<subset_size;j++) {
1008be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1009be83ff47SStefano Zampini               }
1010be83ff47SStefano Zampini             }
1011be83ff47SStefano Zampini           } else { /* copy to workspace */
1012be83ff47SStefano Zampini             PetscInt k;
1013be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1014be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1015be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1016be83ff47SStefano Zampini               }
1017be83ff47SStefano Zampini             }
1018be83ff47SStefano Zampini           }
1019be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
10205a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1021be83ff47SStefano Zampini           cum += subset_size;
1022be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1023883469d8SStefano Zampini         }
1024683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
10254a6c6b0dSStefano Zampini         if (mumps_S) {
10266dba178dSStefano Zampini           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
10275db18549SStefano Zampini         }
1028683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1029683d3df6SStefano Zampini       }
1030683d3df6SStefano Zampini 
1031683d3df6SStefano Zampini       /* move back factors to Schur data into F */
1032683d3df6SStefano Zampini       if (factor_workaround) {
1033683d3df6SStefano Zampini         Mat S_tmp;
1034683d3df6SStefano Zampini         PetscInt nv,nd=0;
1035683d3df6SStefano Zampini 
1036683d3df6SStefano Zampini         if (!mumps_S) {
1037683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1038683d3df6SStefano Zampini         }
1039683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
10406dba178dSStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1041683d3df6SStefano Zampini         if (sub_schurs->is_posdef) {
1042683d3df6SStefano Zampini           PetscScalar *data;
1043683d3df6SStefano Zampini 
1044683d3df6SStefano Zampini           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1045683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1046683d3df6SStefano Zampini 
1047683d3df6SStefano Zampini           if (S_lower_triangular) {
1048683d3df6SStefano Zampini             cum = 0;
1049683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1050683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1051683d3df6SStefano Zampini               cum += size_active_schur-i;
1052683d3df6SStefano Zampini 	    }
1053683d3df6SStefano Zampini           } else {
1054683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1055683d3df6SStefano Zampini           }
1056683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1057683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1058683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1059683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1060683d3df6SStefano Zampini 	    }
1061683d3df6SStefano Zampini           }
10626dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1063683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1064683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
1065683d3df6SStefano Zampini             data[i*(size_schur+1)] = PETSC_MAX_REAL;
1066683d3df6SStefano Zampini           }
1067683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1068683d3df6SStefano Zampini         } else {
1069683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1070683d3df6SStefano Zampini         }
10716dba178dSStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
10729087bf02SStefano Zampini       }
10734a6c6b0dSStefano Zampini     }
10744a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1075683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
10764a6c6b0dSStefano Zampini   }
10775a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1078be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1079a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1080a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1081a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1082a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1083a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
10845db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10855db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10865a95e1ceSStefano Zampini   if (compute_Stilda) {
10875a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10885a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10895a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10905a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109108122e43SStefano Zampini   }
1092a1337663SStefano Zampini 
10935db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
10945db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
10953927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
10963927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
10975a95e1ceSStefano Zampini 
10985db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
10995a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
1100dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
11015a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
11025a94c880SStefano Zampini   } else {
11035a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
11045a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
1105dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
11065a94c880SStefano Zampini   }
110708122e43SStefano Zampini 
1108ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
1109ac632422SStefano Zampini   if (faster_deluxe) {
11105a95e1ceSStefano Zampini     Mat         tmpmat;
11115a95e1ceSStefano Zampini     PetscScalar *array;
11125a95e1ceSStefano Zampini     PetscInt    cum;
11135a95e1ceSStefano Zampini 
11145a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
11155a95e1ceSStefano Zampini     cum = 0;
11165a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
11175a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
11185a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
11195a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
112006a4e24aSStefano Zampini       if (!use_getr) {
11215a95e1ceSStefano Zampini         PetscInt j,k;
11225a95e1ceSStefano Zampini 
11235a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
11245a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
11255a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
11265a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
11275a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
11285a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
11295a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
11305a95e1ceSStefano Zampini           }
11315a95e1ceSStefano Zampini         }
1132d6462365SStefano Zampini       } else {
1133d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1134d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1135d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1136d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
11375a95e1ceSStefano Zampini       }
11385a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
11395a95e1ceSStefano Zampini       cum += subset_size*subset_size;
11405a95e1ceSStefano Zampini     }
11415a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
11425a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
11435a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1144ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
11455a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
11465a95e1ceSStefano Zampini   }
11475a95e1ceSStefano Zampini 
1148f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
11495a95e1ceSStefano Zampini   if (compute_Stilda) {
11505a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1151a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
11525a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1153dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
11545a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
115508122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
11565a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_inv_all;
1157dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
115808122e43SStefano Zampini   }
11593202ece2SStefano Zampini 
11605a95e1ceSStefano Zampini   /* free workspace */
11615a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
116206a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1163a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
11643202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1165dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
11665a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1167b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1168b1b3d7a2SStefano Zampini }
1169b1b3d7a2SStefano Zampini 
1170b1b3d7a2SStefano Zampini #undef __FUNCT__
1171b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
1172a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1173b1b3d7a2SStefano Zampini {
11749bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
11755a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1176b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1177b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1178b1b3d7a2SStefano Zampini 
1179b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1180b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1181b1b3d7a2SStefano Zampini   if (!is_sorted) {
1182b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1183b1b3d7a2SStefano Zampini   }
1184b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1185b1b3d7a2SStefano Zampini   if (!is_sorted) {
1186b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1187b1b3d7a2SStefano Zampini   }
1188b1b3d7a2SStefano Zampini 
1189b1b3d7a2SStefano Zampini   /* reset any previous data */
1190b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1191b1b3d7a2SStefano Zampini 
11925a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
11939bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1194b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
119508122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1196b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1197b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
1198b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
1199b1b3d7a2SStefano Zampini   }
1200b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
1201b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
120208122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1203b1b3d7a2SStefano Zampini   }
1204b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
1205b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
1206d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1207d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1208b1b3d7a2SStefano Zampini 
12096816873aSStefano Zampini   /* Determine if MUMPS can be used */
1210883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
1211883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1212a64f4aa4SStefano Zampini   sub_schurs->use_mumps = PETSC_TRUE;
1213883469d8SStefano Zampini #endif
1214b1b3d7a2SStefano Zampini 
1215b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1216b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1217b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1218b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
12195db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
12205db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
12215db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
12225db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
12235a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1224b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1225a64f4aa4SStefano Zampini   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1226b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1227b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1228b96c3477SStefano Zampini     }
12299bb4a8caSStefano Zampini   }
1230d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
1231b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1232b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
123308122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1234b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1235b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1236b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1237b1b3d7a2SStefano Zampini }
1238b1b3d7a2SStefano Zampini 
1239b1b3d7a2SStefano Zampini #undef __FUNCT__
124034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
124134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
124234a97f8cSStefano Zampini {
124334a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
124434a97f8cSStefano Zampini   PetscErrorCode  ierr;
124534a97f8cSStefano Zampini 
124634a97f8cSStefano Zampini   PetscFunctionBegin;
124734a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
12485ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
124934a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
125034a97f8cSStefano Zampini   PetscFunctionReturn(0);
125134a97f8cSStefano Zampini }
125234a97f8cSStefano Zampini 
125334a97f8cSStefano Zampini #undef __FUNCT__
125434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
125534a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
125634a97f8cSStefano Zampini {
125734a97f8cSStefano Zampini   PetscErrorCode ierr;
125834a97f8cSStefano Zampini 
125934a97f8cSStefano Zampini   PetscFunctionBegin;
126034a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
126134a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
126234a97f8cSStefano Zampini   PetscFunctionReturn(0);
126334a97f8cSStefano Zampini }
126434a97f8cSStefano Zampini 
126534a97f8cSStefano Zampini #undef __FUNCT__
126634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
126734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
126834a97f8cSStefano Zampini {
126934a97f8cSStefano Zampini   PetscInt       i;
127034a97f8cSStefano Zampini   PetscErrorCode ierr;
127134a97f8cSStefano Zampini 
127234a97f8cSStefano Zampini   PetscFunctionBegin;
12731e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1274b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1275b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1276b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
12775db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
12785db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
127941c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
128041c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
128108122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1282a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
12835db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1284d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1285d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
128608122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
128734a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1288b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
128934a97f8cSStefano Zampini   }
12905ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1291b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
12923dc780c3SStefano Zampini   }
1293d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps) {
1294d62866d3SStefano Zampini     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1295d62866d3SStefano Zampini   }
1296d62866d3SStefano Zampini   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
129734a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
129834a97f8cSStefano Zampini   PetscFunctionReturn(0);
129934a97f8cSStefano Zampini }
130034a97f8cSStefano Zampini 
130134a97f8cSStefano Zampini #undef __FUNCT__
130234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
13032a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
130434a97f8cSStefano Zampini {
130534a97f8cSStefano Zampini   PetscInt       i,j,n;
130634a97f8cSStefano Zampini   PetscErrorCode ierr;
130734a97f8cSStefano Zampini 
130834a97f8cSStefano Zampini   PetscFunctionBegin;
130934a97f8cSStefano Zampini   n = 0;
131034a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
131134a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
131234a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
131334a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
131434a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
131534a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
131634a97f8cSStefano Zampini         queue_tip[n] = dof;
131734a97f8cSStefano Zampini         n++;
131834a97f8cSStefano Zampini       }
131934a97f8cSStefano Zampini     }
132034a97f8cSStefano Zampini   }
132134a97f8cSStefano Zampini   *n_added = n;
132234a97f8cSStefano Zampini   PetscFunctionReturn(0);
132334a97f8cSStefano Zampini }
1324