xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 6c4ed00291fe44f94936dd2f04c02ab3c442e77c)
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);
158*6c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
159f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
160f11841e3SStefano Zampini     PetscBool Sdense;
161f11841e3SStefano Zampini 
162f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
163*6c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
164f11841e3SStefano Zampini   }
1653202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
1663202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
1673202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1683202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
1693202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
1703202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
1713202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
1723202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
173f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
174f11841e3SStefano Zampini   if (n_I) {
1753202ece2SStefano Zampini     if (!Bdense) {
1763202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
1773202ece2SStefano Zampini     } else {
1783202ece2SStefano Zampini       Bd = B;
1793202ece2SStefano Zampini     }
1803202ece2SStefano Zampini 
1813202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
1823202ece2SStefano Zampini       Mat fact;
1833202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
1843202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
1853202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
1863202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
1873202ece2SStefano Zampini     } else {
18807b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
18907b1e237SStefano Zampini 
19007b1e237SStefano Zampini       if (ex) {
1913202ece2SStefano Zampini         Mat Ainvd;
1923202ece2SStefano Zampini 
1933202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
1943202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
1953202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
19607b1e237SStefano Zampini       } else {
19707b1e237SStefano Zampini         Vec         sol,rhs;
19807b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
19907b1e237SStefano Zampini         PetscInt    i,nrhs,n;
20007b1e237SStefano Zampini 
20107b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
20207b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
20307b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
20407b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
20507b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
20607b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
20707b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
20807b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
20907b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
21007b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
21107b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
21207b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
21307b1e237SStefano Zampini         }
21407b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
21507b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
21607b1e237SStefano Zampini       }
2173202ece2SStefano Zampini     }
2185ec10c6aSStefano Zampini     if (!Bdense & !issym) {
2193202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2203202ece2SStefano Zampini     }
2215ec10c6aSStefano Zampini 
2225ec10c6aSStefano Zampini     if (!issym) {
2233202ece2SStefano Zampini       if (!Cdense) {
2243202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
2253202ece2SStefano Zampini       } else {
2263202ece2SStefano Zampini         Cd = C;
2273202ece2SStefano Zampini       }
2285ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2293202ece2SStefano Zampini       if (!Cdense) {
2303202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
2313202ece2SStefano Zampini       }
2325ec10c6aSStefano Zampini     } else {
2335ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
2345ec10c6aSStefano Zampini       if (!Bdense) {
2355ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
2365ec10c6aSStefano Zampini       }
2375ec10c6aSStefano Zampini     }
2385ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
239f11841e3SStefano Zampini   }
2403202ece2SStefano Zampini 
2413202ece2SStefano Zampini   if (D) {
2423202ece2SStefano Zampini     Mat       Dd;
2433202ece2SStefano Zampini     PetscBool Ddense;
2443202ece2SStefano Zampini 
2453202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
2463202ece2SStefano Zampini     if (!Ddense) {
2473202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
2483202ece2SStefano Zampini     } else {
2493202ece2SStefano Zampini       Dd = D;
2503202ece2SStefano Zampini     }
251f11841e3SStefano Zampini     if (n_I) {
2523202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
253f11841e3SStefano Zampini     } else {
254f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
255f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
256f11841e3SStefano Zampini       } else {
257f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
258f11841e3SStefano Zampini       }
259f11841e3SStefano Zampini     }
2603202ece2SStefano Zampini     if (!Ddense) {
2613202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
2623202ece2SStefano Zampini     }
2633202ece2SStefano Zampini   } else {
2643202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
2653202ece2SStefano Zampini   }
2663202ece2SStefano Zampini   PetscFunctionReturn(0);
2673202ece2SStefano Zampini }
26834a97f8cSStefano Zampini 
26934a97f8cSStefano Zampini #undef __FUNCT__
2701580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
271862806e4SStefano 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)
272b1b3d7a2SStefano Zampini {
273be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
274be83ff47SStefano Zampini   Mat                    S_all;
2755a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
2765db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
277dc456d91SStefano Zampini   IS                     is_I,is_I_layer;
278dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
279dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
280dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
2815a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
282883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
28308122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
28406a4b1faSStefano Zampini   PetscScalar            *Bwork;
2855a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
2865a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
2875a95e1ceSStefano Zampini   MPI_Comm               comm_n;
288b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
289b1b3d7a2SStefano Zampini 
290b1b3d7a2SStefano Zampini   PetscFunctionBegin;
291a64f4aa4SStefano Zampini   /* update info in sub_schurs */
292a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
293a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
2943301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2953301b35fSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
296a64f4aa4SStefano Zampini   if (Ain) {
297a64f4aa4SStefano Zampini     PetscBool isseqaij;
2983301b35fSStefano Zampini     /* determine if we are dealing with hermitian positive definite problems */
2993301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
3003301b35fSStefano Zampini     if (Ain->symmetric_set) {
3013301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->symmetric;
3023301b35fSStefano Zampini     }
3033301b35fSStefano Zampini #else
3043301b35fSStefano Zampini     if (Ain->hermitian_set) {
3053301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->hermitian;
3063301b35fSStefano Zampini     }
3073301b35fSStefano Zampini #endif
3083301b35fSStefano Zampini     if (Ain->spd_set) {
3093301b35fSStefano Zampini       sub_schurs->is_posdef = Ain->spd;
3103301b35fSStefano Zampini     }
311a64f4aa4SStefano Zampini 
3123301b35fSStefano Zampini     /* check */
313a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
3143301b35fSStefano Zampini     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
3153301b35fSStefano Zampini       PetscInt lsize;
3163301b35fSStefano Zampini 
3173301b35fSStefano Zampini       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
3183301b35fSStefano Zampini       if (lsize) {
3193301b35fSStefano Zampini         PetscScalar val;
3203301b35fSStefano Zampini         PetscReal   norm;
3213301b35fSStefano Zampini         Vec         vec1,vec2,vec3;
3223301b35fSStefano Zampini 
3233301b35fSStefano Zampini         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
3243301b35fSStefano Zampini         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
3253301b35fSStefano Zampini         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
3263301b35fSStefano Zampini         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
3273301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
3283301b35fSStefano Zampini         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3293301b35fSStefano Zampini #else
3303301b35fSStefano Zampini         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
3313301b35fSStefano Zampini #endif
3323301b35fSStefano Zampini         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
3333301b35fSStefano Zampini         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
334862806e4SStefano Zampini         if (norm > PetscSqrtReal(PETSC_SMALL)) {
3353301b35fSStefano Zampini           sub_schurs->is_hermitian = PETSC_FALSE;
3363301b35fSStefano Zampini         } else {
3373301b35fSStefano Zampini           sub_schurs->is_hermitian = PETSC_TRUE;
3383301b35fSStefano Zampini         }
3393301b35fSStefano Zampini         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
3403301b35fSStefano Zampini         if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
3413301b35fSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
3423301b35fSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
3433301b35fSStefano Zampini         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
3443301b35fSStefano Zampini       } else {
3453301b35fSStefano Zampini         sub_schurs->is_hermitian = PETSC_TRUE;
3463301b35fSStefano Zampini         sub_schurs->is_posdef = PETSC_TRUE;
3473301b35fSStefano Zampini       }
3483301b35fSStefano Zampini     }
349a64f4aa4SStefano Zampini     if (isseqaij) {
350a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
351a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
3523301b35fSStefano Zampini     } else {
353a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
354a64f4aa4SStefano Zampini     }
355a64f4aa4SStefano Zampini   }
356*6c4ed002SBarry Smith   if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"General matrix pencils are not currently supported (%D,%D)",sub_schurs->is_hermitian,sub_schurs->is_posdef);
3573301b35fSStefano Zampini 
358a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
359a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
360a64f4aa4SStefano Zampini   if (sub_schurs->use_mumps) {
361a64f4aa4SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
362a64f4aa4SStefano Zampini   }
363a64f4aa4SStefano Zampini 
3645a95e1ceSStefano Zampini   /* preliminary checks */
365*6c4ed002SBarry Smith   if (!sub_schurs->use_mumps && compute_Stilda) SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
3665a95e1ceSStefano Zampini 
3675a95e1ceSStefano Zampini   /* restrict work on active processes */
3685a95e1ceSStefano Zampini   color = 0;
3695a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
3705a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
3715a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
3725a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
3735a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
3745a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
3755a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
3765a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3775a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
3785a95e1ceSStefano Zampini     PetscFunctionReturn(0);
3795a95e1ceSStefano Zampini   }
38081ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
3815a95e1ceSStefano Zampini 
382b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
383883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
384a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
3853301b35fSStefano Zampini     PetscBool isseqsbaij;
386a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
3873301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
3883301b35fSStefano Zampini     if (isseqsbaij) {
389a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
390a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
391a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
392a64f4aa4SStefano Zampini     } else {
393a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
394a64f4aa4SStefano Zampini       A_BB = tA_BB;
395a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
396a64f4aa4SStefano Zampini       A_IB = tA_IB;
397a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
398a64f4aa4SStefano Zampini       A_BI = tA_BI;
399f11841e3SStefano Zampini     }
400a58a30b4SStefano Zampini   } else {
4015a95e1ceSStefano Zampini     A_II = NULL;
4025a95e1ceSStefano Zampini     A_IB = NULL;
4035a95e1ceSStefano Zampini     A_BI = NULL;
4045a95e1ceSStefano Zampini     A_BB = NULL;
405b1b3d7a2SStefano Zampini   }
4065a95e1ceSStefano Zampini   S_all = NULL;
407b1b3d7a2SStefano Zampini 
408b1b3d7a2SStefano Zampini   /* determine interior problems */
4093dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4103dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
411b1b3d7a2SStefano Zampini     PetscBT                touched;
412b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
413b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
414b1b3d7a2SStefano Zampini 
415*6c4ed002SBarry Smith     if (!xadj || !adjncy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
416b1b3d7a2SStefano Zampini     /* get sizes */
417b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
418b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
419b1b3d7a2SStefano Zampini 
420b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
421b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
422b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
423b1b3d7a2SStefano Zampini 
424b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
425b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
426b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
427b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
428b1b3d7a2SStefano Zampini     }
429b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
430b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
431b1b3d7a2SStefano Zampini 
432b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
433b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
434b1b3d7a2SStefano Zampini     n_prev_added = n_B;
435b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
436b1b3d7a2SStefano Zampini       PetscInt n_added;
437b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
438*6c4ed002SBarry Smith       if (n_local_dofs > n_I+n_B) 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);
439b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
440b1b3d7a2SStefano Zampini       n_prev_added = n_added;
441b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
442b1b3d7a2SStefano Zampini       if (!n_added) break;
443b1b3d7a2SStefano Zampini     }
444b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
445b1b3d7a2SStefano Zampini 
446883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
447a9b99552SStefano 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);
448b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
449a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
450883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
451883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
452b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
453b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
454a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
455b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
456b1b3d7a2SStefano Zampini 
457b1b3d7a2SStefano Zampini       /* II block */
458b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
459b1b3d7a2SStefano Zampini     }
460b1b3d7a2SStefano Zampini   } else {
461b1b3d7a2SStefano Zampini     PetscInt n_I;
462b1b3d7a2SStefano Zampini 
463b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
464b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
465a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
466b1b3d7a2SStefano Zampini 
467b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
468883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
469b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
470b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
471b1b3d7a2SStefano Zampini 
472b1b3d7a2SStefano Zampini       /* II block is the same */
473b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
474b1b3d7a2SStefano Zampini       AE_II = A_II;
475b1b3d7a2SStefano Zampini     }
476b1b3d7a2SStefano Zampini   }
4775a95e1ceSStefano Zampini 
478883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
4795a95e1ceSStefano Zampini   max_subset_size = 0;
480883469d8SStefano Zampini   local_size = 0;
4815a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
4825a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
4835a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
484883469d8SStefano Zampini     local_size += subset_size;
485883469d8SStefano Zampini   }
486883469d8SStefano Zampini 
487883469d8SStefano Zampini   /* Work arrays for local indices */
488883469d8SStefano Zampini   extra = 0;
48981ea8064SStefano Zampini   if (sub_schurs->use_mumps && is_I_layer) {
490a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
491883469d8SStefano Zampini   }
492883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
493883469d8SStefano Zampini   if (extra) {
494883469d8SStefano Zampini     const PetscInt *idxs;
495a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
496883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
497a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
498883469d8SStefano Zampini   }
499883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
500dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
501dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
502883469d8SStefano Zampini 
503883469d8SStefano Zampini   /* Get local indices in local numbering */
504883469d8SStefano Zampini   local_size = 0;
5055a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
506883469d8SStefano Zampini     PetscInt j;
507883469d8SStefano Zampini     const    PetscInt *idxs;
508883469d8SStefano Zampini 
5095a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5105a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
511eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
512eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
513eb595f79SStefano Zampini     auxnum2[i] = subset_size;
514883469d8SStefano Zampini     /* subset indices in local numbering */
515883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5165a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
517883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
518883469d8SStefano Zampini     local_size += subset_size;
519883469d8SStefano Zampini   }
520883469d8SStefano Zampini 
5215a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
522d2627357SStefano Zampini   Bwork = NULL;
523d2627357SStefano Zampini   pivots = NULL;
52481ea8064SStefano Zampini   if (local_size && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
525d2627357SStefano Zampini     PetscScalar lwork;
526d2627357SStefano Zampini 
527d2627357SStefano Zampini     B_lwork = -1;
528d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
529d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
530d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
531d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
532d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
533d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
534d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
535d2627357SStefano Zampini   }
536d2627357SStefano Zampini 
537d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
538dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
539dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
540dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
541dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
542dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
543dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
544dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
545dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
546*6c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
547dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
548e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
549d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
550d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
551d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
552d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
5532972d61bSStefano Zampini 
5545a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
5555a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
5565a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
5575a95e1ceSStefano Zampini 
5585a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
5595a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
560*6c4ed002SBarry Smith     if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D\n",subset_size,local_size);
5615a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
562b1b3d7a2SStefano Zampini   }
563b1b3d7a2SStefano Zampini 
5645a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
5655a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
5665a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
5675a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
5685a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
5695a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
570aa83b6aeSStefano Zampini   }
571b1b3d7a2SStefano Zampini 
5725a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
573be83ff47SStefano Zampini   F = NULL;
5745a95e1ceSStefano Zampini   if (!sub_schurs->use_mumps) {
5755a95e1ceSStefano Zampini     Mat         S_Ej_expl;
5765a95e1ceSStefano Zampini     PetscScalar *work;
5775a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
5785a95e1ceSStefano Zampini     PetscBool   Sdense;
5795a95e1ceSStefano Zampini 
5805a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
5815a95e1ceSStefano Zampini     local_size = 0;
582b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
5835a95e1ceSStefano Zampini       IS  is_subset_B;
5845a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
5855a95e1ceSStefano Zampini 
5865a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
5875a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
5885a95e1ceSStefano Zampini       /* EE block */
5895a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
5905a95e1ceSStefano Zampini       /* IE block */
5915a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
5925a95e1ceSStefano Zampini       /* EI block */
5935a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
5945a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
5955a95e1ceSStefano Zampini       } else {
5965a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
5975a95e1ceSStefano Zampini       }
598a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
5995a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
6005a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
6015a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
6025a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
603b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
604b1b3d7a2SStefano Zampini         KSP ksp;
605b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
6065a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
607b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
608b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
609b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
610b1b3d7a2SStefano Zampini         KSPType   ksp_type;
611b1b3d7a2SStefano Zampini         PetscInt  n_internal;
6125a95e1ceSStefano Zampini         PetscBool ispcnone;
613b1b3d7a2SStefano Zampini 
614b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
6155a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
616b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
617b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
618b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
619b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
6205a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
6215a95e1ceSStefano Zampini         if (!ispcnone) {
6225a95e1ceSStefano Zampini           PCType pc_type;
623b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
624b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
6255a95e1ceSStefano Zampini         } else {
6265a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
6275a95e1ceSStefano Zampini         }
628b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
629b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
630b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
631b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
632b1b3d7a2SStefano Zampini           if (solver) {
633b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
634b1b3d7a2SStefano Zampini           }
635b1b3d7a2SStefano Zampini         }
636b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
637b1b3d7a2SStefano Zampini       }
6385a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6395a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
6405a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
6415a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
6425a95e1ceSStefano Zampini       if (Sdense) {
6435a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
6445a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
645b1b3d7a2SStefano Zampini         }
6465a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
647*6c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
6485a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
649a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
6505a95e1ceSStefano Zampini       local_size += subset_size;
6515a95e1ceSStefano Zampini     }
6525a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
653b1b3d7a2SStefano Zampini     /* free */
654b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
655b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
6565a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
657883469d8SStefano Zampini   } else {
658be83ff47SStefano Zampini     Mat         A;
659883469d8SStefano Zampini     IS          is_A_all;
660be83ff47SStefano Zampini     PetscScalar *work,*S_data;
661be83ff47SStefano Zampini     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2;
6624a6c6b0dSStefano Zampini     PetscBool   mumps_S;
663883469d8SStefano Zampini 
664883469d8SStefano Zampini     /* get working mat */
66581ea8064SStefano Zampini     n_I = 0;
66681ea8064SStefano Zampini     if (is_I_layer) {
667a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
66881ea8064SStefano Zampini     }
669d62866d3SStefano Zampini     if (!sub_schurs->is_dir) {
670883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
671d62866d3SStefano Zampini       size_schur = local_size;
672d62866d3SStefano Zampini     } else {
673d62866d3SStefano Zampini       IS list[2];
674d62866d3SStefano Zampini 
675d62866d3SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr);
676d62866d3SStefano Zampini       list[1] = sub_schurs->is_dir;
677d62866d3SStefano Zampini       ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr);
678d62866d3SStefano Zampini       ierr = ISDestroy(&list[0]);CHKERRQ(ierr);
679d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr);
680d62866d3SStefano Zampini       size_schur += local_size;
681d62866d3SStefano Zampini     }
68253892102SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
68353892102SStefano Zampini     size_active_schur = local_size; /* size active schurs does not count any dirichlet dof on the interface */
6846816873aSStefano Zampini     ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
685a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
6863301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
6873301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
6883301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
689883469d8SStefano Zampini 
69008122e43SStefano Zampini     if (n_I) {
6919ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
692883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
693883469d8SStefano Zampini       } else {
694883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
695883469d8SStefano Zampini       }
696883469d8SStefano Zampini       /* subsets ordered last */
697d62866d3SStefano Zampini       ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr);
698d62866d3SStefano Zampini       for (i=0;i<size_schur;i++) {
699883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
700883469d8SStefano Zampini       }
7015a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
702d62866d3SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr);
7035a95e1ceSStefano Zampini #endif
704883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
705883469d8SStefano Zampini 
706883469d8SStefano Zampini       /* factorization step */
7079ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
708883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
709be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
710be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
711be83ff47SStefano Zampini #endif
712883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
713883469d8SStefano Zampini       } else {
714883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
715be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
716be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
717be83ff47SStefano Zampini #endif
718883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
719883469d8SStefano Zampini       }
720883469d8SStefano Zampini 
721883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
7225a95e1ceSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
723883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
7245a95e1ceSStefano Zampini #endif
725d5574798SStefano Zampini 
726d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
727d5574798SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
728be83ff47SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all));
7294a6c6b0dSStefano Zampini       mumps_S = PETSC_TRUE;
730be83ff47SStefano Zampini     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
731be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
732be83ff47SStefano Zampini       reuse_solvers = PETSC_FALSE;
7334a6c6b0dSStefano Zampini       mumps_S = PETSC_FALSE;
734be83ff47SStefano Zampini     }
735be83ff47SStefano Zampini 
736be83ff47SStefano Zampini     if (reuse_solvers) {
737d62866d3SStefano Zampini       Mat              A_II;
73853892102SStefano Zampini       Vec              vec1_B;
739d62866d3SStefano Zampini       PCBDDCReuseMumps msolv_ctx;
740d5574798SStefano Zampini 
741d62866d3SStefano Zampini       if (sub_schurs->reuse_mumps) {
7426816873aSStefano Zampini         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
743e28d306cSStefano Zampini       } else {
744e28d306cSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
745d62866d3SStefano Zampini       }
746e28d306cSStefano Zampini       msolv_ctx = sub_schurs->reuse_mumps;
747be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
748be83ff47SStefano Zampini       ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr);
749d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
750d5574798SStefano Zampini       msolv_ctx->F = F;
751d5574798SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
752d62866d3SStefano Zampini 
753d62866d3SStefano Zampini       /* interior solver */
754d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr);
755d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
756d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
757d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
758d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
759e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
760d62866d3SStefano Zampini 
761d62866d3SStefano Zampini       /* correction solver */
762d62866d3SStefano Zampini       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr);
763d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
764d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
765d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
766d62866d3SStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
767e28d306cSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
76853892102SStefano Zampini 
76953892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
77053892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
77153892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
77253892102SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
77353892102SStefano Zampini       ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
77453892102SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
77553892102SStefano Zampini       ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
77653892102SStefano Zampini       msolv_ctx->is_R = is_A_all;
777d5574798SStefano Zampini     }
77808122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
77953892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
7805db18549SStefano Zampini 
781be83ff47SStefano Zampini     /* Work arrays */
782be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
7835a94c880SStefano Zampini 
7845a94c880SStefano Zampini     /* matrices for adaptive selection */
78512d906b1SStefano Zampini     if (compute_Stilda) {
7865a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
7875a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
7885a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
7895a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
7905a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
7915a94c880SStefano Zampini       }
7925a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all) {
7935a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
7945a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
7955a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
7965a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
7975a94c880SStefano Zampini       }
79812d906b1SStefano Zampini     }
799d2627357SStefano Zampini 
800be83ff47SStefano Zampini     /* S_Ej_all */
801be83ff47SStefano Zampini     cum = cum2 = 0;
802be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
80365d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
80465d8bf0aSStefano Zampini       PetscInt j;
80565d8bf0aSStefano Zampini 
8065a95e1ceSStefano Zampini       /* get S_E */
807b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
8084a6c6b0dSStefano Zampini       if (mumps_S && sub_schurs->is_hermitian) { /* When using MUMPS data I need to expand to upper triangular (column oriented) */
809be83ff47SStefano Zampini         PetscInt k;
810be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
811be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
812be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
813be83ff47SStefano Zampini             work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
814be83ff47SStefano Zampini           }
815be83ff47SStefano Zampini         }
816be83ff47SStefano Zampini       } else { /* copy to workspace */
817be83ff47SStefano Zampini         PetscInt k;
818be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
819be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
820be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
821be83ff47SStefano Zampini           }
822be83ff47SStefano Zampini         }
8239087bf02SStefano Zampini       }
8245a95e1ceSStefano Zampini       /* insert S_E values */
825be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
8265a95e1ceSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
827a1337663SStefano Zampini 
828be83ff47SStefano Zampini       /* if adaptivity is requested, invert S_E block */
829862806e4SStefano Zampini       if (compute_Stilda) {
83008122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
83108122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
832d6462365SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { /* TODO add sytrf/i for symmetric non hermitian */
8335a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
8342972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
8355a95e1ceSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
8362972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
837d6462365SStefano Zampini         } else {
838d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
839d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
840d6462365SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
841d6462365SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
8422972d61bSStefano Zampini         }
84308122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
8445a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8459087bf02SStefano Zampini       }
846be83ff47SStefano Zampini       cum += subset_size;
847be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
848be83ff47SStefano Zampini     }
849be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
8504a6c6b0dSStefano Zampini 
851be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
8524a6c6b0dSStefano Zampini     if (mumps_S) {
853be83ff47SStefano Zampini       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
8544a6c6b0dSStefano Zampini     }
8554a6c6b0dSStefano Zampini #endif
8564a6c6b0dSStefano Zampini 
85745951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
8584a6c6b0dSStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
8594a6c6b0dSStefano Zampini         PetscInt j;
8604a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
8615a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8624a6c6b0dSStefano Zampini       } else {
8634a6c6b0dSStefano Zampini         if (mumps_S) { /* use MatMumps calls to invert S */
8644a6c6b0dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
865be83ff47SStefano Zampini           ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr);
866be83ff47SStefano Zampini           ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
8674a6c6b0dSStefano Zampini #endif
8684a6c6b0dSStefano Zampini         } else { /* we need to invert explicitly since we are not using MUMPS for S */
869be83ff47SStefano Zampini           ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
870be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
871be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
872d6462365SStefano Zampini           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { /* TODO add sytrf/i for symmetric non hermitian */
873be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
874be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
875be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
876be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
877d6462365SStefano Zampini           } else {
878d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
879d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
880d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
881d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
882be83ff47SStefano Zampini           }
883be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
884be83ff47SStefano Zampini           ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
885be83ff47SStefano Zampini         }
886be83ff47SStefano Zampini         /* S_Ej_tilda_all */
887be83ff47SStefano Zampini         cum = cum2 = 0;
888be83ff47SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
889be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
890be83ff47SStefano Zampini           PetscInt j;
891862806e4SStefano Zampini 
892862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
893be83ff47SStefano Zampini           /* get (St^-1)_E */
894aff50787SStefano Zampini           if (sub_schurs->is_hermitian) { /* Here I don't need to expand to upper triangular (column oriented) */
895be83ff47SStefano Zampini             PetscInt k;
896be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
897be83ff47SStefano Zampini               for (j=k;j<subset_size;j++) {
898be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
899be83ff47SStefano Zampini               }
900be83ff47SStefano Zampini             }
901be83ff47SStefano Zampini           } else { /* copy to workspace */
902be83ff47SStefano Zampini             PetscInt k;
903be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
904be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
905be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
906be83ff47SStefano Zampini               }
907be83ff47SStefano Zampini             }
908be83ff47SStefano Zampini           }
909be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
9105a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
911be83ff47SStefano Zampini           cum += subset_size;
912be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
913883469d8SStefano Zampini         }
914be83ff47SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
915be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
9164a6c6b0dSStefano Zampini         if (mumps_S) {
917be83ff47SStefano Zampini           ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
9185db18549SStefano Zampini         }
919be83ff47SStefano Zampini #endif
9209087bf02SStefano Zampini       }
9214a6c6b0dSStefano Zampini     }
9224a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
9234a6c6b0dSStefano Zampini   }
9245a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
925be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
926a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
927a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
928a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
929a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
930a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
9315db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9325db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9335a95e1ceSStefano Zampini   if (compute_Stilda) {
9345a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9355a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9365a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9375a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93808122e43SStefano Zampini   }
939a1337663SStefano Zampini 
9405db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
9415db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
9423927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
9433927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
9445a95e1ceSStefano Zampini 
9455db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
9465a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
947dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
9485a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
9495a94c880SStefano Zampini   } else {
9505a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
9515a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
952dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
9535a94c880SStefano Zampini   }
95408122e43SStefano Zampini 
955ac632422SStefano Zampini   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
956ac632422SStefano Zampini   if (faster_deluxe) {
9575a95e1ceSStefano Zampini     Mat         tmpmat;
9585a95e1ceSStefano Zampini     PetscScalar *array;
9595a95e1ceSStefano Zampini     PetscInt    cum;
9605a95e1ceSStefano Zampini 
9615a95e1ceSStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
9625a95e1ceSStefano Zampini     cum = 0;
9635a95e1ceSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
9645a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
9655a95e1ceSStefano Zampini       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
9665a95e1ceSStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
967d6462365SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
9685a95e1ceSStefano Zampini         PetscInt j,k;
9695a95e1ceSStefano Zampini 
9705a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
9715a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
9725a95e1ceSStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
9735a95e1ceSStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
9745a95e1ceSStefano Zampini         for (j=0;j<B_N;j++) {
9755a95e1ceSStefano Zampini           for (k=j+1;k<B_N;k++) {
9765a95e1ceSStefano Zampini             array[k*B_N+j+cum] = array[j*B_N+k+cum];
9775a95e1ceSStefano Zampini           }
9785a95e1ceSStefano Zampini         }
979d6462365SStefano Zampini       } else {
980d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
981d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
982d6462365SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
983d6462365SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
9845a95e1ceSStefano Zampini       }
9855a95e1ceSStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
9865a95e1ceSStefano Zampini       cum += subset_size*subset_size;
9875a95e1ceSStefano Zampini     }
9885a95e1ceSStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
9895a95e1ceSStefano Zampini     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
9905a95e1ceSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
991ac632422SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
9925a95e1ceSStefano Zampini     sub_schurs->S_Ej_all = tmpmat;
9935a95e1ceSStefano Zampini   }
9945a95e1ceSStefano Zampini 
995f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
9965a95e1ceSStefano Zampini   if (compute_Stilda) {
9975a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
998a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
9995a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1000dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
10015a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
100208122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
10035a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_inv_all;
1004dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
100508122e43SStefano Zampini   }
10063202ece2SStefano Zampini 
10075a95e1ceSStefano Zampini   /* free workspace */
10085a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
100906a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1010a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
10113202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1012dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
10135a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1014b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1015b1b3d7a2SStefano Zampini }
1016b1b3d7a2SStefano Zampini 
1017b1b3d7a2SStefano Zampini #undef __FUNCT__
1018b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
1019a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1020b1b3d7a2SStefano Zampini {
10219bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
10225a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1023b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1024b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1025b1b3d7a2SStefano Zampini 
1026b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1027b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1028*6c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1029b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1030*6c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1031b1b3d7a2SStefano Zampini 
1032b1b3d7a2SStefano Zampini   /* reset any previous data */
1033b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1034b1b3d7a2SStefano Zampini 
10355a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
10369bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1037b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
103808122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1039b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1040b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
1041b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
1042b1b3d7a2SStefano Zampini   }
1043b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
1044b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
104508122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1046b1b3d7a2SStefano Zampini   }
1047b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
1048b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
1049d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1050d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1051b1b3d7a2SStefano Zampini 
10526816873aSStefano Zampini   /* Determine if MUMPS can be used */
1053883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
1054883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1055a64f4aa4SStefano Zampini   sub_schurs->use_mumps = PETSC_TRUE;
1056883469d8SStefano Zampini #endif
1057b1b3d7a2SStefano Zampini 
1058b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1059b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1060b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1061b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
10625db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
10635db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
10645db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
10655db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
10665a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1067b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1068a64f4aa4SStefano Zampini   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1069b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1070b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1071b96c3477SStefano Zampini     }
10729bb4a8caSStefano Zampini   }
1073d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
1074b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1075b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
107608122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1077b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1078b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1079b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1080b1b3d7a2SStefano Zampini }
1081b1b3d7a2SStefano Zampini 
1082b1b3d7a2SStefano Zampini #undef __FUNCT__
108334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
108434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
108534a97f8cSStefano Zampini {
108634a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
108734a97f8cSStefano Zampini   PetscErrorCode  ierr;
108834a97f8cSStefano Zampini 
108934a97f8cSStefano Zampini   PetscFunctionBegin;
109034a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
10915ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
109234a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
109334a97f8cSStefano Zampini   PetscFunctionReturn(0);
109434a97f8cSStefano Zampini }
109534a97f8cSStefano Zampini 
109634a97f8cSStefano Zampini #undef __FUNCT__
109734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
109834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
109934a97f8cSStefano Zampini {
110034a97f8cSStefano Zampini   PetscErrorCode ierr;
110134a97f8cSStefano Zampini 
110234a97f8cSStefano Zampini   PetscFunctionBegin;
110334a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
110434a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
110534a97f8cSStefano Zampini   PetscFunctionReturn(0);
110634a97f8cSStefano Zampini }
110734a97f8cSStefano Zampini 
110834a97f8cSStefano Zampini #undef __FUNCT__
110934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
111034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
111134a97f8cSStefano Zampini {
111234a97f8cSStefano Zampini   PetscInt       i;
111334a97f8cSStefano Zampini   PetscErrorCode ierr;
111434a97f8cSStefano Zampini 
111534a97f8cSStefano Zampini   PetscFunctionBegin;
11161e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1117b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1118b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1119b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
11205db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
11215db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
112241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
112341c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
112408122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1125a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11265db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1127d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1128d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
112908122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
113034a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1131b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
113234a97f8cSStefano Zampini   }
11335ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1134b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
11353dc780c3SStefano Zampini   }
1136d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps) {
1137d62866d3SStefano Zampini     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1138d62866d3SStefano Zampini   }
1139d62866d3SStefano Zampini   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
114034a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
114134a97f8cSStefano Zampini   PetscFunctionReturn(0);
114234a97f8cSStefano Zampini }
114334a97f8cSStefano Zampini 
114434a97f8cSStefano Zampini #undef __FUNCT__
114534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
11462a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
114734a97f8cSStefano Zampini {
114834a97f8cSStefano Zampini   PetscInt       i,j,n;
114934a97f8cSStefano Zampini   PetscErrorCode ierr;
115034a97f8cSStefano Zampini 
115134a97f8cSStefano Zampini   PetscFunctionBegin;
115234a97f8cSStefano Zampini   n = 0;
115334a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
115434a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
115534a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
115634a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
115734a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
115834a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
115934a97f8cSStefano Zampini         queue_tip[n] = dof;
116034a97f8cSStefano Zampini         n++;
116134a97f8cSStefano Zampini       }
116234a97f8cSStefano Zampini     }
116334a97f8cSStefano Zampini   }
116434a97f8cSStefano Zampini   *n_added = n;
116534a97f8cSStefano Zampini   PetscFunctionReturn(0);
116634a97f8cSStefano Zampini }
1167