xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 06a4e24a4ea34074cd1dc021a3bf2cd6abf1fe12)
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__
11e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve_Private"
12e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose)
13d62866d3SStefano Zampini {
14d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
15ffbec71fSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
16d62866d3SStefano Zampini   PetscInt         ival;
17ffbec71fSStefano Zampini #endif
18d62866d3SStefano Zampini   PetscErrorCode   ierr;
19d62866d3SStefano Zampini 
20d62866d3SStefano Zampini   PetscFunctionBegin;
21d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
22d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
23d62866d3SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
24d62866d3SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
25d62866d3SStefano Zampini #endif
26e28d306cSStefano Zampini   if (transpose) {
27e28d306cSStefano Zampini     ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
28e28d306cSStefano Zampini   } else {
296816873aSStefano Zampini     ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
30e28d306cSStefano Zampini   }
31d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
32d62866d3SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
33d62866d3SStefano Zampini #endif
34d62866d3SStefano Zampini   PetscFunctionReturn(0);
35d62866d3SStefano Zampini }
36d62866d3SStefano Zampini 
37d62866d3SStefano Zampini #undef __FUNCT__
38e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
39e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
40e28d306cSStefano Zampini {
41e28d306cSStefano Zampini   PetscErrorCode   ierr;
42e28d306cSStefano Zampini 
43e28d306cSStefano Zampini   PetscFunctionBegin;
44e28d306cSStefano Zampini   ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr);
45e28d306cSStefano Zampini   PetscFunctionReturn(0);
46e28d306cSStefano Zampini }
47e28d306cSStefano Zampini 
48e28d306cSStefano Zampini #undef __FUNCT__
49e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose"
50e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol)
51e28d306cSStefano Zampini {
52e28d306cSStefano Zampini   PetscErrorCode   ierr;
53e28d306cSStefano Zampini 
54e28d306cSStefano Zampini   PetscFunctionBegin;
55e28d306cSStefano Zampini   ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr);
56e28d306cSStefano Zampini   PetscFunctionReturn(0);
57e28d306cSStefano Zampini }
58e28d306cSStefano Zampini 
59e28d306cSStefano Zampini #undef __FUNCT__
60d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCReuseMumpsReset"
61d62866d3SStefano Zampini static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
62d62866d3SStefano Zampini {
63d62866d3SStefano Zampini   PetscErrorCode ierr;
64d62866d3SStefano Zampini 
65d62866d3SStefano Zampini   PetscFunctionBegin;
66d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
67d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
68d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
69d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
70d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
7153892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
7253892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
73d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
7453892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
7553892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
76d62866d3SStefano Zampini   PetscFunctionReturn(0);
77d62866d3SStefano Zampini }
78d5574798SStefano Zampini 
79d5574798SStefano Zampini #undef __FUNCT__
80e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve_Private"
81e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose)
82d5574798SStefano Zampini {
83d62866d3SStefano Zampini   PCBDDCReuseMumps ctx;
84d5574798SStefano Zampini   PetscScalar      *array,*array_mumps;
85ffbec71fSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
86d5574798SStefano Zampini   PetscInt         ival;
87ffbec71fSStefano Zampini #endif
88d5574798SStefano Zampini   PetscErrorCode   ierr;
89d5574798SStefano Zampini 
90d5574798SStefano Zampini   PetscFunctionBegin;
91d5574798SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
92d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
93d5574798SStefano Zampini   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
94d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
95d62866d3SStefano Zampini #endif
96d5574798SStefano Zampini   /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */
97d5574798SStefano Zampini   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
98d5574798SStefano Zampini   ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
99d62866d3SStefano Zampini   ierr = PetscMemcpy(array_mumps,array,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
100d5574798SStefano Zampini   ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
101d5574798SStefano Zampini   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
102d5574798SStefano Zampini 
103e28d306cSStefano Zampini   if (transpose) {
104e28d306cSStefano Zampini     ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
105e28d306cSStefano Zampini   } else {
106d5574798SStefano Zampini     ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
107e28d306cSStefano Zampini   }
108d5574798SStefano Zampini 
109d5574798SStefano Zampini   /* get back data to caller worskpace */
110d5574798SStefano Zampini   ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
111d5574798SStefano Zampini   ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
112d62866d3SStefano Zampini   ierr = PetscMemcpy(array,array_mumps,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
113d5574798SStefano Zampini   ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
114d5574798SStefano Zampini   ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
115d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
116d5574798SStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
117d62866d3SStefano Zampini #endif
118d5574798SStefano Zampini   PetscFunctionReturn(0);
119d5574798SStefano Zampini }
1203202ece2SStefano Zampini 
1213202ece2SStefano Zampini #undef __FUNCT__
122e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
123e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
124e28d306cSStefano Zampini {
125e28d306cSStefano Zampini   PetscErrorCode   ierr;
126e28d306cSStefano Zampini 
127e28d306cSStefano Zampini   PetscFunctionBegin;
128e28d306cSStefano Zampini   ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr);
129e28d306cSStefano Zampini   PetscFunctionReturn(0);
130e28d306cSStefano Zampini }
131e28d306cSStefano Zampini 
132e28d306cSStefano Zampini #undef __FUNCT__
133e28d306cSStefano Zampini #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose"
134e28d306cSStefano Zampini static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol)
135e28d306cSStefano Zampini {
136e28d306cSStefano Zampini   PetscErrorCode   ierr;
137e28d306cSStefano Zampini 
138e28d306cSStefano Zampini   PetscFunctionBegin;
139e28d306cSStefano Zampini   ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr);
140e28d306cSStefano Zampini   PetscFunctionReturn(0);
141e28d306cSStefano Zampini }
142e28d306cSStefano Zampini 
143e28d306cSStefano Zampini #undef __FUNCT__
1443202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
1455ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
1463202ece2SStefano Zampini {
1473202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
1483202ece2SStefano Zampini   KSP            ksp;
1493202ece2SStefano Zampini   PC             pc;
1503202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
1513202ece2SStefano Zampini   PetscReal      fill = 2.0;
152f11841e3SStefano Zampini   PetscInt       n_I;
1533202ece2SStefano Zampini   PetscMPIInt    size;
1543202ece2SStefano Zampini   PetscErrorCode ierr;
1553202ece2SStefano Zampini 
1563202ece2SStefano Zampini   PetscFunctionBegin;
1573202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
1583202ece2SStefano Zampini   if (size != 1) {
1593202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
1603202ece2SStefano Zampini   }
161f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
162f11841e3SStefano Zampini     PetscBool Sdense;
163f11841e3SStefano Zampini 
164f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
165f11841e3SStefano Zampini     if (!Sdense) {
166f11841e3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
167f11841e3SStefano Zampini     }
168f11841e3SStefano Zampini   }
1693202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
1703202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
1713202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1723202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
1733202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
1743202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
1753202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
1763202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
177f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
178f11841e3SStefano Zampini   if (n_I) {
1793202ece2SStefano Zampini     if (!Bdense) {
1803202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
1813202ece2SStefano Zampini     } else {
1823202ece2SStefano Zampini       Bd = B;
1833202ece2SStefano Zampini     }
1843202ece2SStefano Zampini 
1853202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1863202ece2SStefano Zampini       Mat fact;
1873202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1883202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1893202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1903202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1913202ece2SStefano Zampini     } else {
19207b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
19307b1e237SStefano Zampini 
19407b1e237SStefano Zampini       if (ex) {
1953202ece2SStefano Zampini         Mat Ainvd;
1963202ece2SStefano Zampini 
1973202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1983202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1993202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
20007b1e237SStefano Zampini       } else {
20107b1e237SStefano Zampini         Vec         sol,rhs;
20207b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
20307b1e237SStefano Zampini         PetscInt    i,nrhs,n;
20407b1e237SStefano Zampini 
20507b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
20607b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
20707b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
20807b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
20907b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
21007b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
21107b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
21207b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
21307b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
21407b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
21507b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
21607b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
21707b1e237SStefano Zampini         }
21807b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
21907b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
22007b1e237SStefano Zampini       }
2213202ece2SStefano Zampini     }
2225ec10c6aSStefano Zampini     if (!Bdense & !issym) {
2233202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2243202ece2SStefano Zampini     }
2255ec10c6aSStefano Zampini 
2265ec10c6aSStefano Zampini     if (!issym) {
2273202ece2SStefano Zampini       if (!Cdense) {
2283202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
2293202ece2SStefano Zampini       } else {
2303202ece2SStefano Zampini         Cd = C;
2313202ece2SStefano Zampini       }
2325ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2333202ece2SStefano Zampini       if (!Cdense) {
2343202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
2353202ece2SStefano Zampini       }
2365ec10c6aSStefano Zampini     } else {
2375ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2385ec10c6aSStefano Zampini       if (!Bdense) {
2395ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2405ec10c6aSStefano Zampini       }
2415ec10c6aSStefano Zampini     }
2425ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
243f11841e3SStefano Zampini   }
2443202ece2SStefano Zampini 
2453202ece2SStefano Zampini   if (D) {
2463202ece2SStefano Zampini     Mat       Dd;
2473202ece2SStefano Zampini     PetscBool Ddense;
2483202ece2SStefano Zampini 
2493202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
2503202ece2SStefano Zampini     if (!Ddense) {
2513202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
2523202ece2SStefano Zampini     } else {
2533202ece2SStefano Zampini       Dd = D;
2543202ece2SStefano Zampini     }
255f11841e3SStefano Zampini     if (n_I) {
2563202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
257f11841e3SStefano Zampini     } else {
258f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
259f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
260f11841e3SStefano Zampini       } else {
261f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
262f11841e3SStefano Zampini       }
263f11841e3SStefano Zampini     }
2643202ece2SStefano Zampini     if (!Ddense) {
2653202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
2663202ece2SStefano Zampini     }
2673202ece2SStefano Zampini   } else {
2683202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
2693202ece2SStefano Zampini   }
2703202ece2SStefano Zampini   PetscFunctionReturn(0);
2713202ece2SStefano Zampini }
27234a97f8cSStefano Zampini 
27334a97f8cSStefano Zampini #undef __FUNCT__
2741580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
275*06a4e24aSStefano 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 benign_trick)
276b1b3d7a2SStefano Zampini {
277be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
278be83ff47SStefano Zampini   Mat                    S_all;
2795a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
2805db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
281dc456d91SStefano Zampini   IS                     is_I,is_I_layer;
282dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
283dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
284dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
2855a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
286883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
28708122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
28806a4b1faSStefano Zampini   PetscScalar            *Bwork;
2895a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2905a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2915a95e1ceSStefano Zampini   MPI_Comm               comm_n;
292*06a4e24aSStefano Zampini   PetscBool              use_getr = PETSC_FALSE;
293b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
294b1b3d7a2SStefano Zampini 
295b1b3d7a2SStefano Zampini   PetscFunctionBegin;
296a64f4aa4SStefano Zampini   /* update info in sub_schurs */
297a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
298a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
2993301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
3003301b35fSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
301a64f4aa4SStefano Zampini   if (Ain) {
302a64f4aa4SStefano Zampini     PetscBool isseqaij;
3033301b35fSStefano Zampini     /* determine if we are dealing with hermitian positive definite problems */
3043301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
3053301b35fSStefano Zampini     if (Ain->symmetric_set) {
3063301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->symmetric;
3073301b35fSStefano Zampini     }
3083301b35fSStefano Zampini #else
3093301b35fSStefano Zampini     if (Ain->hermitian_set) {
3103301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->hermitian;
3113301b35fSStefano Zampini     }
3123301b35fSStefano Zampini #endif
3133301b35fSStefano Zampini     if (Ain->spd_set) {
3143301b35fSStefano Zampini       sub_schurs->is_posdef = Ain->spd;
3153301b35fSStefano Zampini     }
316a64f4aa4SStefano Zampini 
3173301b35fSStefano Zampini     /* check */
318a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
319*06a4e24aSStefano Zampini     if (compute_Stilda && (!Ain->hermitian_set || !Ain->spd_set)) {
3203301b35fSStefano Zampini       PetscInt lsize;
3213301b35fSStefano Zampini 
3223301b35fSStefano Zampini       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
3233301b35fSStefano Zampini       if (lsize) {
3243301b35fSStefano Zampini         PetscScalar val;
3253301b35fSStefano Zampini         PetscReal   norm;
3263301b35fSStefano Zampini         Vec         vec1,vec2,vec3;
3273301b35fSStefano Zampini 
3283301b35fSStefano Zampini         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
3293301b35fSStefano Zampini         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
3303301b35fSStefano Zampini         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
3313301b35fSStefano Zampini         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
3323301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
3333301b35fSStefano Zampini         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3343301b35fSStefano Zampini #else
3353301b35fSStefano Zampini         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3363301b35fSStefano Zampini #endif
3373301b35fSStefano Zampini         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
3383301b35fSStefano Zampini         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
339862806e4SStefano Zampini         if (norm > PetscSqrtReal(PETSC_SMALL)) {
3403301b35fSStefano Zampini           sub_schurs->is_hermitian = PETSC_FALSE;
3413301b35fSStefano Zampini         } else {
3423301b35fSStefano Zampini           sub_schurs->is_hermitian = PETSC_TRUE;
3433301b35fSStefano Zampini         }
344*06a4e24aSStefano Zampini         if (!Ain->spd_set) {
3453301b35fSStefano Zampini           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
3463301b35fSStefano Zampini           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
347*06a4e24aSStefano Zampini         }
3483301b35fSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
3493301b35fSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
3503301b35fSStefano Zampini         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
3513301b35fSStefano Zampini       } else {
3523301b35fSStefano Zampini         sub_schurs->is_hermitian = PETSC_TRUE;
3533301b35fSStefano Zampini         sub_schurs->is_posdef = PETSC_TRUE;
3543301b35fSStefano Zampini       }
3553301b35fSStefano Zampini     }
356a64f4aa4SStefano Zampini     if (isseqaij) {
357a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
358a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
3593301b35fSStefano Zampini     } else {
360a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
361a64f4aa4SStefano Zampini     }
362a64f4aa4SStefano Zampini   }
3633301b35fSStefano Zampini 
364a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
365a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
366a64f4aa4SStefano Zampini   if (sub_schurs->use_mumps) {
367a64f4aa4SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
368a64f4aa4SStefano Zampini   }
369a64f4aa4SStefano Zampini 
3705a95e1ceSStefano Zampini   /* preliminary checks */
3715a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps && compute_Stilda) {
3725a95e1ceSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
3735a95e1ceSStefano Zampini   }
3745a95e1ceSStefano Zampini 
3755a95e1ceSStefano Zampini   /* restrict work on active processes */
3765a95e1ceSStefano Zampini   color = 0;
3775a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
3785a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
3795a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
3805a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
3815a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
3825a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
3835a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
3845a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3855a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
3865a95e1ceSStefano Zampini     PetscFunctionReturn(0);
3875a95e1ceSStefano Zampini   }
38881ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
3895a95e1ceSStefano Zampini 
390b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
391883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
392a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
3933301b35fSStefano Zampini     PetscBool isseqsbaij;
394a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
3953301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
3963301b35fSStefano Zampini     if (isseqsbaij) {
397a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
398a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
399a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
400a64f4aa4SStefano Zampini     } else {
401a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
402a64f4aa4SStefano Zampini       A_BB = tA_BB;
403a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
404a64f4aa4SStefano Zampini       A_IB = tA_IB;
405a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
406a64f4aa4SStefano Zampini       A_BI = tA_BI;
407f11841e3SStefano Zampini     }
408a58a30b4SStefano Zampini   } else {
4095a95e1ceSStefano Zampini     A_II = NULL;
4105a95e1ceSStefano Zampini     A_IB = NULL;
4115a95e1ceSStefano Zampini     A_BI = NULL;
4125a95e1ceSStefano Zampini     A_BB = NULL;
413b1b3d7a2SStefano Zampini   }
4145a95e1ceSStefano Zampini   S_all = NULL;
415b1b3d7a2SStefano Zampini 
416b1b3d7a2SStefano Zampini   /* determine interior problems */
4173dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4183dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
419b1b3d7a2SStefano Zampini     PetscBT                touched;
420b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
421b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
422b1b3d7a2SStefano Zampini 
4233dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
4243dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
4253dc780c3SStefano Zampini     }
426b1b3d7a2SStefano Zampini     /* get sizes */
427b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
428b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
429b1b3d7a2SStefano Zampini 
430b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
431b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
432b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
433b1b3d7a2SStefano Zampini 
434b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
435b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
436b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
437b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
438b1b3d7a2SStefano Zampini     }
439b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
440b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
441b1b3d7a2SStefano Zampini 
442b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
443b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
444b1b3d7a2SStefano Zampini     n_prev_added = n_B;
445b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
446b1b3d7a2SStefano Zampini       PetscInt n_added;
447b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
448b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
449b1b3d7a2SStefano 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);
450b1b3d7a2SStefano Zampini       }
451b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
452b1b3d7a2SStefano Zampini       n_prev_added = n_added;
453b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
454b1b3d7a2SStefano Zampini       if (!n_added) break;
455b1b3d7a2SStefano Zampini     }
456b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
457b1b3d7a2SStefano Zampini 
458883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
459a9b99552SStefano 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);
460b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
461a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
462883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
463883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
464b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
465b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
466a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
467b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
468b1b3d7a2SStefano Zampini 
469b1b3d7a2SStefano Zampini       /* II block */
470b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
471b1b3d7a2SStefano Zampini     }
472b1b3d7a2SStefano Zampini   } else {
473b1b3d7a2SStefano Zampini     PetscInt n_I;
474b1b3d7a2SStefano Zampini 
475b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
476b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
477a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
478b1b3d7a2SStefano Zampini 
479b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
480883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
481b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
482b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
483b1b3d7a2SStefano Zampini 
484b1b3d7a2SStefano Zampini       /* II block is the same */
485b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
486b1b3d7a2SStefano Zampini       AE_II = A_II;
487b1b3d7a2SStefano Zampini     }
488b1b3d7a2SStefano Zampini   }
4895a95e1ceSStefano Zampini 
490883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
4915a95e1ceSStefano Zampini   max_subset_size = 0;
492883469d8SStefano Zampini   local_size = 0;
4935a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4945a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4955a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
496883469d8SStefano Zampini     local_size += subset_size;
497883469d8SStefano Zampini   }
498883469d8SStefano Zampini 
499883469d8SStefano Zampini   /* Work arrays for local indices */
500883469d8SStefano Zampini   extra = 0;
50181ea8064SStefano Zampini   if (sub_schurs->use_mumps && is_I_layer) {
502a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
503883469d8SStefano Zampini   }
504883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
505883469d8SStefano Zampini   if (extra) {
506883469d8SStefano Zampini     const PetscInt *idxs;
507a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
508883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
509a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
510883469d8SStefano Zampini   }
511883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
512dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
513dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
514883469d8SStefano Zampini 
515883469d8SStefano Zampini   /* Get local indices in local numbering */
516883469d8SStefano Zampini   local_size = 0;
5175a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
518883469d8SStefano Zampini     PetscInt j;
519883469d8SStefano Zampini     const    PetscInt *idxs;
520883469d8SStefano Zampini 
5215a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5225a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
523eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
524eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
525eb595f79SStefano Zampini     auxnum2[i] = subset_size;
526883469d8SStefano Zampini     /* subset indices in local numbering */
527883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5285a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
529883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
530883469d8SStefano Zampini     local_size += subset_size;
531883469d8SStefano Zampini   }
532883469d8SStefano Zampini 
5335a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
534d2627357SStefano Zampini   Bwork = NULL;
535d2627357SStefano Zampini   pivots = NULL;
536*06a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
537d2627357SStefano Zampini     PetscScalar lwork;
538d2627357SStefano Zampini 
539*06a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
540d2627357SStefano Zampini     B_lwork = -1;
541d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
542d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
543d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
544d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
545d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
546d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
547d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
548d2627357SStefano Zampini   }
549d2627357SStefano Zampini 
550d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
551dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
552dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
553dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
554dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
555dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
556dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
557dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
558dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
559dc456d91SStefano Zampini   if (i != local_size) {
560dc456d91SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size);
561eb595f79SStefano Zampini   }
562dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
563e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
564d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
565d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
566d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
567d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
5682972d61bSStefano Zampini 
5695a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
5705a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
5715a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
5725a95e1ceSStefano Zampini 
5735a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
5745a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
5755a95e1ceSStefano Zampini     if (subset_size != local_size) {
5765a95e1ceSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
5775a95e1ceSStefano Zampini     }
5785a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
579b1b3d7a2SStefano Zampini   }
580b1b3d7a2SStefano Zampini 
5815a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
5825a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
5835a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
5845a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
5855a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
5865a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
587aa83b6aeSStefano Zampini   }
588b1b3d7a2SStefano Zampini 
5895a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
590be83ff47SStefano Zampini   F = NULL;
5915a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps) {
5925a95e1ceSStefano Zampini     Mat         S_Ej_expl;
5935a95e1ceSStefano Zampini     PetscScalar *work;
5945a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
5955a95e1ceSStefano Zampini     PetscBool   Sdense;
5965a95e1ceSStefano Zampini 
5975a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
5985a95e1ceSStefano Zampini     local_size = 0;
599b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
6005a95e1ceSStefano Zampini       IS  is_subset_B;
6015a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
6025a95e1ceSStefano Zampini 
6035a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
6045a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
6055a95e1ceSStefano Zampini       /* EE block */
6065a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
6075a95e1ceSStefano Zampini       /* IE block */
6085a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
6095a95e1ceSStefano Zampini       /* EI block */
6105a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
6115a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
6125a95e1ceSStefano Zampini       } else {
6135a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
6145a95e1ceSStefano Zampini       }
615a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
6165a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
6175a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
6185a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
6195a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
620b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
621b1b3d7a2SStefano Zampini         KSP ksp;
622b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
6235a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
624b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
625b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
626b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
627b1b3d7a2SStefano Zampini         KSPType   ksp_type;
628b1b3d7a2SStefano Zampini         PetscInt  n_internal;
6295a95e1ceSStefano Zampini         PetscBool ispcnone;
630b1b3d7a2SStefano Zampini 
631b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
6325a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
633b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
634b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
635b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
636b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
6375a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
6385a95e1ceSStefano Zampini         if (!ispcnone) {
6395a95e1ceSStefano Zampini           PCType pc_type;
640b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
641b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
6425a95e1ceSStefano Zampini         } else {
6435a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
6445a95e1ceSStefano Zampini         }
645b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
646b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
647b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
648b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
649b1b3d7a2SStefano Zampini           if (solver) {
650b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
651b1b3d7a2SStefano Zampini           }
652b1b3d7a2SStefano Zampini         }
653b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
654b1b3d7a2SStefano Zampini       }
6555a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6565a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
6575a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
6585a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
6595a95e1ceSStefano Zampini       if (Sdense) {
6605a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
6615a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
662b1b3d7a2SStefano Zampini         }
6635a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
6645a95e1ceSStefano Zampini       } else {
6655a95e1ceSStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
6665a95e1ceSStefano Zampini       }
6675a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
668a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
6695a95e1ceSStefano Zampini       local_size += subset_size;
6705a95e1ceSStefano Zampini     }
6715a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
672b1b3d7a2SStefano Zampini     /* free */
673b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
674b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
6755a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
676883469d8SStefano Zampini   } else {
677be83ff47SStefano Zampini     Mat         A;
678883469d8SStefano Zampini     IS          is_A_all;
679be83ff47SStefano Zampini     PetscScalar *work,*S_data;
680be83ff47SStefano Zampini     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2;
681*06a4e24aSStefano Zampini     PetscBool   mumps_S,S_upper_triangular = PETSC_FALSE;
682883469d8SStefano Zampini 
683883469d8SStefano Zampini     /* get working mat */
68481ea8064SStefano Zampini     n_I = 0;
68581ea8064SStefano Zampini     if (is_I_layer) {
686a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
68781ea8064SStefano Zampini     }
688d62866d3SStefano Zampini     if (!sub_schurs->is_dir) {
689883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
690d62866d3SStefano Zampini       size_schur = local_size;
691d62866d3SStefano Zampini     } else {
692d62866d3SStefano Zampini       IS list[2];
693d62866d3SStefano Zampini 
694d62866d3SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr);
695d62866d3SStefano Zampini       list[1] = sub_schurs->is_dir;
696d62866d3SStefano Zampini       ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr);
697d62866d3SStefano Zampini       ierr = ISDestroy(&list[0]);CHKERRQ(ierr);
698d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr);
699d62866d3SStefano Zampini       size_schur += local_size;
700d62866d3SStefano Zampini     }
70153892102SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
70253892102SStefano Zampini     size_active_schur = local_size; /* size active schurs does not count any dirichlet dof on the interface */
7036816873aSStefano Zampini     ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
704a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
7053301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
7063301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
7073301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
708883469d8SStefano Zampini 
70908122e43SStefano Zampini     if (n_I) {
7109ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
711883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
712883469d8SStefano Zampini       } else {
713883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
714883469d8SStefano Zampini       }
715883469d8SStefano Zampini       /* subsets ordered last */
716d62866d3SStefano Zampini       ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr);
717d62866d3SStefano Zampini       for (i=0;i<size_schur;i++) {
718883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
719883469d8SStefano Zampini       }
7205a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
721d62866d3SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr);
7225a95e1ceSStefano Zampini #endif
723883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
724883469d8SStefano Zampini 
725883469d8SStefano Zampini       /* factorization step */
7269ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
727883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
728be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
729be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
730be83ff47SStefano Zampini #endif
731883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
732*06a4e24aSStefano Zampini         S_upper_triangular = PETSC_TRUE;
733883469d8SStefano Zampini       } else {
734883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
735be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
736be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
737be83ff47SStefano Zampini #endif
738883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
739883469d8SStefano Zampini       }
740883469d8SStefano Zampini 
741883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
7425a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
743883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
7445a95e1ceSStefano Zampini #endif
745d5574798SStefano Zampini 
746d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
747d5574798SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
748be83ff47SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all));
7494a6c6b0dSStefano Zampini       mumps_S = PETSC_TRUE;
750be83ff47SStefano Zampini     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
751be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
752be83ff47SStefano Zampini       reuse_solvers = PETSC_FALSE;
7534a6c6b0dSStefano Zampini       mumps_S = PETSC_FALSE;
754be83ff47SStefano Zampini     }
755be83ff47SStefano Zampini 
756*06a4e24aSStefano Zampini     /* when using te benign subspace trick, the local Schur complements are SPD */
757*06a4e24aSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
758*06a4e24aSStefano Zampini 
759be83ff47SStefano Zampini     if (reuse_solvers) {
760d62866d3SStefano Zampini       Mat              A_II;
76153892102SStefano Zampini       Vec              vec1_B;
762d62866d3SStefano Zampini       PCBDDCReuseMumps msolv_ctx;
763d5574798SStefano Zampini 
764d62866d3SStefano Zampini       if (sub_schurs->reuse_mumps) {
7656816873aSStefano Zampini         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
766e28d306cSStefano Zampini       } else {
767e28d306cSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
768d62866d3SStefano Zampini       }
769e28d306cSStefano Zampini       msolv_ctx = sub_schurs->reuse_mumps;
770be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
771be83ff47SStefano Zampini       ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr);
772d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
773d5574798SStefano Zampini       msolv_ctx->F = F;
774d5574798SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
775d62866d3SStefano Zampini 
776d62866d3SStefano Zampini       /* interior solver */
777d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr);
778d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
779d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
780d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
781d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
782e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
783d62866d3SStefano Zampini 
784d62866d3SStefano Zampini       /* correction solver */
785d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr);
786d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
787d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
788d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
789d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
790e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
79153892102SStefano Zampini 
79253892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
79353892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
79453892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
79553892102SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
79653892102SStefano Zampini       ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
79753892102SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
79853892102SStefano Zampini       ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
79953892102SStefano Zampini       msolv_ctx->is_R = is_A_all;
800d5574798SStefano Zampini     }
80108122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
80253892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
8035db18549SStefano Zampini 
804be83ff47SStefano Zampini     /* Work arrays */
805be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
8065a94c880SStefano Zampini 
8075a94c880SStefano Zampini     /* matrices for adaptive selection */
80812d906b1SStefano Zampini     if (compute_Stilda) {
8095a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
8105a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
8115a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
8125a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
8135a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
8145a94c880SStefano Zampini       }
8155a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all) {
8165a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
8175a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
8185a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
8195a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
8205a94c880SStefano Zampini       }
82112d906b1SStefano Zampini     }
822d2627357SStefano Zampini 
823be83ff47SStefano Zampini     /* S_Ej_all */
824be83ff47SStefano Zampini     cum = cum2 = 0;
825be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
82665d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
82765d8bf0aSStefano Zampini       PetscInt j;
82865d8bf0aSStefano Zampini 
8295a95e1ceSStefano Zampini       /* get S_E */
830b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
831*06a4e24aSStefano Zampini       if (mumps_S && S_upper_triangular) { /* I need to expand the upper triangular data (column oriented) */
832be83ff47SStefano Zampini         PetscInt k;
833be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
834be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
835be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
836be83ff47SStefano Zampini             work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
837be83ff47SStefano Zampini           }
838be83ff47SStefano Zampini         }
839*06a4e24aSStefano Zampini       } else { /* just copy to workspace */
840be83ff47SStefano Zampini         PetscInt k;
841be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
842be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
843be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
844be83ff47SStefano Zampini           }
845be83ff47SStefano Zampini         }
8469087bf02SStefano Zampini       }
8475a95e1ceSStefano Zampini       /* insert S_E values */
848be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
8495a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
850a1337663SStefano Zampini 
851be83ff47SStefano Zampini       /* if adaptivity is requested, invert S_E block */
852862806e4SStefano Zampini       if (compute_Stilda) {
85308122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
85408122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
855*06a4e24aSStefano Zampini         if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
8565a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
8572972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
8585a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
8592972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
860d6462365SStefano Zampini         } else {
861d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
862d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
863d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
864d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
8652972d61bSStefano Zampini         }
86608122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
8675a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8689087bf02SStefano Zampini       }
869be83ff47SStefano Zampini       cum += subset_size;
870be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
871be83ff47SStefano Zampini     }
872be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
8734a6c6b0dSStefano Zampini 
874be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
8754a6c6b0dSStefano Zampini     if (mumps_S) {
876be83ff47SStefano Zampini       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
8774a6c6b0dSStefano Zampini     }
8784a6c6b0dSStefano Zampini #endif
8794a6c6b0dSStefano Zampini 
88045951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
8814a6c6b0dSStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
8824a6c6b0dSStefano Zampini         PetscInt j;
8834a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
8845a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8854a6c6b0dSStefano Zampini       } else {
8864a6c6b0dSStefano Zampini         if (mumps_S) { /* use MatMumps calls to invert S */
8874a6c6b0dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
888be83ff47SStefano Zampini           ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr);
889be83ff47SStefano Zampini           ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
8904a6c6b0dSStefano Zampini #endif
8914a6c6b0dSStefano Zampini         } else { /* we need to invert explicitly since we are not using MUMPS for S */
892be83ff47SStefano Zampini           ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
893be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
894be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
895*06a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
896be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
897be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
898be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
899be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
900d6462365SStefano Zampini           } else {
901d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
902d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
903d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
904d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
905be83ff47SStefano Zampini           }
906be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
907be83ff47SStefano Zampini           ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
908be83ff47SStefano Zampini         }
909be83ff47SStefano Zampini         /* S_Ej_tilda_all */
910be83ff47SStefano Zampini         cum = cum2 = 0;
911be83ff47SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
912be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
913be83ff47SStefano Zampini           PetscInt j;
914862806e4SStefano Zampini 
915862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
916be83ff47SStefano Zampini           /* get (St^-1)_E */
917*06a4e24aSStefano Zampini           /* Here I don't need to expand to upper triangular since St^-1
918*06a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
919*06a4e24aSStefano Zampini           if (S_upper_triangular) {
920be83ff47SStefano Zampini             PetscInt k;
921be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
922be83ff47SStefano Zampini               for (j=k;j<subset_size;j++) {
923be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
924be83ff47SStefano Zampini               }
925be83ff47SStefano Zampini             }
926be83ff47SStefano Zampini           } else { /* copy to workspace */
927be83ff47SStefano Zampini             PetscInt k;
928be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
929be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
930be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
931be83ff47SStefano Zampini               }
932be83ff47SStefano Zampini             }
933be83ff47SStefano Zampini           }
934be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
9355a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
936be83ff47SStefano Zampini           cum += subset_size;
937be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
938883469d8SStefano Zampini         }
939be83ff47SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
940be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
9414a6c6b0dSStefano Zampini         if (mumps_S) {
942be83ff47SStefano Zampini           ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
9435db18549SStefano Zampini         }
944be83ff47SStefano Zampini #endif
9459087bf02SStefano Zampini       }
9464a6c6b0dSStefano Zampini     }
9474a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
9484a6c6b0dSStefano Zampini   }
9495a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
950be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
951a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
952a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
953a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
954a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
955a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
9565db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9575db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9585a95e1ceSStefano Zampini   if (compute_Stilda) {
9595a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9605a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9615a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9625a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
96308122e43SStefano Zampini   }
964a1337663SStefano Zampini 
9655db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
9665db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
9673927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
9683927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
9695a95e1ceSStefano Zampini 
9705db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
9715a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
972dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
9735a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
9745a94c880SStefano Zampini   } else {
9755a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
9765a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
977dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
9785a94c880SStefano Zampini   }
97908122e43SStefano Zampini 
980ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
981ac632422SStefano Zampini   if (faster_deluxe) {
9825a95e1ceSStefano Zampini     Mat         tmpmat;
9835a95e1ceSStefano Zampini     PetscScalar *array;
9845a95e1ceSStefano Zampini     PetscInt    cum;
9855a95e1ceSStefano Zampini 
9865a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
9875a95e1ceSStefano Zampini     cum = 0;
9885a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
9895a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
9905a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
9915a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
992*06a4e24aSStefano Zampini       if (!use_getr) {
9935a95e1ceSStefano Zampini         PetscInt j,k;
9945a95e1ceSStefano Zampini 
9955a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
9965a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
9975a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
9985a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
9995a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
10005a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
10015a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
10025a95e1ceSStefano Zampini           }
10035a95e1ceSStefano Zampini         }
1004d6462365SStefano Zampini       } else {
1005d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1006d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1007d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1008d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
10095a95e1ceSStefano Zampini       }
10105a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
10115a95e1ceSStefano Zampini       cum += subset_size*subset_size;
10125a95e1ceSStefano Zampini     }
10135a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
10145a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
10155a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1016ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
10175a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
10185a95e1ceSStefano Zampini   }
10195a95e1ceSStefano Zampini 
1020f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
10215a95e1ceSStefano Zampini   if (compute_Stilda) {
10225a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1023a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
10245a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1025dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
10265a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
102708122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
10285a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_inv_all;
1029dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
103008122e43SStefano Zampini   }
10313202ece2SStefano Zampini 
10325a95e1ceSStefano Zampini   /* free workspace */
10335a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
103406a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1035a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
10363202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1037dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
10385a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1039b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1040b1b3d7a2SStefano Zampini }
1041b1b3d7a2SStefano Zampini 
1042b1b3d7a2SStefano Zampini #undef __FUNCT__
1043b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
1044a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1045b1b3d7a2SStefano Zampini {
10469bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
10475a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1048b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1049b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1050b1b3d7a2SStefano Zampini 
1051b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1052b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1053b1b3d7a2SStefano Zampini   if (!is_sorted) {
1054b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1055b1b3d7a2SStefano Zampini   }
1056b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1057b1b3d7a2SStefano Zampini   if (!is_sorted) {
1058b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1059b1b3d7a2SStefano Zampini   }
1060b1b3d7a2SStefano Zampini 
1061b1b3d7a2SStefano Zampini   /* reset any previous data */
1062b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1063b1b3d7a2SStefano Zampini 
10645a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
10659bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1066b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
106708122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1068b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1069b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
1070b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
1071b1b3d7a2SStefano Zampini   }
1072b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
1073b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
107408122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1075b1b3d7a2SStefano Zampini   }
1076b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
1077b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
1078d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1079d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1080b1b3d7a2SStefano Zampini 
10816816873aSStefano Zampini   /* Determine if MUMPS can be used */
1082883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
1083883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1084a64f4aa4SStefano Zampini   sub_schurs->use_mumps = PETSC_TRUE;
1085883469d8SStefano Zampini #endif
1086b1b3d7a2SStefano Zampini 
1087b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1088b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1089b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1090b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
10915db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
10925db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
10935db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
10945db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
10955a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1096b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1097a64f4aa4SStefano Zampini   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1098b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1099b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1100b96c3477SStefano Zampini     }
11019bb4a8caSStefano Zampini   }
1102d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
1103b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1104b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
110508122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1106b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1107b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1108b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1109b1b3d7a2SStefano Zampini }
1110b1b3d7a2SStefano Zampini 
1111b1b3d7a2SStefano Zampini #undef __FUNCT__
111234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
111334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
111434a97f8cSStefano Zampini {
111534a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
111634a97f8cSStefano Zampini   PetscErrorCode  ierr;
111734a97f8cSStefano Zampini 
111834a97f8cSStefano Zampini   PetscFunctionBegin;
111934a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
11205ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
112134a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
112234a97f8cSStefano Zampini   PetscFunctionReturn(0);
112334a97f8cSStefano Zampini }
112434a97f8cSStefano Zampini 
112534a97f8cSStefano Zampini #undef __FUNCT__
112634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
112734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
112834a97f8cSStefano Zampini {
112934a97f8cSStefano Zampini   PetscErrorCode ierr;
113034a97f8cSStefano Zampini 
113134a97f8cSStefano Zampini   PetscFunctionBegin;
113234a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
113334a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
113434a97f8cSStefano Zampini   PetscFunctionReturn(0);
113534a97f8cSStefano Zampini }
113634a97f8cSStefano Zampini 
113734a97f8cSStefano Zampini #undef __FUNCT__
113834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
113934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
114034a97f8cSStefano Zampini {
114134a97f8cSStefano Zampini   PetscInt       i;
114234a97f8cSStefano Zampini   PetscErrorCode ierr;
114334a97f8cSStefano Zampini 
114434a97f8cSStefano Zampini   PetscFunctionBegin;
11451e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1146b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1147b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1148b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
11495db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
11505db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
115141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
115241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
115308122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1154a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11555db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1156d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1157d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
115808122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
115934a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1160b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
116134a97f8cSStefano Zampini   }
11625ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1163b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
11643dc780c3SStefano Zampini   }
1165d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps) {
1166d62866d3SStefano Zampini     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1167d62866d3SStefano Zampini   }
1168d62866d3SStefano Zampini   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
116934a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
117034a97f8cSStefano Zampini   PetscFunctionReturn(0);
117134a97f8cSStefano Zampini }
117234a97f8cSStefano Zampini 
117334a97f8cSStefano Zampini #undef __FUNCT__
117434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
11752a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
117634a97f8cSStefano Zampini {
117734a97f8cSStefano Zampini   PetscInt       i,j,n;
117834a97f8cSStefano Zampini   PetscErrorCode ierr;
117934a97f8cSStefano Zampini 
118034a97f8cSStefano Zampini   PetscFunctionBegin;
118134a97f8cSStefano Zampini   n = 0;
118234a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
118334a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
118434a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
118534a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
118634a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
118734a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
118834a97f8cSStefano Zampini         queue_tip[n] = dof;
118934a97f8cSStefano Zampini         n++;
119034a97f8cSStefano Zampini       }
119134a97f8cSStefano Zampini     }
119234a97f8cSStefano Zampini   }
119334a97f8cSStefano Zampini   *n_added = n;
119434a97f8cSStefano Zampini   PetscFunctionReturn(0);
119534a97f8cSStefano Zampini }
1196