xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision b3cb21dd70d8f158736d96e11785e5dd22772d64)
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*);
7df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
8df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
9d62866d3SStefano Zampini 
10ca92afb2SStefano Zampini /* if v2 is not present, correction is done in-place */
115cbda25cSStefano Zampini PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
12ca92afb2SStefano Zampini {
13ca92afb2SStefano Zampini   PetscScalar    *array;
14ca92afb2SStefano Zampini   PetscScalar    *array2;
15ca92afb2SStefano Zampini   PetscErrorCode ierr;
16ca92afb2SStefano Zampini 
17ca92afb2SStefano Zampini   PetscFunctionBegin;
18ca92afb2SStefano Zampini   if (!ctx->benign_n) PetscFunctionReturn(0);
195cbda25cSStefano Zampini   if (sol && full) {
205cbda25cSStefano Zampini     PetscInt n_I,size_schur;
215cbda25cSStefano Zampini 
225cbda25cSStefano Zampini     /* get sizes */
23282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
245cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
255cbda25cSStefano Zampini     n_I = n_I - size_schur;
265cbda25cSStefano Zampini     /* get schur sol from array */
275cbda25cSStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
285cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
295cbda25cSStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
305cbda25cSStefano Zampini     /* apply interior sol correction */
31282d6408SStefano Zampini     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
325cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
335cbda25cSStefano Zampini     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
345cbda25cSStefano Zampini   }
35ca92afb2SStefano Zampini   if (v2) {
36ca92afb2SStefano Zampini     PetscInt nl;
37ca92afb2SStefano Zampini 
38ca92afb2SStefano Zampini     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
39ca92afb2SStefano Zampini     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
40ca92afb2SStefano Zampini     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
41ca92afb2SStefano Zampini     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
42ca92afb2SStefano Zampini   } else {
43ca92afb2SStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
44ca92afb2SStefano Zampini     array2 = array;
45ca92afb2SStefano Zampini   }
46ca92afb2SStefano Zampini   if (!sol) { /* change rhs */
47ca92afb2SStefano Zampini     PetscInt n;
48ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
49ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
50ca92afb2SStefano Zampini       const PetscInt *cols;
51ca92afb2SStefano Zampini       PetscInt       nz,i;
52ca92afb2SStefano Zampini 
53ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
54ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
55ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
5622db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
5722db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
5822db5ddcSStefano Zampini #else
59ca92afb2SStefano Zampini       sum = -sum/nz;
6022db5ddcSStefano Zampini #endif
61ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
62ca92afb2SStefano Zampini       ctx->benign_save_vals[n] = array2[cols[nz-1]];
63ca92afb2SStefano Zampini       array2[cols[nz-1]] = sum;
64ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
65ca92afb2SStefano Zampini     }
66ca92afb2SStefano Zampini   } else {
67ca92afb2SStefano Zampini     PetscInt n;
68ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
69ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
70ca92afb2SStefano Zampini       const PetscInt *cols;
71ca92afb2SStefano Zampini       PetscInt       nz,i;
72ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
73ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
74ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
7522db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
7622db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
7722db5ddcSStefano Zampini #else
78ca92afb2SStefano Zampini       sum = -sum/nz;
7922db5ddcSStefano Zampini #endif
80ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
81ca92afb2SStefano Zampini       array2[cols[nz-1]] = ctx->benign_save_vals[n];
82ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
83ca92afb2SStefano Zampini     }
84ca92afb2SStefano Zampini   }
85ca92afb2SStefano Zampini   if (v2) {
86ca92afb2SStefano Zampini     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
8772b8c272SStefano Zampini     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
88ca92afb2SStefano Zampini   } else {
89ca92afb2SStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
90ca92afb2SStefano Zampini   }
915cbda25cSStefano Zampini   if (!sol && full) {
925cbda25cSStefano Zampini     Vec      usedv;
935cbda25cSStefano Zampini     PetscInt n_I,size_schur;
945cbda25cSStefano Zampini 
955cbda25cSStefano Zampini     /* get sizes */
96282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
975cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
985cbda25cSStefano Zampini     n_I = n_I - size_schur;
995cbda25cSStefano Zampini     /* compute schur rhs correction */
1005cbda25cSStefano Zampini     if (v2) {
1015cbda25cSStefano Zampini       usedv = v2;
1025cbda25cSStefano Zampini     } else {
1035cbda25cSStefano Zampini       usedv = v;
1045cbda25cSStefano Zampini     }
1055cbda25cSStefano Zampini     /* apply schur rhs correction */
1065cbda25cSStefano Zampini     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
1075cbda25cSStefano Zampini     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1085cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
1095cbda25cSStefano Zampini     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
110282d6408SStefano Zampini     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1115cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1125cbda25cSStefano Zampini   }
113ca92afb2SStefano Zampini   PetscFunctionReturn(0);
114ca92afb2SStefano Zampini }
115ca92afb2SStefano Zampini 
116df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
117d62866d3SStefano Zampini {
118df4d28bfSStefano Zampini   PCBDDCReuseSolvers ctx;
119683d3df6SStefano Zampini   PetscBool          copy = PETSC_FALSE;
120d62866d3SStefano Zampini   PetscErrorCode     ierr;
121d62866d3SStefano Zampini 
122d62866d3SStefano Zampini   PetscFunctionBegin;
123d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
124683d3df6SStefano Zampini   if (full) {
125d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
126d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
127d62866d3SStefano Zampini #endif
1285cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1295cbda25cSStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
1305cbda25cSStefano Zampini #endif
131683d3df6SStefano Zampini     copy = ctx->has_vertices;
132d4933d67SStefano Zampini   } else { /* interior solver */
1336dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
134683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
1356dba178dSStefano Zampini #endif
136d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
137d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
138d4933d67SStefano Zampini #endif
139683d3df6SStefano Zampini     copy = PETSC_TRUE;
140683d3df6SStefano Zampini   }
141683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
142683d3df6SStefano Zampini   if (copy) {
143ca92afb2SStefano Zampini     PetscInt    n;
144df4d28bfSStefano Zampini     PetscScalar *array,*array_solver;
145ca92afb2SStefano Zampini 
146ca92afb2SStefano Zampini     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
147683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
148df4d28bfSStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
149df4d28bfSStefano Zampini     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
150df4d28bfSStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
151683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
152683d3df6SStefano Zampini 
1535cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
154683d3df6SStefano Zampini     if (transpose) {
155683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
156683d3df6SStefano Zampini     } else {
157683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
158683d3df6SStefano Zampini     }
1595cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
160683d3df6SStefano Zampini 
161683d3df6SStefano Zampini     /* get back data to caller worskpace */
162df4d28bfSStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
163683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
164df4d28bfSStefano Zampini     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
165683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
166df4d28bfSStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
167683d3df6SStefano Zampini   } else {
168ca92afb2SStefano Zampini     if (ctx->benign_n) {
1695cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
170ca92afb2SStefano Zampini       if (transpose) {
171ca92afb2SStefano Zampini         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
172ca92afb2SStefano Zampini       } else {
173ca92afb2SStefano Zampini         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
174ca92afb2SStefano Zampini       }
1755cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
176ca92afb2SStefano Zampini     } else {
177e28d306cSStefano Zampini       if (transpose) {
178e28d306cSStefano Zampini         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
179e28d306cSStefano Zampini       } else {
1806816873aSStefano Zampini         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
181e28d306cSStefano Zampini       }
182683d3df6SStefano Zampini     }
183ca92afb2SStefano Zampini   }
1845cbda25cSStefano Zampini   /* restore defaults */
1855cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1865cbda25cSStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
1875cbda25cSStefano Zampini #endif
188d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
189d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
190d4933d67SStefano Zampini #endif
191d62866d3SStefano Zampini   PetscFunctionReturn(0);
192d62866d3SStefano Zampini }
193d62866d3SStefano Zampini 
194df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
195e28d306cSStefano Zampini {
196e28d306cSStefano Zampini   PetscErrorCode   ierr;
197e28d306cSStefano Zampini 
198e28d306cSStefano Zampini   PetscFunctionBegin;
199df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
200e28d306cSStefano Zampini   PetscFunctionReturn(0);
201e28d306cSStefano Zampini }
202e28d306cSStefano Zampini 
203df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
204e28d306cSStefano Zampini {
205e28d306cSStefano Zampini   PetscErrorCode   ierr;
206e28d306cSStefano Zampini 
207e28d306cSStefano Zampini   PetscFunctionBegin;
208df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
209683d3df6SStefano Zampini   PetscFunctionReturn(0);
210683d3df6SStefano Zampini }
211683d3df6SStefano Zampini 
212df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
213683d3df6SStefano Zampini {
214683d3df6SStefano Zampini   PetscErrorCode   ierr;
215683d3df6SStefano Zampini 
216683d3df6SStefano Zampini   PetscFunctionBegin;
217df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
218683d3df6SStefano Zampini   PetscFunctionReturn(0);
219683d3df6SStefano Zampini }
220683d3df6SStefano Zampini 
221df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
222683d3df6SStefano Zampini {
223683d3df6SStefano Zampini   PetscErrorCode   ierr;
224683d3df6SStefano Zampini 
225683d3df6SStefano Zampini   PetscFunctionBegin;
226df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
227e28d306cSStefano Zampini   PetscFunctionReturn(0);
228e28d306cSStefano Zampini }
229e28d306cSStefano Zampini 
230df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
231d62866d3SStefano Zampini {
232ca92afb2SStefano Zampini   PetscInt       i;
233d62866d3SStefano Zampini   PetscErrorCode ierr;
234d62866d3SStefano Zampini 
235d62866d3SStefano Zampini   PetscFunctionBegin;
236d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
237d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
238d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
239d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
240d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
24153892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
24253892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
243d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
24453892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
24553892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
246ca92afb2SStefano Zampini   for (i=0;i<reuse->benign_n;i++) {
247ca92afb2SStefano Zampini     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
248ca92afb2SStefano Zampini   }
249ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
250ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
2515cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
2525cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
2535cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
2545cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
255d62866d3SStefano Zampini   PetscFunctionReturn(0);
256d62866d3SStefano Zampini }
257d5574798SStefano Zampini 
2585ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
2593202ece2SStefano Zampini {
2603202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
2613202ece2SStefano Zampini   KSP            ksp;
2623202ece2SStefano Zampini   PC             pc;
2633202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
2643202ece2SStefano Zampini   PetscReal      fill = 2.0;
265f11841e3SStefano Zampini   PetscInt       n_I;
2663202ece2SStefano Zampini   PetscMPIInt    size;
2673202ece2SStefano Zampini   PetscErrorCode ierr;
2683202ece2SStefano Zampini 
2693202ece2SStefano Zampini   PetscFunctionBegin;
2703202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
2716c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
272f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
273f11841e3SStefano Zampini     PetscBool Sdense;
274f11841e3SStefano Zampini 
275f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
2766c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
277f11841e3SStefano Zampini   }
2783202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
2793202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
2803202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
2813202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
2823202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
2833202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
2843202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
2853202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
286f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
287f11841e3SStefano Zampini   if (n_I) {
2883202ece2SStefano Zampini     if (!Bdense) {
2893202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
2903202ece2SStefano Zampini     } else {
2913202ece2SStefano Zampini       Bd = B;
2923202ece2SStefano Zampini     }
2933202ece2SStefano Zampini 
2943202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
2953202ece2SStefano Zampini       Mat fact;
2963202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
2973202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
2983202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
2993202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
3003202ece2SStefano Zampini     } else {
30107b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
30207b1e237SStefano Zampini 
30307b1e237SStefano Zampini       if (ex) {
3043202ece2SStefano Zampini         Mat Ainvd;
3053202ece2SStefano Zampini 
3063202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
3073202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
3083202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
30907b1e237SStefano Zampini       } else {
31007b1e237SStefano Zampini         Vec         sol,rhs;
31107b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
31207b1e237SStefano Zampini         PetscInt    i,nrhs,n;
31307b1e237SStefano Zampini 
31407b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
31507b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
31607b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
31707b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
31807b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
31907b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
32007b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
32107b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
32207b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
32307b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
32407b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
32507b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
32607b1e237SStefano Zampini         }
32707b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
32807b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
32907b1e237SStefano Zampini       }
3303202ece2SStefano Zampini     }
3315ec10c6aSStefano Zampini     if (!Bdense & !issym) {
3323202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3333202ece2SStefano Zampini     }
3345ec10c6aSStefano Zampini 
3355ec10c6aSStefano Zampini     if (!issym) {
3363202ece2SStefano Zampini       if (!Cdense) {
3373202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
3383202ece2SStefano Zampini       } else {
3393202ece2SStefano Zampini         Cd = C;
3403202ece2SStefano Zampini       }
3415ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3423202ece2SStefano Zampini       if (!Cdense) {
3433202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
3443202ece2SStefano Zampini       }
3455ec10c6aSStefano Zampini     } else {
3465ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3475ec10c6aSStefano Zampini       if (!Bdense) {
3485ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3495ec10c6aSStefano Zampini       }
3505ec10c6aSStefano Zampini     }
3515ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
352f11841e3SStefano Zampini   }
3533202ece2SStefano Zampini 
3543202ece2SStefano Zampini   if (D) {
3553202ece2SStefano Zampini     Mat       Dd;
3563202ece2SStefano Zampini     PetscBool Ddense;
3573202ece2SStefano Zampini 
3583202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
3593202ece2SStefano Zampini     if (!Ddense) {
3603202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
3613202ece2SStefano Zampini     } else {
3623202ece2SStefano Zampini       Dd = D;
3633202ece2SStefano Zampini     }
364f11841e3SStefano Zampini     if (n_I) {
3653202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
366f11841e3SStefano Zampini     } else {
367f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
368f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
369f11841e3SStefano Zampini       } else {
370f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
371f11841e3SStefano Zampini       }
372f11841e3SStefano Zampini     }
3733202ece2SStefano Zampini     if (!Ddense) {
3743202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
3753202ece2SStefano Zampini     }
3763202ece2SStefano Zampini   } else {
3773202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
3783202ece2SStefano Zampini   }
3793202ece2SStefano Zampini   PetscFunctionReturn(0);
3803202ece2SStefano Zampini }
38134a97f8cSStefano Zampini 
38291af6908SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal)
383b1b3d7a2SStefano Zampini {
384be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
385be83ff47SStefano Zampini   Mat                    S_all;
3865a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
3875db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
388b7ab4a40SStefano Zampini   IS                     is_I,is_I_layer;
389dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
390dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
391dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
3925a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
393683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
39408122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
39506a4b1faSStefano Zampini   PetscScalar            *Bwork;
3965a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
3975a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
3985a95e1ceSStefano Zampini   MPI_Comm               comm_n;
3999d54b7f4SStefano Zampini   PetscBool              deluxe = PETSC_TRUE, use_getr = PETSC_FALSE;
400b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
401b1b3d7a2SStefano Zampini 
402b1b3d7a2SStefano Zampini   PetscFunctionBegin;
403a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
404a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
405a64f4aa4SStefano Zampini 
4063301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_TRUE;
4073301b35fSStefano Zampini   sub_schurs->is_posdef    = PETSC_TRUE;
408e62b6521Sstefano_zampini   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
409e62b6521Sstefano_zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
410e62b6521Sstefano_zampini   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
411e62b6521Sstefano_zampini   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
412e62b6521Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
413e62b6521Sstefano_zampini 
414e62b6521Sstefano_zampini   /* convert matrix if needed */
415e62b6521Sstefano_zampini   if (Ain) {
416e62b6521Sstefano_zampini     PetscBool isseqaij;
417e62b6521Sstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
418a64f4aa4SStefano Zampini     if (isseqaij) {
419a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
420a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4213301b35fSStefano Zampini     } else {
422a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
423a64f4aa4SStefano Zampini     }
424a64f4aa4SStefano Zampini   }
4253301b35fSStefano Zampini 
426a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
427a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
428df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
429df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
430a64f4aa4SStefano Zampini   }
431a64f4aa4SStefano Zampini 
4325a95e1ceSStefano Zampini   /* preliminary checks */
433af25d912SStefano Zampini   if (!sub_schurs->schur_explicit && compute_Stilda) SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO");
4345a95e1ceSStefano Zampini 
4355a95e1ceSStefano Zampini   /* restrict work on active processes */
4365a95e1ceSStefano Zampini   color = 0;
4375a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4385a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4395a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4405a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4415a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4425a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
4435a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
4445a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4455a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
4465a95e1ceSStefano Zampini     PetscFunctionReturn(0);
4475a95e1ceSStefano Zampini   }
44881ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
4495a95e1ceSStefano Zampini 
450b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
451df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
452a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
4533301b35fSStefano Zampini     PetscBool isseqsbaij;
454a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
4553301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
4563301b35fSStefano Zampini     if (isseqsbaij) {
457a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
458a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
459a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
460a64f4aa4SStefano Zampini     } else {
461a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
462a64f4aa4SStefano Zampini       A_BB = tA_BB;
463a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
464a64f4aa4SStefano Zampini       A_IB = tA_IB;
465a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
466a64f4aa4SStefano Zampini       A_BI = tA_BI;
467f11841e3SStefano Zampini     }
468a58a30b4SStefano Zampini   } else {
4695a95e1ceSStefano Zampini     A_II = NULL;
4705a95e1ceSStefano Zampini     A_IB = NULL;
4715a95e1ceSStefano Zampini     A_BI = NULL;
4725a95e1ceSStefano Zampini     A_BB = NULL;
473b1b3d7a2SStefano Zampini   }
4745a95e1ceSStefano Zampini   S_all = NULL;
475b1b3d7a2SStefano Zampini 
476b1b3d7a2SStefano Zampini   /* determine interior problems */
4773dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4783dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
479b1b3d7a2SStefano Zampini     PetscBT                touched;
480b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
481b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
482b1b3d7a2SStefano Zampini 
4831633d1f0SStefano Zampini     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
484b1b3d7a2SStefano Zampini     /* get sizes */
485b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
486b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
487b1b3d7a2SStefano Zampini 
488b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
489b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
490b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
491b1b3d7a2SStefano Zampini 
492b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
493b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
494b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
495b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
496b1b3d7a2SStefano Zampini     }
497b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
498b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
499b1b3d7a2SStefano Zampini 
500b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
501b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
502b1b3d7a2SStefano Zampini     n_prev_added = n_B;
503b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
504b1b3d7a2SStefano Zampini       PetscInt n_added;
505b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5066c4ed002SBarry 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);
507b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
508b1b3d7a2SStefano Zampini       n_prev_added = n_added;
509b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
510b1b3d7a2SStefano Zampini       if (!n_added) break;
511b1b3d7a2SStefano Zampini     }
512b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
513b1b3d7a2SStefano Zampini 
514883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
515a9b99552SStefano 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);
516b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
517a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
518883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
519df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
520b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
521b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
522a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
523b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
524b1b3d7a2SStefano Zampini 
525b1b3d7a2SStefano Zampini       /* II block */
5267dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
527b1b3d7a2SStefano Zampini     }
528b1b3d7a2SStefano Zampini   } else {
529b1b3d7a2SStefano Zampini     PetscInt n_I;
530b1b3d7a2SStefano Zampini 
531b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
532b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
533a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
534b1b3d7a2SStefano Zampini 
535b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
536df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
537b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
538b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
539b1b3d7a2SStefano Zampini 
540b1b3d7a2SStefano Zampini       /* II block is the same */
541b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
542b1b3d7a2SStefano Zampini       AE_II = A_II;
543b1b3d7a2SStefano Zampini     }
544b1b3d7a2SStefano Zampini   }
5455a95e1ceSStefano Zampini 
546883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
5475a95e1ceSStefano Zampini   max_subset_size = 0;
548883469d8SStefano Zampini   local_size = 0;
5495a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
5505a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5515a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
552883469d8SStefano Zampini     local_size += subset_size;
553883469d8SStefano Zampini   }
554883469d8SStefano Zampini 
555883469d8SStefano Zampini   /* Work arrays for local indices */
556883469d8SStefano Zampini   extra = 0;
557683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
558df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
559a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
560883469d8SStefano Zampini   }
561683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
562883469d8SStefano Zampini   if (extra) {
563883469d8SStefano Zampini     const PetscInt *idxs;
564a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
565883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
566a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
567883469d8SStefano Zampini   }
568883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
569dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
570dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
571883469d8SStefano Zampini 
572883469d8SStefano Zampini   /* Get local indices in local numbering */
573883469d8SStefano Zampini   local_size = 0;
5745a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
575883469d8SStefano Zampini     PetscInt j;
576883469d8SStefano Zampini     const    PetscInt *idxs;
577883469d8SStefano Zampini 
5785a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5795a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
580eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
581eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
582eb595f79SStefano Zampini     auxnum2[i] = subset_size;
583883469d8SStefano Zampini     /* subset indices in local numbering */
584883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5855a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
586883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
587883469d8SStefano Zampini     local_size += subset_size;
588883469d8SStefano Zampini   }
589883469d8SStefano Zampini 
5905a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
59106a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
59259ac4de7SStefano Zampini     PetscScalar  lwork,dummyscalar = 0.;
59359ac4de7SStefano Zampini     PetscBLASInt dummyint = 0;
594d2627357SStefano Zampini 
59506a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
596d2627357SStefano Zampini     B_lwork = -1;
597d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
598d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
59959ac4de7SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
600d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
601d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
602d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
603d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
604056290a2SStefano Zampini   } else {
605056290a2SStefano Zampini     Bwork = NULL;
606056290a2SStefano Zampini     pivots = NULL;
607d2627357SStefano Zampini   }
608d2627357SStefano Zampini 
609d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
610dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
611dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
612dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
613dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
6146583bcc1SStefano Zampini   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
615dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
616dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
617dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6186c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
619dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
620e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
621d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
622d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
623d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
624d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6252972d61bSStefano Zampini 
6265a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6275a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6285a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6295a95e1ceSStefano Zampini 
6305a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6315a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
6326c4ed002SBarry 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);
6335a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
634b1b3d7a2SStefano Zampini   }
635b1b3d7a2SStefano Zampini 
63672b8c272SStefano Zampini   if (change) {
63772b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
63872b8c272SStefano Zampini     IS                     change_primal_B;
63972b8c272SStefano Zampini     IS                     change_primal_all;
64072b8c272SStefano Zampini 
641b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
642b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
643b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
64472b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
64572b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
64672b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
647b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
64872b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
64972b8c272SStefano Zampini     }
65072b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
65172b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
65272b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
65372b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
65472b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
655b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
65672b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
65772b8c272SStefano Zampini       Mat change_sub;
65872b8c272SStefano Zampini 
659b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
66072b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
66172b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
66272b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
6637dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
66472b8c272SStefano Zampini       } else {
66572b8c272SStefano Zampini         Mat change_subt;
6667dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
66772b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
66872b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
66972b8c272SStefano Zampini       }
67072b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
67172b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
672e62b6521Sstefano_zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
673e62b6521Sstefano_zampini       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
67472b8c272SStefano Zampini     }
67572b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
67672b8c272SStefano Zampini   }
67772b8c272SStefano Zampini 
6785a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
6795a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
6805a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
6815a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
6825a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
6835a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
684aa83b6aeSStefano Zampini   }
685b1b3d7a2SStefano Zampini 
6865a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
687be83ff47SStefano Zampini   F = NULL;
688df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) { /* this code branch is used when MatFactor with Schur complemnent support is not present; it is not very efficient, unless the economic version of the scaling is required */
6895a95e1ceSStefano Zampini     Mat         S_Ej_expl;
6905a95e1ceSStefano Zampini     PetscScalar *work;
6915a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
6925a95e1ceSStefano Zampini     PetscBool   Sdense;
6935a95e1ceSStefano Zampini 
6945a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
6955a95e1ceSStefano Zampini     local_size = 0;
696b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
6975a95e1ceSStefano Zampini       IS  is_subset_B;
6985a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
6995a95e1ceSStefano Zampini 
7005a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7015a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7025a95e1ceSStefano Zampini       /* EE block */
7037dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7045a95e1ceSStefano Zampini       /* IE block */
7057dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7065a95e1ceSStefano Zampini       /* EI block */
7075a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
7085a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7095a95e1ceSStefano Zampini       } else {
7107dae84e0SHong Zhang         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7115a95e1ceSStefano Zampini       }
712a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7135a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7145a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7155a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7165a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
717b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
718b1b3d7a2SStefano Zampini         KSP ksp;
719b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7205a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
721b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
722b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
723b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
724b1b3d7a2SStefano Zampini         KSPType   ksp_type;
725b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7265a95e1ceSStefano Zampini         PetscBool ispcnone;
727b1b3d7a2SStefano Zampini 
728b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7295a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
730b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
731b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
732b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
733b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7345a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7355a95e1ceSStefano Zampini         if (!ispcnone) {
7365a95e1ceSStefano Zampini           PCType pc_type;
737b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
738b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7395a95e1ceSStefano Zampini         } else {
7405a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7415a95e1ceSStefano Zampini         }
742b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
743b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
744b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
745b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
746b1b3d7a2SStefano Zampini           if (solver) {
747b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
748b1b3d7a2SStefano Zampini           }
749b1b3d7a2SStefano Zampini         }
750b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
751b1b3d7a2SStefano Zampini       }
7525a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7535a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
7545a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
7555a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
7565a95e1ceSStefano Zampini       if (Sdense) {
7575a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
7585a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
759b1b3d7a2SStefano Zampini         }
7605a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
7616c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
7625a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
763a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
7645a95e1ceSStefano Zampini       local_size += subset_size;
7655a95e1ceSStefano Zampini     }
7665a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
767b1b3d7a2SStefano Zampini     /* free */
768b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
769b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
7705a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
771883469d8SStefano Zampini   } else {
7725cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
7739d54b7f4SStefano Zampini     Vec         Dall = NULL;
774ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
7755cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
7766dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
777df4d28bfSStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
778d4933d67SStefano Zampini     char        solver[256];
779883469d8SStefano Zampini 
780683d3df6SStefano Zampini     /* get sizes */
78181ea8064SStefano Zampini     n_I = 0;
78281ea8064SStefano Zampini     if (is_I_layer) {
783a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
78481ea8064SStefano Zampini     }
785683d3df6SStefano Zampini     economic = PETSC_FALSE;
786683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
787683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
788683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
7899d54b7f4SStefano Zampini     size_active_schur = local_size;
7909d54b7f4SStefano Zampini 
791f17d2ae1SStefano Zampini     /* import scaling vector (wrong formulation if we have 3D edges) */
7929d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
7939d54b7f4SStefano Zampini       const PetscScalar *array;
7949d54b7f4SStefano Zampini       PetscScalar       *array2;
7959d54b7f4SStefano Zampini       const PetscInt    *idxs;
7969d54b7f4SStefano Zampini       PetscInt          i;
7979d54b7f4SStefano Zampini 
7989d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
7999d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8009d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8019d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8029d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8039d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8049d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8059d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8069d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8079d54b7f4SStefano Zampini     }
808d62866d3SStefano Zampini 
809683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
810683d3df6SStefano Zampini     cum = n_I+size_active_schur;
811683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
812683d3df6SStefano Zampini       const PetscInt* idxs;
813683d3df6SStefano Zampini       PetscInt        n_dir;
814683d3df6SStefano Zampini 
815683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
816683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
817683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
818683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
819683d3df6SStefano Zampini       cum += n_dir;
820d62866d3SStefano Zampini     }
821683d3df6SStefano Zampini     factor_workaround = PETSC_FALSE;
822683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
823367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
824683d3df6SStefano Zampini       PetscInt n_v;
825683d3df6SStefano Zampini 
826683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
827683d3df6SStefano Zampini       if (n_v) {
828683d3df6SStefano Zampini         const PetscInt* idxs;
829683d3df6SStefano Zampini 
830683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
831683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
832683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
833683d3df6SStefano Zampini         cum += n_v;
834683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
835683d3df6SStefano Zampini       }
836683d3df6SStefano Zampini     }
837683d3df6SStefano Zampini     size_schur = cum - n_I;
838683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
839683d3df6SStefano Zampini     /* get working mat (TODO: factorize without actually permuting)  */
840683d3df6SStefano Zampini     if (cum == n) {
841683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
842683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
843683d3df6SStefano Zampini     } else {
8447dae84e0SHong Zhang       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
845683d3df6SStefano Zampini     }
846e62b6521Sstefano_zampini     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
847e62b6521Sstefano_zampini     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
848ca92afb2SStefano Zampini 
849ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
850367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
851ca92afb2SStefano Zampini     if (benign_n) {
852ca92afb2SStefano Zampini       Vec                    v;
853ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
854ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
8555cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
8565cbda25cSStefano Zampini       PetscInt               sizeA;
857ca92afb2SStefano Zampini 
858ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
859ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
860ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
861ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
862ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
8635cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
8645cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
865282d6408SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
8665cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
8675cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
868ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
869ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
870ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
8715cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
872ca92afb2SStefano Zampini         PetscScalar    *array;
873ca92afb2SStefano Zampini         const PetscInt *idxs;
874ca92afb2SStefano Zampini         PetscInt       j,nz;
875ca92afb2SStefano Zampini 
876ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
877ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
878ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8795cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
880ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8815cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
8825cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
8835cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
884ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
88522db5ddcSStefano Zampini         for (j=0;j<size_schur;j++) {
88622db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
88722db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
88822db5ddcSStefano Zampini #else
88922db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
89022db5ddcSStefano Zampini #endif
89122db5ddcSStefano Zampini         }
892ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
893ca92afb2SStefano Zampini       }
8945cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
8955cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
896ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
897ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
898ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
899ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
900ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
901ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
902ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
903ca92afb2SStefano Zampini     }
9043301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
9053301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9063301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
907df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
908d4933d67SStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
909df4d28bfSStefano Zampini #else
910df4d28bfSStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr);
911df4d28bfSStefano Zampini #endif
912e62b6521Sstefano_zampini     ierr = PetscOptionsGetString(NULL,((PetscObject)A)->prefix,"-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
913883469d8SStefano Zampini 
914683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
915d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
916d47842beSStefano Zampini 
917683d3df6SStefano Zampini     if (n_I) { /* TODO add ordering from options */
9185a05ddb0SStefano Zampini       IS is_schur;
9195a05ddb0SStefano Zampini 
9209ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
921d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
922883469d8SStefano Zampini       } else {
923d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
924883469d8SStefano Zampini       }
9258c09ecd8Sstefano_zampini       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
9268c09ecd8Sstefano_zampini 
927883469d8SStefano Zampini       /* subsets ordered last */
9285a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9295a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9305a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
931883469d8SStefano Zampini 
932883469d8SStefano Zampini       /* factorization step */
9339ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
934883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
935be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
936be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
937be83ff47SStefano Zampini #endif
938883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
939a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
940883469d8SStefano Zampini       } else {
941883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
942be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
943be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
944be83ff47SStefano Zampini #endif
945883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
946883469d8SStefano Zampini       }
947ec968b6dSstefano_zampini       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
948883469d8SStefano Zampini 
949883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
9507c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
951*b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
952*b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
953*b3cb21ddSStefano Zampini 
954d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
955683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
956683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
957df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
958ca92afb2SStefano Zampini 
95972b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
960ca92afb2SStefano Zampini       if (benign_n) {
9615cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
962ca92afb2SStefano Zampini         Vec         v;
9635cbda25cSStefano Zampini         PetscInt    sizeA;
964ca92afb2SStefano Zampini 
965ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
9665cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
9675cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
968ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
969ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
970ca92afb2SStefano Zampini #endif
971ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
972ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
973ca92afb2SStefano Zampini #endif
9745cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9755cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
976ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
9775cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
97847484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
979ca92afb2SStefano Zampini           const PetscInt *idxs;
980ca92afb2SStefano Zampini           PetscInt       j,nz;
98147484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
982ca92afb2SStefano Zampini 
9835cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
9845cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
9855cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
986ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
987ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
988282d6408SStefano Zampini           /* p0 dof (eliminated) is excluded from the sum */
989282d6408SStefano Zampini           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
990ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9915cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
992ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
99347484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
99447484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
99547484b83SStefano Zampini           B_k = 1;
99647484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
99747484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
9985cbda25cSStefano Zampini           sum  = 1.;
99947484b83SStefano Zampini           PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1000ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10015cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10025cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1003282d6408SStefano Zampini           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
10045cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10055cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1006ca92afb2SStefano Zampini         }
10075cbda25cSStefano Zampini   /* restore defaults */
10085cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
10095cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
10105cbda25cSStefano Zampini #endif
10115cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
10125cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
10135cbda25cSStefano Zampini #endif
10145cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10155cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1016ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1017ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1018ca92afb2SStefano Zampini       }
1019a3df083aSStefano Zampini       if (!reuse_solvers) {
1020a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1021a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1022a3df083aSStefano Zampini         }
1023a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
10245cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
10255cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1026a3df083aSStefano Zampini       }
1027df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1028be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1029683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1030166598c1SStefano Zampini       factor_workaround = PETSC_FALSE;
1031df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1032be83ff47SStefano Zampini     }
1033be83ff47SStefano Zampini 
1034be83ff47SStefano Zampini     if (reuse_solvers) {
1035a00504b5SStefano Zampini       Mat                A_II,Afake;
103653892102SStefano Zampini       Vec                vec1_B;
1037df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
10383462e049SStefano Zampini       PetscInt           n_R;
1039d5574798SStefano Zampini 
1040df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1041df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1042e28d306cSStefano Zampini       } else {
1043df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1044d62866d3SStefano Zampini       }
1045df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1046be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1047d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1048d5574798SStefano Zampini       msolv_ctx->F = F;
1049683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1050683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1051683d3df6SStefano Zampini       {
1052683d3df6SStefano Zampini         PetscScalar *array;
1053683d3df6SStefano Zampini         PetscInt    n;
1054683d3df6SStefano Zampini 
1055683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1056683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1057683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1058683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1059683d3df6SStefano Zampini       }
1060683d3df6SStefano Zampini       msolv_ctx->has_vertices = factor_workaround;
1061d62866d3SStefano Zampini 
1062d62866d3SStefano Zampini       /* interior solver */
1063683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1064d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1065d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1066d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1067df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1068df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1069d62866d3SStefano Zampini 
1070d62866d3SStefano Zampini       /* correction solver */
10713462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1072d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1073d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1074df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1075df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
107653892102SStefano Zampini 
107753892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
107853892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
107953892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1080683d3df6SStefano Zampini       if (!factor_workaround) {
108153892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
108253892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
108353892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
108453892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1085683d3df6SStefano Zampini       } else {
1086683d3df6SStefano Zampini         IS              is_B_all;
1087683d3df6SStefano Zampini         const PetscInt* idxs;
1088683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1089683d3df6SStefano Zampini 
1090683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1091683d3df6SStefano Zampini         dual = size_schur - n_v;
1092683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1093683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1094683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1095683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1096683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1097683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1098683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1099683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1100683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1101683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1102683d3df6SStefano Zampini       }
11033462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
11043462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
11053462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
11063462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1107683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1108ca92afb2SStefano Zampini 
1109ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1110ca92afb2SStefano Zampini       if (benign_n) {
11115cbda25cSStefano Zampini         PetscScalar *array;
11125cbda25cSStefano Zampini 
1113ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1114ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1115ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
11165cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
1117282d6408SStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
11185cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11195cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
11205cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11215cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1122ca92afb2SStefano Zampini       }
1123d5574798SStefano Zampini     }
112408122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
112553892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
11265db18549SStefano Zampini 
1127be83ff47SStefano Zampini     /* Work arrays */
1128be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
11295a94c880SStefano Zampini 
11305a94c880SStefano Zampini     /* matrices for adaptive selection */
113112d906b1SStefano Zampini     if (compute_Stilda) {
11325a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
11335a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11345a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
11355a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
11365a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
11375a94c880SStefano Zampini       }
11389d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
11395a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
11405a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
11415a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
11425a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
11435a94c880SStefano Zampini       }
114412d906b1SStefano Zampini     }
1145d2627357SStefano Zampini 
1146be83ff47SStefano Zampini     /* S_Ej_all */
1147be83ff47SStefano Zampini     cum = cum2 = 0;
1148be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
114965d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
115065d8bf0aSStefano Zampini       PetscInt j;
115165d8bf0aSStefano Zampini 
11525a95e1ceSStefano Zampini       /* get S_E */
1153b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1154683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1155be83ff47SStefano Zampini         PetscInt k;
1156be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1157be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1158be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1159683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1160be83ff47SStefano Zampini           }
1161be83ff47SStefano Zampini         }
116206a4e24aSStefano Zampini       } else { /* just copy to workspace */
1163be83ff47SStefano Zampini         PetscInt k;
1164be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1165be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1166be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1167be83ff47SStefano Zampini           }
1168be83ff47SStefano Zampini         }
11699087bf02SStefano Zampini       }
11705a95e1ceSStefano Zampini       /* insert S_E values */
1171be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1172b7ab4a40SStefano Zampini       if (sub_schurs->change) {
11738760537fSStefano Zampini         Mat            change_sub,SEj,T;
117472b8c272SStefano Zampini 
117572b8c272SStefano Zampini         /* change basis */
117672b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
117772b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
11788760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
11798760537fSStefano Zampini           Mat T2;
11808760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
11818760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
11828760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
118372b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
11848760537fSStefano Zampini         } else {
11858760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
118672b8c272SStefano Zampini         }
118772b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
118872b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
1189afa0f562SStefano Zampini         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
119072b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
119172b8c272SStefano Zampini       }
11929d54b7f4SStefano Zampini       if (deluxe) {
11935a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1194683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1195862806e4SStefano Zampini         if (compute_Stilda) {
119608122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
119708122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
119806a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
11996c3e6151SStefano Zampini             PetscInt k;
12005a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12012972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
12025a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
12032972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
12046c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
12056c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
12066c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
12076c3e6151SStefano Zampini               }
12086c3e6151SStefano Zampini             }
1209d6462365SStefano Zampini           } else {
1210d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1211d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1212d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1213d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
12142972d61bSStefano Zampini           }
121508122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
12165a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12179087bf02SStefano Zampini         }
12189d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
12199d54b7f4SStefano Zampini         Mat         SEj;
12209d54b7f4SStefano Zampini         Vec         D;
12219d54b7f4SStefano Zampini         PetscScalar *array;
12229d54b7f4SStefano Zampini 
12239d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12249d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
12259d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
12269d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1227f17d2ae1SStefano Zampini         ierr = VecShift(D,-1.);CHKERRQ(ierr);
12289d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
12299d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
12309d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
12319d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12329d54b7f4SStefano Zampini       }
1233be83ff47SStefano Zampini       cum += subset_size;
1234be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1235be83ff47SStefano Zampini     }
1236be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
12374a6c6b0dSStefano Zampini 
1238df4d28bfSStefano Zampini     if (solver_S) {
12397c2f51b8SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
12404a6c6b0dSStefano Zampini     }
1241683d3df6SStefano Zampini 
1242683d3df6SStefano Zampini     schur_factor = NULL;
124345951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1244683d3df6SStefano Zampini 
12459d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
12464a6c6b0dSStefano Zampini         PetscInt j;
12474a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
12485a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12494a6c6b0dSStefano Zampini       } else {
1250683d3df6SStefano Zampini         Mat S_all_inv=NULL;
1251df4d28bfSStefano Zampini         if (solver_S) { /* use MatFactor calls to invert S */
125272b8c272SStefano Zampini             /* let MatFactor factor the Schur complement */
12536dba178dSStefano Zampini           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
1254683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1255683d3df6SStefano Zampini              The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */
1256683d3df6SStefano Zampini           if (factor_workaround) {
1257683d3df6SStefano Zampini             PetscScalar *data;
1258683d3df6SStefano Zampini             PetscInt nd = 0;
12596dba178dSStefano Zampini 
1260683d3df6SStefano Zampini             if (!sub_schurs->is_posdef) {
1261683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1262683d3df6SStefano Zampini             }
12637c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1264683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1265683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1266683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1267683d3df6SStefano Zampini             }
1268683d3df6SStefano Zampini             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
1269683d3df6SStefano Zampini             cum = 0;
1270683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1271683d3df6SStefano Zampini               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1272683d3df6SStefano Zampini               cum += size_active_schur-i;
1273683d3df6SStefano Zampini             }
1274683d3df6SStefano Zampini             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
12756dba178dSStefano Zampini             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1276683d3df6SStefano Zampini             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1277683d3df6SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1278683d3df6SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
1279683d3df6SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1280683d3df6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1281683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1282683d3df6SStefano Zampini           } else {
12836dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
12847c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1285683d3df6SStefano Zampini           }
1286df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1287683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1288683d3df6SStefano Zampini           S_all_inv = S_all;
1289683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1290be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1291be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
129206a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1293be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1294be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1295be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1296be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1297d6462365SStefano Zampini           } else {
1298d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1299d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1300d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1301d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1302be83ff47SStefano Zampini           }
1303be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1304683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1305be83ff47SStefano Zampini         }
1306be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1307be83ff47SStefano Zampini         cum = cum2 = 0;
1308683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1309be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1310be83ff47SStefano Zampini           PetscInt j;
1311862806e4SStefano Zampini 
1312862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1313be83ff47SStefano Zampini           /* get (St^-1)_E */
131472b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
131506a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1316a0b0af32SStefano Zampini           if (S_lower_triangular) {
1317be83ff47SStefano Zampini             PetscInt k;
1318b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1319be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1320be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1321be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
13226c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1323be83ff47SStefano Zampini                 }
1324be83ff47SStefano Zampini               }
132572b8c272SStefano Zampini             } else {
132672b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
132772b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
132872b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
132972b8c272SStefano Zampini                 }
133072b8c272SStefano Zampini               }
133172b8c272SStefano Zampini             }
133272b8c272SStefano Zampini           } else {
1333be83ff47SStefano Zampini             PetscInt k;
1334be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1335be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1336be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1337be83ff47SStefano Zampini               }
1338be83ff47SStefano Zampini             }
1339be83ff47SStefano Zampini           }
1340b7ab4a40SStefano Zampini           if (sub_schurs->change) {
13418760537fSStefano Zampini             Mat change_sub,SEj,T;
134272b8c272SStefano Zampini 
134372b8c272SStefano Zampini             /* change basis */
134472b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
134572b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
13468760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
13478760537fSStefano Zampini               Mat T2;
13488760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
13498760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
135072b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
13518760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
13528760537fSStefano Zampini             } else {
13538760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
135472b8c272SStefano Zampini             }
13558760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
13568760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
135772b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1358afa0f562SStefano Zampini             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
135972b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
136072b8c272SStefano Zampini           }
1361be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
13625a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1363be83ff47SStefano Zampini           cum += subset_size;
1364be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1365883469d8SStefano Zampini         }
1366683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1367df4d28bfSStefano Zampini         if (solver_S) {
13687c2f51b8SStefano Zampini           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
13695db18549SStefano Zampini         }
1370683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1371683d3df6SStefano Zampini       }
1372683d3df6SStefano Zampini 
1373683d3df6SStefano Zampini       /* move back factors to Schur data into F */
1374683d3df6SStefano Zampini       if (factor_workaround) {
1375683d3df6SStefano Zampini         Mat S_tmp;
1376683d3df6SStefano Zampini         PetscInt nv,nd=0;
1377683d3df6SStefano Zampini 
1378df4d28bfSStefano Zampini         if (!solver_S) {
1379683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1380683d3df6SStefano Zampini         }
1381683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
13827c2f51b8SStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1383683d3df6SStefano Zampini         if (sub_schurs->is_posdef) {
1384683d3df6SStefano Zampini           PetscScalar *data;
1385683d3df6SStefano Zampini 
1386683d3df6SStefano Zampini           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1387683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1388683d3df6SStefano Zampini 
1389683d3df6SStefano Zampini           if (S_lower_triangular) {
1390683d3df6SStefano Zampini             cum = 0;
1391683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1392683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1393683d3df6SStefano Zampini               cum += size_active_schur-i;
1394683d3df6SStefano Zampini 	    }
1395683d3df6SStefano Zampini           } else {
1396683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1397683d3df6SStefano Zampini           }
1398683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1399683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1400683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1401683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1402683d3df6SStefano Zampini 	    }
1403683d3df6SStefano Zampini           }
14046dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1405683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1406683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
14076c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1408683d3df6SStefano Zampini           }
1409683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1410683d3df6SStefano Zampini         } else {
1411683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1412683d3df6SStefano Zampini         }
14137c2f51b8SStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
14149087bf02SStefano Zampini       }
1415367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1416367aa537SStefano Zampini       PetscScalar *data;
1417367aa537SStefano Zampini       PetscInt    nd = 0;
1418367aa537SStefano Zampini 
1419367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1420367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1421367aa537SStefano Zampini       }
14227c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1423367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1424367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1425367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1426367aa537SStefano Zampini       }
1427367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1428367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
14296c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1430367aa537SStefano Zampini       }
1431367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
14327c2f51b8SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
14334a6c6b0dSStefano Zampini     }
14344a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1435683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
14369d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
14374a6c6b0dSStefano Zampini   }
14385a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1439be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1440a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1441a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1442a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1443a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1444a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
14455db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14465db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14475a95e1ceSStefano Zampini   if (compute_Stilda) {
14485a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14495a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14509d54b7f4SStefano Zampini     if (deluxe) {
14515a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14525a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145308122e43SStefano Zampini     }
14549d54b7f4SStefano Zampini   }
1455a1337663SStefano Zampini 
14565db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
14575db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
14583927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
14593927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
14605a95e1ceSStefano Zampini 
14615db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
14625a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
14637dae84e0SHong Zhang     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
14645a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
14655a94c880SStefano Zampini   } else {
14665a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
14675a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
14687dae84e0SHong Zhang     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
14695a94c880SStefano Zampini   }
147008122e43SStefano Zampini 
1471f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
14725a95e1ceSStefano Zampini   if (compute_Stilda) {
14735a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1474a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
147553153037Sstefano_zampini     submats[0] = NULL;
14767dae84e0SHong Zhang     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
147753153037Sstefano_zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
147853153037Sstefano_zampini     sub_schurs->sum_S_Ej_tilda_all = submats[0];
14799d54b7f4SStefano Zampini     if (deluxe) {
14805a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
148108122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
148253153037Sstefano_zampini       submats[0] = NULL;
14837dae84e0SHong Zhang       ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
148453153037Sstefano_zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
148553153037Sstefano_zampini       sub_schurs->sum_S_Ej_inv_all = submats[0];
14869d54b7f4SStefano Zampini     } else {
14879d54b7f4SStefano Zampini       PetscScalar *array;
14889d54b7f4SStefano Zampini       PetscInt    cum;
14899d54b7f4SStefano Zampini 
14909d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
14919d54b7f4SStefano Zampini       cum = 0;
14929d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
14939d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
14949d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
14959d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
14969d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
14979d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
14989d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
14999d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
15009d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
15019d54b7f4SStefano Zampini         cum += subset_size*subset_size;
15029d54b7f4SStefano Zampini       }
15039d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15049d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
15059d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
15069d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
15079d54b7f4SStefano Zampini     }
150808122e43SStefano Zampini   }
15093202ece2SStefano Zampini 
15105a95e1ceSStefano Zampini   /* free workspace */
15115a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
151206a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1513a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
15143202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1515dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
15165a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1517b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1518b1b3d7a2SStefano Zampini }
1519b1b3d7a2SStefano Zampini 
15208b6046baSStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1521b1b3d7a2SStefano Zampini {
15229bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
15235a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1524b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1525b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1526b1b3d7a2SStefano Zampini 
1527b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1528b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
15296c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1530b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
15316c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1532b1b3d7a2SStefano Zampini 
1533b1b3d7a2SStefano Zampini   /* reset any previous data */
1534b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1535b1b3d7a2SStefano Zampini 
15365a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
15379bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1538b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
153908122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1540b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1541b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
15428b6046baSStefano Zampini     if (copycc) {
15438b6046baSStefano Zampini       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
15448b6046baSStefano Zampini     } else {
1545c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1546b1b3d7a2SStefano Zampini       all_cc[i] = faces[i];
1547b1b3d7a2SStefano Zampini     }
15488b6046baSStefano Zampini   }
1549b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
15508b6046baSStefano Zampini     if (copycc) {
15518b6046baSStefano Zampini       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
15528b6046baSStefano Zampini     } else {
1553c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1554b1b3d7a2SStefano Zampini       all_cc[n_faces+i] = edges[i];
15558b6046baSStefano Zampini     }
155608122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1557b1b3d7a2SStefano Zampini   }
1558c8272957SStefano Zampini   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1559c8272957SStefano Zampini   sub_schurs->is_vertices = vertices;
1560c8272957SStefano Zampini   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1561d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1562d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1563b1b3d7a2SStefano Zampini 
1564df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
1565df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_FALSE;
1566883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1567df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1568df4d28bfSStefano Zampini #endif
1569df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1570df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1571883469d8SStefano Zampini #endif
1572b1b3d7a2SStefano Zampini 
1573b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1574b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1575b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1576b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
15775db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
15785db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
15795db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
15805db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
15815a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1582b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1583df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1584b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1585b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1586b96c3477SStefano Zampini     }
15879bb4a8caSStefano Zampini   }
1588b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1589b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
159008122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1591b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1592b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1593b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1594b1b3d7a2SStefano Zampini }
1595b1b3d7a2SStefano Zampini 
159634a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
159734a97f8cSStefano Zampini {
159834a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
159934a97f8cSStefano Zampini   PetscErrorCode  ierr;
160034a97f8cSStefano Zampini 
160134a97f8cSStefano Zampini   PetscFunctionBegin;
160234a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
16035ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
160434a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
160534a97f8cSStefano Zampini   PetscFunctionReturn(0);
160634a97f8cSStefano Zampini }
160734a97f8cSStefano Zampini 
160834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
160934a97f8cSStefano Zampini {
161034a97f8cSStefano Zampini   PetscInt       i;
161134a97f8cSStefano Zampini   PetscErrorCode ierr;
161234a97f8cSStefano Zampini 
161334a97f8cSStefano Zampini   PetscFunctionBegin;
1614aea80f77Sstefano_zampini   if (!sub_schurs) PetscFunctionReturn(0);
16151e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1616b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1617b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1618b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
16195db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
16205db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
162141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
162241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
162308122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1624a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
16255db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1626d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1627d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
162808122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
162934a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1630b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
163134a97f8cSStefano Zampini   }
16325ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1633b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
16343dc780c3SStefano Zampini   }
1635df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1636df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1637d62866d3SStefano Zampini   }
1638df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
163972b8c272SStefano Zampini   if (sub_schurs->change) {
164072b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
164172b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1642b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
164372b8c272SStefano Zampini     }
164472b8c272SStefano Zampini   }
164572b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1646b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
164734a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
164834a97f8cSStefano Zampini   PetscFunctionReturn(0);
164934a97f8cSStefano Zampini }
165034a97f8cSStefano Zampini 
1651aea80f77Sstefano_zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1652aea80f77Sstefano_zampini {
1653aea80f77Sstefano_zampini   PetscErrorCode ierr;
1654aea80f77Sstefano_zampini 
1655aea80f77Sstefano_zampini   PetscFunctionBegin;
1656aea80f77Sstefano_zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1657aea80f77Sstefano_zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1658aea80f77Sstefano_zampini   PetscFunctionReturn(0);
1659aea80f77Sstefano_zampini }
1660aea80f77Sstefano_zampini 
16612a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
166234a97f8cSStefano Zampini {
166334a97f8cSStefano Zampini   PetscInt       i,j,n;
166434a97f8cSStefano Zampini   PetscErrorCode ierr;
166534a97f8cSStefano Zampini 
166634a97f8cSStefano Zampini   PetscFunctionBegin;
166734a97f8cSStefano Zampini   n = 0;
166834a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
166934a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
167034a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
167134a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
167234a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
167334a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
167434a97f8cSStefano Zampini         queue_tip[n] = dof;
167534a97f8cSStefano Zampini         n++;
167634a97f8cSStefano Zampini       }
167734a97f8cSStefano Zampini     }
167834a97f8cSStefano Zampini   }
167934a97f8cSStefano Zampini   *n_added = n;
168034a97f8cSStefano Zampini   PetscFunctionReturn(0);
168134a97f8cSStefano Zampini }
1682