xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision be83ff473a1423cc2ab5775993c40ed29465416d)
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__
11d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
12d62866d3SStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
13d62866d3SStefano Zampini {
14d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
15d62866d3SStefano Zampini   PetscInt         ival;
16d62866d3SStefano Zampini   PetscErrorCode   ierr;
17d62866d3SStefano Zampini 
18d62866d3SStefano Zampini   PetscFunctionBegin;
19d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
20d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
21d62866d3SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
22d62866d3SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
23d62866d3SStefano Zampini #endif
246816873aSStefano Zampini   ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
25d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
26d62866d3SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
27d62866d3SStefano Zampini #endif
28d62866d3SStefano Zampini   PetscFunctionReturn(0);
29d62866d3SStefano Zampini }
30d62866d3SStefano Zampini 
31d62866d3SStefano Zampini #undef __FUNCT__
32d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCReuseMumpsReset"
33d62866d3SStefano Zampini static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
34d62866d3SStefano Zampini {
35d62866d3SStefano Zampini   PetscErrorCode ierr;
36d62866d3SStefano Zampini 
37d62866d3SStefano Zampini   PetscFunctionBegin;
38d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
39d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
40d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
41d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
42d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
43d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
446816873aSStefano Zampini   ierr = VecDestroy(&reuse->solB);CHKERRQ(ierr);
456816873aSStefano Zampini   ierr = VecDestroy(&reuse->rhsB);CHKERRQ(ierr);
46d62866d3SStefano Zampini   PetscFunctionReturn(0);
47d62866d3SStefano Zampini }
48d5574798SStefano Zampini 
49d5574798SStefano Zampini #undef __FUNCT__
50d5574798SStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
51d5574798SStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
52d5574798SStefano Zampini {
53d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
54d5574798SStefano Zampini   PetscScalar      *array,*array_mumps;
55d5574798SStefano Zampini   PetscInt         ival;
56d5574798SStefano Zampini   PetscErrorCode   ierr;
57d5574798SStefano Zampini 
58d5574798SStefano Zampini   PetscFunctionBegin;
59d5574798SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
60d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
61d5574798SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
62d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
63d62866d3SStefano Zampini #endif
64d5574798SStefano Zampini   /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */
65d5574798SStefano Zampini   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
66d5574798SStefano Zampini   ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
67d62866d3SStefano Zampini   ierr = PetscMemcpy(array_mumps,array,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
68d5574798SStefano Zampini   ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
69d5574798SStefano Zampini   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
70d5574798SStefano Zampini 
71d5574798SStefano Zampini   ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
72d5574798SStefano Zampini 
73d5574798SStefano Zampini   /* get back data to caller worskpace */
74d5574798SStefano Zampini   ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
75d5574798SStefano Zampini   ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
76d62866d3SStefano Zampini   ierr = PetscMemcpy(array,array_mumps,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
77d5574798SStefano Zampini   ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
78d5574798SStefano Zampini   ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
79d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
80d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
81d62866d3SStefano Zampini #endif
82d5574798SStefano Zampini   PetscFunctionReturn(0);
83d5574798SStefano Zampini }
843202ece2SStefano Zampini 
853202ece2SStefano Zampini #undef __FUNCT__
863202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
875ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
883202ece2SStefano Zampini {
893202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
903202ece2SStefano Zampini   KSP            ksp;
913202ece2SStefano Zampini   PC             pc;
923202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
933202ece2SStefano Zampini   PetscReal      fill = 2.0;
94f11841e3SStefano Zampini   PetscInt       n_I;
953202ece2SStefano Zampini   PetscMPIInt    size;
963202ece2SStefano Zampini   PetscErrorCode ierr;
973202ece2SStefano Zampini 
983202ece2SStefano Zampini   PetscFunctionBegin;
993202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
1003202ece2SStefano Zampini   if (size != 1) {
1013202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
1023202ece2SStefano Zampini   }
103f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
104f11841e3SStefano Zampini     PetscBool Sdense;
105f11841e3SStefano Zampini 
106f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
107f11841e3SStefano Zampini     if (!Sdense) {
108f11841e3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
109f11841e3SStefano Zampini     }
110f11841e3SStefano Zampini   }
1113202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
1123202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
1133202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1143202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
1153202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
1163202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
1173202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
1183202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
119f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
120f11841e3SStefano Zampini   if (n_I) {
1213202ece2SStefano Zampini     if (!Bdense) {
1223202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
1233202ece2SStefano Zampini     } else {
1243202ece2SStefano Zampini       Bd = B;
1253202ece2SStefano Zampini     }
1263202ece2SStefano Zampini 
1273202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1283202ece2SStefano Zampini       Mat fact;
1293202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1303202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1313202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1323202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1333202ece2SStefano Zampini     } else {
13407b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
13507b1e237SStefano Zampini 
13607b1e237SStefano Zampini       if (ex) {
1373202ece2SStefano Zampini         Mat Ainvd;
1383202ece2SStefano Zampini 
1393202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1403202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1413202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
14207b1e237SStefano Zampini       } else {
14307b1e237SStefano Zampini         Vec         sol,rhs;
14407b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
14507b1e237SStefano Zampini         PetscInt    i,nrhs,n;
14607b1e237SStefano Zampini 
14707b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
14807b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
14907b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
15007b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
15107b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
15207b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
15307b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
15407b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
15507b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
15607b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
15707b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
15807b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
15907b1e237SStefano Zampini         }
16007b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
16107b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
16207b1e237SStefano Zampini       }
1633202ece2SStefano Zampini     }
1645ec10c6aSStefano Zampini     if (!Bdense & !issym) {
1653202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
1663202ece2SStefano Zampini     }
1675ec10c6aSStefano Zampini 
1685ec10c6aSStefano Zampini     if (!issym) {
1693202ece2SStefano Zampini       if (!Cdense) {
1703202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
1713202ece2SStefano Zampini       } else {
1723202ece2SStefano Zampini         Cd = C;
1733202ece2SStefano Zampini       }
1745ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
1753202ece2SStefano Zampini       if (!Cdense) {
1763202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
1773202ece2SStefano Zampini       }
1785ec10c6aSStefano Zampini     } else {
1795ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
1805ec10c6aSStefano Zampini       if (!Bdense) {
1815ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
1825ec10c6aSStefano Zampini       }
1835ec10c6aSStefano Zampini     }
1845ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
185f11841e3SStefano Zampini   }
1863202ece2SStefano Zampini 
1873202ece2SStefano Zampini   if (D) {
1883202ece2SStefano Zampini     Mat       Dd;
1893202ece2SStefano Zampini     PetscBool Ddense;
1903202ece2SStefano Zampini 
1913202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
1923202ece2SStefano Zampini     if (!Ddense) {
1933202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
1943202ece2SStefano Zampini     } else {
1953202ece2SStefano Zampini       Dd = D;
1963202ece2SStefano Zampini     }
197f11841e3SStefano Zampini     if (n_I) {
1983202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
199f11841e3SStefano Zampini     } else {
200f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
201f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
202f11841e3SStefano Zampini       } else {
203f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
204f11841e3SStefano Zampini       }
205f11841e3SStefano Zampini     }
2063202ece2SStefano Zampini     if (!Ddense) {
2073202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
2083202ece2SStefano Zampini     }
2093202ece2SStefano Zampini   } else {
2103202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
2113202ece2SStefano Zampini   }
2123202ece2SStefano Zampini   PetscFunctionReturn(0);
2133202ece2SStefano Zampini }
21434a97f8cSStefano Zampini 
21534a97f8cSStefano Zampini #undef __FUNCT__
2161580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
2176816873aSStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers,PetscBool use_edges, PetscBool use_faces)
218b1b3d7a2SStefano Zampini {
219*be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
220*be83ff47SStefano Zampini   Mat                    S_all;
221d2627357SStefano Zampini   Mat                    global_schur_subsets,work_mat;
22208122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
2235db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
224a9b99552SStefano Zampini   IS                     is_I,is_I_layer,temp_is;
225d648f858SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
226eb595f79SStefano Zampini   PetscInt               *auxnum1,*auxnum2,*all_local_idx_G_rep;
2275a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
228883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
22908122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
23006a4b1faSStefano Zampini   PetscScalar            *Bwork;
2315a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2325a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2335a95e1ceSStefano Zampini   MPI_Comm               comm_n;
234b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
235b1b3d7a2SStefano Zampini 
236b1b3d7a2SStefano Zampini   PetscFunctionBegin;
237a64f4aa4SStefano Zampini   /* update info in sub_schurs */
238a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
239a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
240a64f4aa4SStefano Zampini   if (Ain) {
241a64f4aa4SStefano Zampini     PetscBool isseqaij;
242a64f4aa4SStefano Zampini 
243a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
244a64f4aa4SStefano Zampini     if (isseqaij) {
245a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
246a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
247a64f4aa4SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
248a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
249a64f4aa4SStefano Zampini     }
250a64f4aa4SStefano Zampini   }
251a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
252a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
253a64f4aa4SStefano Zampini   if (sub_schurs->use_mumps) {
254a64f4aa4SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
255a64f4aa4SStefano Zampini   }
256a64f4aa4SStefano Zampini 
2575a95e1ceSStefano Zampini   /* preliminary checks */
2585a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
2595a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
2605a95e1ceSStefano Zampini   }
2615a95e1ceSStefano Zampini   /* determine if we are dealing with hermitian positive definite problems */
2625a95e1ceSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2635a95e1ceSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
2645a95e1ceSStefano Zampini   if (sub_schurs->A) {
2655a95e1ceSStefano Zampini     PetscInt lsize;
2665a95e1ceSStefano Zampini 
2675a95e1ceSStefano Zampini     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
2685a95e1ceSStefano Zampini     if (lsize) {
2695a95e1ceSStefano Zampini       ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
2705a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
2715a95e1ceSStefano Zampini         PetscScalar val;
2725a95e1ceSStefano Zampini         Vec         vec1,vec2;
2735a95e1ceSStefano Zampini 
2745a95e1ceSStefano Zampini         ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
2755a95e1ceSStefano Zampini         ierr = VecSetRandom(vec1,NULL);
2765a95e1ceSStefano Zampini         ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
2775a95e1ceSStefano Zampini         ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
2785a95e1ceSStefano Zampini         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
2795a95e1ceSStefano Zampini         if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
2805a95e1ceSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
2815a95e1ceSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
2825a95e1ceSStefano Zampini       }
2835a95e1ceSStefano Zampini     } else {
2845a95e1ceSStefano Zampini       sub_schurs->is_hermitian = PETSC_TRUE;
2855a95e1ceSStefano Zampini       sub_schurs->is_posdef = PETSC_TRUE;
2865a95e1ceSStefano Zampini     }
2875a95e1ceSStefano Zampini     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
2885a95e1ceSStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported");
2895a95e1ceSStefano Zampini     }
2905a95e1ceSStefano Zampini   }
2915a95e1ceSStefano Zampini   /* restrict work on active processes */
2925a95e1ceSStefano Zampini   color = 0;
2935a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
2945a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
2955a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
2965a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
2975a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
2985a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
2995a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
3005a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3015a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
3025a95e1ceSStefano Zampini     PetscFunctionReturn(0);
3035a95e1ceSStefano Zampini   }
3045a95e1ceSStefano Zampini 
305b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
306883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
307a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
308f11841e3SStefano Zampini     PetscBool isseqaij;
309a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
310a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
311f11841e3SStefano Zampini     if (!isseqaij) {
312a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
313a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
314a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
315a64f4aa4SStefano Zampini     } else {
316a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
317a64f4aa4SStefano Zampini       A_BB = tA_BB;
318a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
319a64f4aa4SStefano Zampini       A_IB = tA_IB;
320a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
321a64f4aa4SStefano Zampini       A_BI = tA_BI;
322f11841e3SStefano Zampini     }
323a58a30b4SStefano Zampini   } else {
3245a95e1ceSStefano Zampini     A_II = NULL;
3255a95e1ceSStefano Zampini     A_IB = NULL;
3265a95e1ceSStefano Zampini     A_BI = NULL;
3275a95e1ceSStefano Zampini     A_BB = NULL;
328b1b3d7a2SStefano Zampini   }
3295a95e1ceSStefano Zampini   S_all = NULL;
3305a95e1ceSStefano Zampini   S_Ej_tilda_all = NULL;
3315a95e1ceSStefano Zampini   S_Ej_inv_all = NULL;
332b1b3d7a2SStefano Zampini 
333b1b3d7a2SStefano Zampini   /* determine interior problems */
3343dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
3353dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
336b1b3d7a2SStefano Zampini     PetscBT                touched;
337b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
338b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
339b1b3d7a2SStefano Zampini 
3403dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
3413dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
3423dc780c3SStefano Zampini     }
343b1b3d7a2SStefano Zampini     /* get sizes */
344b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
345b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
346b1b3d7a2SStefano Zampini 
347b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
348b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
349b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
350b1b3d7a2SStefano Zampini 
351b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
352b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
353b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
354b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
355b1b3d7a2SStefano Zampini     }
356b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
357b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
358b1b3d7a2SStefano Zampini 
359b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
360b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
361b1b3d7a2SStefano Zampini     n_prev_added = n_B;
362b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
363b1b3d7a2SStefano Zampini       PetscInt n_added;
364b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
365b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
366b1b3d7a2SStefano 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);
367b1b3d7a2SStefano Zampini       }
368b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
369b1b3d7a2SStefano Zampini       n_prev_added = n_added;
370b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
371b1b3d7a2SStefano Zampini       if (!n_added) break;
372b1b3d7a2SStefano Zampini     }
373b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
374b1b3d7a2SStefano Zampini 
375883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
376a9b99552SStefano 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);
377b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
378a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
379883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
380883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
381b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
382b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
383a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
384b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
385b1b3d7a2SStefano Zampini 
386b1b3d7a2SStefano Zampini       /* II block */
387b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
388b1b3d7a2SStefano Zampini     }
389b1b3d7a2SStefano Zampini   } else {
390b1b3d7a2SStefano Zampini     PetscInt n_I;
391b1b3d7a2SStefano Zampini 
392b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
393b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
394a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
395b1b3d7a2SStefano Zampini 
396b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
397883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
398b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
399b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
400b1b3d7a2SStefano Zampini 
401b1b3d7a2SStefano Zampini       /* II block is the same */
402b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
403b1b3d7a2SStefano Zampini       AE_II = A_II;
404b1b3d7a2SStefano Zampini     }
405b1b3d7a2SStefano Zampini   }
4065a95e1ceSStefano Zampini 
407883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
4085a95e1ceSStefano Zampini   max_subset_size = 0;
409883469d8SStefano Zampini   local_size = 0;
4105a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4115a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4125a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
413883469d8SStefano Zampini     local_size += subset_size;
414883469d8SStefano Zampini   }
415883469d8SStefano Zampini 
416883469d8SStefano Zampini   /* Work arrays for local indices */
417883469d8SStefano Zampini   extra = 0;
418883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
419a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
420883469d8SStefano Zampini   }
421883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
422883469d8SStefano Zampini   if (extra) {
423883469d8SStefano Zampini     const PetscInt *idxs;
424a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
425883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
426a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
427883469d8SStefano Zampini   }
428883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
429eb595f79SStefano Zampini   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
430883469d8SStefano Zampini 
431883469d8SStefano Zampini   /* Get local indices in local numbering */
432883469d8SStefano Zampini   local_size = 0;
4335a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
434883469d8SStefano Zampini     PetscInt j;
435883469d8SStefano Zampini     const    PetscInt *idxs;
436883469d8SStefano Zampini 
4375a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4385a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
439eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
440eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
441eb595f79SStefano Zampini     auxnum2[i] = subset_size;
442883469d8SStefano Zampini     /* subset indices in local numbering */
443883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
4445a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
445883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
446883469d8SStefano Zampini     local_size += subset_size;
447883469d8SStefano Zampini   }
448883469d8SStefano Zampini 
4495a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
450d2627357SStefano Zampini   Bwork = NULL;
451d2627357SStefano Zampini   pivots = NULL;
4525a95e1ceSStefano Zampini   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
453d2627357SStefano Zampini     PetscScalar lwork;
454d2627357SStefano Zampini 
455d2627357SStefano Zampini     B_lwork = -1;
456d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
457d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
458d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
459d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
460d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
461d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
462d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
463d2627357SStefano Zampini   }
464d2627357SStefano Zampini 
465d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
466eb595f79SStefano Zampini   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
467eb595f79SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
468eb595f79SStefano Zampini   local_size = 0;
469eb595f79SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
470eb595f79SStefano Zampini     PetscInt j;
471eb595f79SStefano Zampini     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
472eb595f79SStefano Zampini   }
473eb595f79SStefano Zampini   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
474eb595f79SStefano Zampini   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
4755a95e1ceSStefano Zampini   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
4765a95e1ceSStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
477d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
478d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
479d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
480d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
4812972d61bSStefano Zampini 
4825a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
4835a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
4845a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
4855a95e1ceSStefano Zampini 
4865a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
4875a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
4885a95e1ceSStefano Zampini     if (subset_size != local_size) {
4895a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
4905a95e1ceSStefano Zampini     }
4915a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
492b1b3d7a2SStefano Zampini   }
493b1b3d7a2SStefano Zampini 
4945a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
4955a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
4965a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
4975a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
4985a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
4995a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
500aa83b6aeSStefano Zampini   } else {
5015a95e1ceSStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
502aa83b6aeSStefano Zampini   }
503b1b3d7a2SStefano Zampini 
5045a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
505*be83ff47SStefano Zampini   F = NULL;
5065a95e1ceSStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
5075a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps) {
5085a95e1ceSStefano Zampini     Mat         S_Ej_expl;
5095a95e1ceSStefano Zampini     PetscScalar *work;
5105a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
5115a95e1ceSStefano Zampini     PetscBool   Sdense;
5125a95e1ceSStefano Zampini 
5135a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
5145a95e1ceSStefano Zampini     local_size = 0;
515b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
5165a95e1ceSStefano Zampini       IS  is_subset_B;
5175a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
5185a95e1ceSStefano Zampini 
5195a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
5205a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
5215a95e1ceSStefano Zampini       /* EE block */
5225a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
5235a95e1ceSStefano Zampini       /* IE block */
5245a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
5255a95e1ceSStefano Zampini       /* EI block */
5265a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
5275a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
5285a95e1ceSStefano Zampini       } else {
5295a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
5305a95e1ceSStefano Zampini       }
531a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
5325a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
5335a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
5345a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
5355a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
536b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
537b1b3d7a2SStefano Zampini         KSP ksp;
538b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
5395a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
540b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
541b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
542b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
543b1b3d7a2SStefano Zampini         KSPType   ksp_type;
544b1b3d7a2SStefano Zampini         PetscInt  n_internal;
5455a95e1ceSStefano Zampini         PetscBool ispcnone;
546b1b3d7a2SStefano Zampini 
547b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
5485a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
549b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
550b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
551b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
552b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
5535a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
5545a95e1ceSStefano Zampini         if (!ispcnone) {
5555a95e1ceSStefano Zampini           PCType pc_type;
556b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
557b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
5585a95e1ceSStefano Zampini         } else {
5595a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
5605a95e1ceSStefano Zampini         }
561b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
562b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
563b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
564b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
565b1b3d7a2SStefano Zampini           if (solver) {
566b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
567b1b3d7a2SStefano Zampini           }
568b1b3d7a2SStefano Zampini         }
569b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
570b1b3d7a2SStefano Zampini       }
5715a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5725a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
5735a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
5745a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
5755a95e1ceSStefano Zampini       if (Sdense) {
5765a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
5775a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
578b1b3d7a2SStefano Zampini         }
5795a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
5805a95e1ceSStefano Zampini       } else {
5815a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
5825a95e1ceSStefano Zampini       }
5835a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
584a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
5855a95e1ceSStefano Zampini       local_size += subset_size;
5865a95e1ceSStefano Zampini     }
5875a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
588b1b3d7a2SStefano Zampini     /* free */
589b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
590b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
5915a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
592883469d8SStefano Zampini   } else {
593*be83ff47SStefano Zampini     Mat         A;
594883469d8SStefano Zampini     IS          is_A_all;
595*be83ff47SStefano Zampini     PetscScalar *work,*S_data;
596*be83ff47SStefano Zampini     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2;
597*be83ff47SStefano Zampini     PetscBool   restore_S;
598883469d8SStefano Zampini 
599883469d8SStefano Zampini     /* get working mat */
600a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
601d62866d3SStefano Zampini     if (!sub_schurs->is_dir) {
602883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
603d62866d3SStefano Zampini       size_schur = local_size;
604d62866d3SStefano Zampini     } else {
605d62866d3SStefano Zampini       IS list[2];
606d62866d3SStefano Zampini 
607d62866d3SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr);
608d62866d3SStefano Zampini       list[1] = sub_schurs->is_dir;
609d62866d3SStefano Zampini       ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr);
610d62866d3SStefano Zampini       ierr = ISDestroy(&list[0]);CHKERRQ(ierr);
611d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr);
612d62866d3SStefano Zampini       size_schur += local_size;
613d62866d3SStefano Zampini     }
614d62866d3SStefano Zampini     size_active_schur = local_size;
6156816873aSStefano Zampini     ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
616a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
617883469d8SStefano Zampini 
61808122e43SStefano Zampini     if (n_I) {
6199ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
620d62866d3SStefano Zampini         ierr = MatSetOption(A,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
621883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
622883469d8SStefano Zampini       } else {
623d62866d3SStefano Zampini         ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
624883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
625883469d8SStefano Zampini       }
626883469d8SStefano Zampini 
627883469d8SStefano Zampini       /* subsets ordered last */
628d62866d3SStefano Zampini       ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr);
629d62866d3SStefano Zampini       for (i=0;i<size_schur;i++) {
630883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
631883469d8SStefano Zampini       }
6325a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
633d62866d3SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr);
6345a95e1ceSStefano Zampini #endif
635883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
636883469d8SStefano Zampini 
637883469d8SStefano Zampini       /* factorization step */
6389ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
639883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
640*be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
641*be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
642*be83ff47SStefano Zampini #endif
643883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
644883469d8SStefano Zampini       } else {
645883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
646*be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
647*be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
648*be83ff47SStefano Zampini #endif
649883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
650883469d8SStefano Zampini       }
651883469d8SStefano Zampini 
652883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
6535a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
654883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
6555a95e1ceSStefano Zampini #endif
656d5574798SStefano Zampini 
657d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
658d5574798SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
659*be83ff47SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all));
660*be83ff47SStefano Zampini       restore_S = PETSC_TRUE;
661*be83ff47SStefano Zampini     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
662*be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
663*be83ff47SStefano Zampini       reuse_solvers = PETSC_FALSE;
664*be83ff47SStefano Zampini       restore_S = PETSC_FALSE;
665*be83ff47SStefano Zampini     }
666*be83ff47SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
667*be83ff47SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
668*be83ff47SStefano Zampini 
669*be83ff47SStefano Zampini     if (reuse_solvers) {
670d62866d3SStefano Zampini       Mat              A_II;
671d62866d3SStefano Zampini       PCBDDCReuseMumps msolv_ctx;
672d5574798SStefano Zampini 
673d62866d3SStefano Zampini       if (sub_schurs->reuse_mumps) {
6746816873aSStefano Zampini         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
6756816873aSStefano Zampini         ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
676d62866d3SStefano Zampini       }
677*be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
678d5574798SStefano Zampini       ierr = PetscNew(&msolv_ctx);CHKERRQ(ierr);
679*be83ff47SStefano Zampini       ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr);
680d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
681d5574798SStefano Zampini       msolv_ctx->F = F;
682d5574798SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
683d62866d3SStefano Zampini 
684d62866d3SStefano Zampini       /* interior solver */
685d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr);
686d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
687d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
688d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
689d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
690d62866d3SStefano Zampini 
691d62866d3SStefano Zampini       /* correction solver */
692d62866d3SStefano Zampini       /* auxiliary scatters are needed and are created in PCBDDCSetUpLocalScatters */
693d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr);
694d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
695d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
696d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
697d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
698*be83ff47SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->solB,&msolv_ctx->rhsB);CHKERRQ(ierr);
699d62866d3SStefano Zampini       sub_schurs->reuse_mumps = msolv_ctx;
700d5574798SStefano Zampini     }
70108122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
7025db18549SStefano Zampini 
703*be83ff47SStefano Zampini     /* Work arrays */
704*be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
70512d906b1SStefano Zampini     if (compute_Stilda) {
706a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
707d62866d3SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
708a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
709a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
71008122e43SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
711d62866d3SStefano Zampini       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
71208122e43SStefano Zampini       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
71308122e43SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
71412d906b1SStefano Zampini     }
715d2627357SStefano Zampini 
716*be83ff47SStefano Zampini     /* S_Ej_all */
717*be83ff47SStefano Zampini     cum = cum2 = 0;
718*be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
71965d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
72065d8bf0aSStefano Zampini       PetscInt j;
72165d8bf0aSStefano Zampini 
7225a95e1ceSStefano Zampini       /* get S_E */
723b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
724*be83ff47SStefano Zampini       if (restore_S && sub_schurs->is_hermitian) { /* When using MUMPS data I need to expand to upper triangular (column oriented) */
725*be83ff47SStefano Zampini         PetscInt k;
726*be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
727*be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
728*be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
729*be83ff47SStefano Zampini             work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
730*be83ff47SStefano Zampini           }
731*be83ff47SStefano Zampini         }
732*be83ff47SStefano Zampini       } else { /* copy to workspace */
733*be83ff47SStefano Zampini         PetscInt k;
734*be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
735*be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
736*be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
737*be83ff47SStefano Zampini           }
738*be83ff47SStefano Zampini         }
7399087bf02SStefano Zampini       }
7405a95e1ceSStefano Zampini       /* insert S_E values */
741*be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
7425a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
743a1337663SStefano Zampini 
744*be83ff47SStefano Zampini       /* if adaptivity is requested, invert S_E block */
745d2627357SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
74608122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
74708122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
748*be83ff47SStefano Zampini         if (!sub_schurs->is_hermitian) { /* TODO add sytrf/i for symmetric non hermitian */
7495a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
75008122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
7515a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
75208122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
7532972d61bSStefano Zampini         } else {
7545a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
7552972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7565a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
7572972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7582972d61bSStefano Zampini         }
75908122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
7605a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
7619087bf02SStefano Zampini       }
762*be83ff47SStefano Zampini       cum += subset_size;
763*be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
764*be83ff47SStefano Zampini     }
765*be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
766*be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
767*be83ff47SStefano Zampini     if (restore_S) {
768*be83ff47SStefano Zampini       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
769*be83ff47SStefano Zampini       ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr);
770*be83ff47SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
771*be83ff47SStefano Zampini     } else if (compute_Stilda) { /* we need to invert explicitly since we are not using MUMPS */
772*be83ff47SStefano Zampini       ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
773*be83ff47SStefano Zampini       ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
774*be83ff47SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
775*be83ff47SStefano Zampini       if (!sub_schurs->is_hermitian) { /* TODO add sytrf/i for symmetric non hermitian */
776*be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
777*be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
778*be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
779*be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
780*be83ff47SStefano Zampini       } else {
781*be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
782*be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
783*be83ff47SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
784*be83ff47SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
785*be83ff47SStefano Zampini       }
786*be83ff47SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
787*be83ff47SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
788*be83ff47SStefano Zampini     }
789*be83ff47SStefano Zampini #endif
790*be83ff47SStefano Zampini     /* S_Ej_tilda_all */
791*be83ff47SStefano Zampini     cum = cum2 = 0;
792*be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
793*be83ff47SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
794*be83ff47SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
795*be83ff47SStefano Zampini       /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
796*be83ff47SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
797*be83ff47SStefano Zampini         PetscInt j;
798*be83ff47SStefano Zampini         /* get (St^-1)_E */
799*be83ff47SStefano Zampini         if (sub_schurs->is_hermitian) { /* I need to expand to upper triangular (column oriented) */
800*be83ff47SStefano Zampini           PetscInt k;
801*be83ff47SStefano Zampini           for (k=0;k<subset_size;k++) {
802*be83ff47SStefano Zampini             for (j=k;j<subset_size;j++) {
803*be83ff47SStefano Zampini               work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
804*be83ff47SStefano Zampini               work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
805*be83ff47SStefano Zampini             }
806*be83ff47SStefano Zampini           }
807*be83ff47SStefano Zampini         } else { /* copy to workspace */
808*be83ff47SStefano Zampini           PetscInt k;
809*be83ff47SStefano Zampini           for (k=0;k<subset_size;k++) {
810*be83ff47SStefano Zampini             for (j=0;j<subset_size;j++) {
811*be83ff47SStefano Zampini               work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
812*be83ff47SStefano Zampini             }
813*be83ff47SStefano Zampini           }
814*be83ff47SStefano Zampini         }
815*be83ff47SStefano Zampini         for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
8165a95e1ceSStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8179087bf02SStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
81808122e43SStefano Zampini       }
819*be83ff47SStefano Zampini       cum += subset_size;
820*be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
821883469d8SStefano Zampini     }
822*be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
8235ec10c6aSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
824*be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
825*be83ff47SStefano Zampini     if (restore_S) {
826*be83ff47SStefano Zampini       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
8275db18549SStefano Zampini     }
828*be83ff47SStefano Zampini #endif
8299087bf02SStefano Zampini   }
830*be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
831a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
832a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
833a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
834a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
835a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
836a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
8375db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8385db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8395a95e1ceSStefano Zampini   if (compute_Stilda) {
840a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
841a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
84208122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
84308122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
84408122e43SStefano Zampini   }
845a1337663SStefano Zampini 
8465db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
8475db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
8483927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
8493927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
8505a95e1ceSStefano Zampini 
8515db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
8525a95e1ceSStefano Zampini   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
853d648f858SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
854d648f858SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
85508122e43SStefano Zampini 
856ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
857ac632422SStefano Zampini   if (faster_deluxe) {
8585a95e1ceSStefano Zampini     Mat         tmpmat;
8595a95e1ceSStefano Zampini     PetscScalar *array;
8605a95e1ceSStefano Zampini     PetscInt    cum;
8615a95e1ceSStefano Zampini 
8625a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
8635a95e1ceSStefano Zampini     cum = 0;
8645a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
8655a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
8665a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
8675a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8685a95e1ceSStefano Zampini       if (!sub_schurs->is_hermitian) {
8695a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
8705a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
8715a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
8725a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
8735a95e1ceSStefano Zampini       } else {
8745a95e1ceSStefano Zampini         PetscInt j,k;
8755a95e1ceSStefano Zampini 
8765a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
8775a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
8785a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
8795a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
8805a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
8815a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
8825a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
8835a95e1ceSStefano Zampini           }
8845a95e1ceSStefano Zampini         }
8855a95e1ceSStefano Zampini       }
8865a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
8875a95e1ceSStefano Zampini       cum += subset_size*subset_size;
8885a95e1ceSStefano Zampini     }
8895a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
8905a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
8915a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
892ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
8935a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
8945a95e1ceSStefano Zampini   }
8955a95e1ceSStefano Zampini 
896f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
8975a95e1ceSStefano Zampini   if (compute_Stilda) {
898a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
899a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
900d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
901d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
90208122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
90308122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
904d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
905d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
90608122e43SStefano Zampini   }
9073202ece2SStefano Zampini 
9085a95e1ceSStefano Zampini   /* free workspace */
90906a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
910a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
911a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
91208122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
9133202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
9145db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
9155a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
916b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
917b1b3d7a2SStefano Zampini }
918b1b3d7a2SStefano Zampini 
919b1b3d7a2SStefano Zampini #undef __FUNCT__
920b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
921a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
922b1b3d7a2SStefano Zampini {
9239bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
9245a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
925b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
926b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
927b1b3d7a2SStefano Zampini 
928b1b3d7a2SStefano Zampini   PetscFunctionBegin;
929b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
930b1b3d7a2SStefano Zampini   if (!is_sorted) {
931b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
932b1b3d7a2SStefano Zampini   }
933b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
934b1b3d7a2SStefano Zampini   if (!is_sorted) {
935b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
936b1b3d7a2SStefano Zampini   }
937b1b3d7a2SStefano Zampini 
938b1b3d7a2SStefano Zampini   /* reset any previous data */
939b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
940b1b3d7a2SStefano Zampini 
9415a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
9429bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
943b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
94408122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
94508122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
946b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
947b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
948b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
949b1b3d7a2SStefano Zampini   }
950b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
951b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
95208122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
953b1b3d7a2SStefano Zampini   }
954b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
955b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
956d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
957d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
958b1b3d7a2SStefano Zampini 
9596816873aSStefano Zampini   /* Determine if MUMPS can be used */
960883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
961883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
962a64f4aa4SStefano Zampini   sub_schurs->use_mumps = PETSC_TRUE;
963883469d8SStefano Zampini #endif
964b1b3d7a2SStefano Zampini 
965b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
966b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
967b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
968b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
9695db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
9705db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
9715db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
9725db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
9735a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
974b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
975a64f4aa4SStefano Zampini   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
976b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
977b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
978b96c3477SStefano Zampini     }
9799bb4a8caSStefano Zampini   }
980d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
981b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
982b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
98308122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
984b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
985b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
986b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
987b1b3d7a2SStefano Zampini }
988b1b3d7a2SStefano Zampini 
989b1b3d7a2SStefano Zampini #undef __FUNCT__
99034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
99134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
99234a97f8cSStefano Zampini {
99334a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
99434a97f8cSStefano Zampini   PetscErrorCode  ierr;
99534a97f8cSStefano Zampini 
99634a97f8cSStefano Zampini   PetscFunctionBegin;
99734a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
9985ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
99934a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
100034a97f8cSStefano Zampini   PetscFunctionReturn(0);
100134a97f8cSStefano Zampini }
100234a97f8cSStefano Zampini 
100334a97f8cSStefano Zampini #undef __FUNCT__
100434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
100534a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
100634a97f8cSStefano Zampini {
100734a97f8cSStefano Zampini   PetscErrorCode ierr;
100834a97f8cSStefano Zampini 
100934a97f8cSStefano Zampini   PetscFunctionBegin;
101034a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
101134a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
101234a97f8cSStefano Zampini   PetscFunctionReturn(0);
101334a97f8cSStefano Zampini }
101434a97f8cSStefano Zampini 
101534a97f8cSStefano Zampini #undef __FUNCT__
101634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
101734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
101834a97f8cSStefano Zampini {
101934a97f8cSStefano Zampini   PetscInt       i;
102034a97f8cSStefano Zampini   PetscErrorCode ierr;
102134a97f8cSStefano Zampini 
102234a97f8cSStefano Zampini   PetscFunctionBegin;
10231e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1024b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1025b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1026b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
10275db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
10285db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
102941c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
103041c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
103108122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1032a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
10335db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1034d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1035d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
103608122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
103708122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
103834a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1039b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
104034a97f8cSStefano Zampini   }
10415ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1042b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
10433dc780c3SStefano Zampini   }
1044d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps) {
1045d62866d3SStefano Zampini     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1046d62866d3SStefano Zampini   }
1047d62866d3SStefano Zampini   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
104834a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
104934a97f8cSStefano Zampini   PetscFunctionReturn(0);
105034a97f8cSStefano Zampini }
105134a97f8cSStefano Zampini 
105234a97f8cSStefano Zampini #undef __FUNCT__
105334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
10542a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
105534a97f8cSStefano Zampini {
105634a97f8cSStefano Zampini   PetscInt       i,j,n;
105734a97f8cSStefano Zampini   PetscErrorCode ierr;
105834a97f8cSStefano Zampini 
105934a97f8cSStefano Zampini   PetscFunctionBegin;
106034a97f8cSStefano Zampini   n = 0;
106134a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
106234a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
106334a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
106434a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
106534a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
106634a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
106734a97f8cSStefano Zampini         queue_tip[n] = dof;
106834a97f8cSStefano Zampini         n++;
106934a97f8cSStefano Zampini       }
107034a97f8cSStefano Zampini     }
107134a97f8cSStefano Zampini   }
107234a97f8cSStefano Zampini   *n_added = n;
107334a97f8cSStefano Zampini   PetscFunctionReturn(0);
107434a97f8cSStefano Zampini }
1075