xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision d943a64230c34c465b2a9e7d554d914b82517bfe)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
305709791SSatish Balay #include <../src/mat/impls/dense/seq/dense.h>
408122e43SStefano Zampini #include <petscblaslapack.h>
534a97f8cSStefano Zampini 
63202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
75ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
8df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
9df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
10d62866d3SStefano Zampini 
11ca92afb2SStefano Zampini /* if v2 is not present, correction is done in-place */
125cbda25cSStefano Zampini PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
13ca92afb2SStefano Zampini {
14ca92afb2SStefano Zampini   PetscScalar    *array;
15ca92afb2SStefano Zampini   PetscScalar    *array2;
16ca92afb2SStefano Zampini   PetscErrorCode ierr;
17ca92afb2SStefano Zampini 
18ca92afb2SStefano Zampini   PetscFunctionBegin;
19ca92afb2SStefano Zampini   if (!ctx->benign_n) PetscFunctionReturn(0);
205cbda25cSStefano Zampini   if (sol && full) {
215cbda25cSStefano Zampini     PetscInt n_I,size_schur;
225cbda25cSStefano Zampini 
235cbda25cSStefano Zampini     /* get sizes */
24282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
255cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
265cbda25cSStefano Zampini     n_I = n_I - size_schur;
275cbda25cSStefano Zampini     /* get schur sol from array */
285cbda25cSStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
295cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
305cbda25cSStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
315cbda25cSStefano Zampini     /* apply interior sol correction */
32282d6408SStefano Zampini     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
335cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
345cbda25cSStefano Zampini     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
355cbda25cSStefano Zampini   }
36ca92afb2SStefano Zampini   if (v2) {
37ca92afb2SStefano Zampini     PetscInt nl;
38ca92afb2SStefano Zampini 
39ca92afb2SStefano Zampini     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
40ca92afb2SStefano Zampini     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
41ca92afb2SStefano Zampini     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
42ca92afb2SStefano Zampini     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
43ca92afb2SStefano Zampini   } else {
44ca92afb2SStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
45ca92afb2SStefano Zampini     array2 = array;
46ca92afb2SStefano Zampini   }
47ca92afb2SStefano Zampini   if (!sol) { /* change rhs */
48ca92afb2SStefano Zampini     PetscInt n;
49ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
50ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
51ca92afb2SStefano Zampini       const PetscInt *cols;
52ca92afb2SStefano Zampini       PetscInt       nz,i;
53ca92afb2SStefano Zampini 
54ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
55ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
56ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
5722db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
5822db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
5922db5ddcSStefano Zampini #else
60ca92afb2SStefano Zampini       sum = -sum/nz;
6122db5ddcSStefano Zampini #endif
62ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
63ca92afb2SStefano Zampini       ctx->benign_save_vals[n] = array2[cols[nz-1]];
64ca92afb2SStefano Zampini       array2[cols[nz-1]] = sum;
65ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
66ca92afb2SStefano Zampini     }
67ca92afb2SStefano Zampini   } else {
68ca92afb2SStefano Zampini     PetscInt n;
69ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
70ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
71ca92afb2SStefano Zampini       const PetscInt *cols;
72ca92afb2SStefano Zampini       PetscInt       nz,i;
73ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
74ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
75ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
7622db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
7722db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
7822db5ddcSStefano Zampini #else
79ca92afb2SStefano Zampini       sum = -sum/nz;
8022db5ddcSStefano Zampini #endif
81ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
82ca92afb2SStefano Zampini       array2[cols[nz-1]] = ctx->benign_save_vals[n];
83ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
84ca92afb2SStefano Zampini     }
85ca92afb2SStefano Zampini   }
86ca92afb2SStefano Zampini   if (v2) {
87ca92afb2SStefano Zampini     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
8872b8c272SStefano Zampini     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
89ca92afb2SStefano Zampini   } else {
90ca92afb2SStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
91ca92afb2SStefano Zampini   }
925cbda25cSStefano Zampini   if (!sol && full) {
935cbda25cSStefano Zampini     Vec      usedv;
945cbda25cSStefano Zampini     PetscInt n_I,size_schur;
955cbda25cSStefano Zampini 
965cbda25cSStefano Zampini     /* get sizes */
97282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
985cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
995cbda25cSStefano Zampini     n_I = n_I - size_schur;
1005cbda25cSStefano Zampini     /* compute schur rhs correction */
1015cbda25cSStefano Zampini     if (v2) {
1025cbda25cSStefano Zampini       usedv = v2;
1035cbda25cSStefano Zampini     } else {
1045cbda25cSStefano Zampini       usedv = v;
1055cbda25cSStefano Zampini     }
1065cbda25cSStefano Zampini     /* apply schur rhs correction */
1075cbda25cSStefano Zampini     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
1085cbda25cSStefano Zampini     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1095cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
1105cbda25cSStefano Zampini     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
111282d6408SStefano Zampini     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1125cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1135cbda25cSStefano Zampini   }
114ca92afb2SStefano Zampini   PetscFunctionReturn(0);
115ca92afb2SStefano Zampini }
116ca92afb2SStefano Zampini 
117df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
118d62866d3SStefano Zampini {
119df4d28bfSStefano Zampini   PCBDDCReuseSolvers ctx;
120683d3df6SStefano Zampini   PetscBool          copy = PETSC_FALSE;
121d62866d3SStefano Zampini   PetscErrorCode     ierr;
122d62866d3SStefano Zampini 
123d62866d3SStefano Zampini   PetscFunctionBegin;
124d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
125683d3df6SStefano Zampini   if (full) {
126d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
127d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
128d62866d3SStefano Zampini #endif
1295cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1305cbda25cSStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
1315cbda25cSStefano Zampini #endif
132683d3df6SStefano Zampini     copy = ctx->has_vertices;
133d4933d67SStefano Zampini   } else { /* interior solver */
1346dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
135683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
1366dba178dSStefano Zampini #endif
137d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
138d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
139d4933d67SStefano Zampini #endif
140683d3df6SStefano Zampini     copy = PETSC_TRUE;
141683d3df6SStefano Zampini   }
142683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
143683d3df6SStefano Zampini   if (copy) {
144ca92afb2SStefano Zampini     PetscInt    n;
145df4d28bfSStefano Zampini     PetscScalar *array,*array_solver;
146ca92afb2SStefano Zampini 
147ca92afb2SStefano Zampini     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
148683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
149df4d28bfSStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
150df4d28bfSStefano Zampini     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
151df4d28bfSStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
152683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
153683d3df6SStefano Zampini 
1545cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
155683d3df6SStefano Zampini     if (transpose) {
156683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
157683d3df6SStefano Zampini     } else {
158683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
159683d3df6SStefano Zampini     }
1605cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
161683d3df6SStefano Zampini 
162683d3df6SStefano Zampini     /* get back data to caller worskpace */
163df4d28bfSStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
164683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
165df4d28bfSStefano Zampini     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
166683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
167df4d28bfSStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
168683d3df6SStefano Zampini   } else {
169ca92afb2SStefano Zampini     if (ctx->benign_n) {
1705cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
171ca92afb2SStefano Zampini       if (transpose) {
172ca92afb2SStefano Zampini         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
173ca92afb2SStefano Zampini       } else {
174ca92afb2SStefano Zampini         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
175ca92afb2SStefano Zampini       }
1765cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
177ca92afb2SStefano Zampini     } else {
178e28d306cSStefano Zampini       if (transpose) {
179e28d306cSStefano Zampini         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
180e28d306cSStefano Zampini       } else {
1816816873aSStefano Zampini         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
182e28d306cSStefano Zampini       }
183683d3df6SStefano Zampini     }
184ca92afb2SStefano Zampini   }
1855cbda25cSStefano Zampini   /* restore defaults */
1865cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1875cbda25cSStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
1885cbda25cSStefano Zampini #endif
189d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
190d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
191d4933d67SStefano Zampini #endif
192d62866d3SStefano Zampini   PetscFunctionReturn(0);
193d62866d3SStefano Zampini }
194d62866d3SStefano Zampini 
195df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
196e28d306cSStefano Zampini {
197e28d306cSStefano Zampini   PetscErrorCode   ierr;
198e28d306cSStefano Zampini 
199e28d306cSStefano Zampini   PetscFunctionBegin;
200df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
201e28d306cSStefano Zampini   PetscFunctionReturn(0);
202e28d306cSStefano Zampini }
203e28d306cSStefano Zampini 
204df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
205e28d306cSStefano Zampini {
206e28d306cSStefano Zampini   PetscErrorCode   ierr;
207e28d306cSStefano Zampini 
208e28d306cSStefano Zampini   PetscFunctionBegin;
209df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
210683d3df6SStefano Zampini   PetscFunctionReturn(0);
211683d3df6SStefano Zampini }
212683d3df6SStefano Zampini 
213df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
214683d3df6SStefano Zampini {
215683d3df6SStefano Zampini   PetscErrorCode   ierr;
216683d3df6SStefano Zampini 
217683d3df6SStefano Zampini   PetscFunctionBegin;
218df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
219683d3df6SStefano Zampini   PetscFunctionReturn(0);
220683d3df6SStefano Zampini }
221683d3df6SStefano Zampini 
222df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
223683d3df6SStefano Zampini {
224683d3df6SStefano Zampini   PetscErrorCode   ierr;
225683d3df6SStefano Zampini 
226683d3df6SStefano Zampini   PetscFunctionBegin;
227df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
228e28d306cSStefano Zampini   PetscFunctionReturn(0);
229e28d306cSStefano Zampini }
230e28d306cSStefano Zampini 
231df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
232d62866d3SStefano Zampini {
233ca92afb2SStefano Zampini   PetscInt       i;
234d62866d3SStefano Zampini   PetscErrorCode ierr;
235d62866d3SStefano Zampini 
236d62866d3SStefano Zampini   PetscFunctionBegin;
237d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
238d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
239d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
240d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
241d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
24253892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
24353892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
244d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
24553892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
24653892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
247ca92afb2SStefano Zampini   for (i=0;i<reuse->benign_n;i++) {
248ca92afb2SStefano Zampini     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
249ca92afb2SStefano Zampini   }
250ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
251ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
2525cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
2535cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
2545cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
2555cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
256d62866d3SStefano Zampini   PetscFunctionReturn(0);
257d62866d3SStefano Zampini }
258d5574798SStefano Zampini 
2595ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
2603202ece2SStefano Zampini {
2613202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
2623202ece2SStefano Zampini   KSP            ksp;
2633202ece2SStefano Zampini   PC             pc;
2643202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
2653202ece2SStefano Zampini   PetscReal      fill = 2.0;
266f11841e3SStefano Zampini   PetscInt       n_I;
2673202ece2SStefano Zampini   PetscMPIInt    size;
2683202ece2SStefano Zampini   PetscErrorCode ierr;
2693202ece2SStefano Zampini 
2703202ece2SStefano Zampini   PetscFunctionBegin;
2713202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
2726c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
273f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
274f11841e3SStefano Zampini     PetscBool Sdense;
275f11841e3SStefano Zampini 
276f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
2776c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
278f11841e3SStefano Zampini   }
2793202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
2803202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
2813202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
2823202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
2833202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
2843202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
2853202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
2863202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
287f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
288f11841e3SStefano Zampini   if (n_I) {
2893202ece2SStefano Zampini     if (!Bdense) {
2903202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
2913202ece2SStefano Zampini     } else {
2923202ece2SStefano Zampini       Bd = B;
2933202ece2SStefano Zampini     }
2943202ece2SStefano Zampini 
2953202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
2963202ece2SStefano Zampini       Mat fact;
2973202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
2983202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
2993202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
3003202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
3013202ece2SStefano Zampini     } else {
30207b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
30307b1e237SStefano Zampini 
30407b1e237SStefano Zampini       if (ex) {
3053202ece2SStefano Zampini         Mat Ainvd;
3063202ece2SStefano Zampini 
3073202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
3083202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
3093202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
31007b1e237SStefano Zampini       } else {
31107b1e237SStefano Zampini         Vec         sol,rhs;
31207b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
31307b1e237SStefano Zampini         PetscInt    i,nrhs,n;
31407b1e237SStefano Zampini 
31507b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
31607b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
31707b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
31807b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
31907b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
32007b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
32107b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
32207b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
32307b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
32407b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
32507b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
32607b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
32707b1e237SStefano Zampini         }
32807b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
32907b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
33007b1e237SStefano Zampini       }
3313202ece2SStefano Zampini     }
3325ec10c6aSStefano Zampini     if (!Bdense & !issym) {
3333202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3343202ece2SStefano Zampini     }
3355ec10c6aSStefano Zampini 
3365ec10c6aSStefano Zampini     if (!issym) {
3373202ece2SStefano Zampini       if (!Cdense) {
3383202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
3393202ece2SStefano Zampini       } else {
3403202ece2SStefano Zampini         Cd = C;
3413202ece2SStefano Zampini       }
3425ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3433202ece2SStefano Zampini       if (!Cdense) {
3443202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
3453202ece2SStefano Zampini       }
3465ec10c6aSStefano Zampini     } else {
3475ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3485ec10c6aSStefano Zampini       if (!Bdense) {
3495ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3505ec10c6aSStefano Zampini       }
3515ec10c6aSStefano Zampini     }
3525ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
353f11841e3SStefano Zampini   }
3543202ece2SStefano Zampini 
3553202ece2SStefano Zampini   if (D) {
3563202ece2SStefano Zampini     Mat       Dd;
3573202ece2SStefano Zampini     PetscBool Ddense;
3583202ece2SStefano Zampini 
3593202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
3603202ece2SStefano Zampini     if (!Ddense) {
3613202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
3623202ece2SStefano Zampini     } else {
3633202ece2SStefano Zampini       Dd = D;
3643202ece2SStefano Zampini     }
365f11841e3SStefano Zampini     if (n_I) {
3663202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
367f11841e3SStefano Zampini     } else {
368f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
369f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
370f11841e3SStefano Zampini       } else {
371f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
372f11841e3SStefano Zampini       }
373f11841e3SStefano Zampini     }
3743202ece2SStefano Zampini     if (!Ddense) {
3753202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
3763202ece2SStefano Zampini     }
3773202ece2SStefano Zampini   } else {
3783202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
3793202ece2SStefano Zampini   }
3803202ece2SStefano Zampini   PetscFunctionReturn(0);
3813202ece2SStefano Zampini }
38234a97f8cSStefano Zampini 
38391af6908SStefano 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)
384b1b3d7a2SStefano Zampini {
385be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
386be83ff47SStefano Zampini   Mat                    S_all;
3875a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
3885db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
389b7ab4a40SStefano Zampini   IS                     is_I,is_I_layer;
390dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
391dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
392dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
3935a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
394683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
39508122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
39606a4b1faSStefano Zampini   PetscScalar            *Bwork;
3975a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
3985a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
3995a95e1ceSStefano Zampini   MPI_Comm               comm_n;
400f4f7d9d6SStefano Zampini   PetscBool              deluxe = PETSC_TRUE;
401f4f7d9d6SStefano Zampini   PetscBool              use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
402b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
403b1b3d7a2SStefano Zampini 
404b1b3d7a2SStefano Zampini   PetscFunctionBegin;
405a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
406a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
407e62b6521Sstefano_zampini   /* convert matrix if needed */
408e62b6521Sstefano_zampini   if (Ain) {
409e62b6521Sstefano_zampini     PetscBool isseqaij;
410e62b6521Sstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
411a64f4aa4SStefano Zampini     if (isseqaij) {
412a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
413a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4143301b35fSStefano Zampini     } else {
415a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
416a64f4aa4SStefano Zampini     }
417a64f4aa4SStefano Zampini   }
4183301b35fSStefano Zampini 
419a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
420a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
421df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
422df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
423a64f4aa4SStefano Zampini   }
424a64f4aa4SStefano Zampini 
4255a95e1ceSStefano Zampini   /* preliminary checks */
426af25d912SStefano 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");
4275a95e1ceSStefano Zampini 
42888113c35SStefano Zampini   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
42988113c35SStefano Zampini 
4305a95e1ceSStefano Zampini   /* restrict work on active processes */
4315a95e1ceSStefano Zampini   color = 0;
4325a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4335a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4345a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4355a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4365a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4375a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
4385a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
4395a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4405a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
4415a95e1ceSStefano Zampini     PetscFunctionReturn(0);
4425a95e1ceSStefano Zampini   }
44381ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
4445a95e1ceSStefano Zampini 
445b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
446df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
447a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
4483301b35fSStefano Zampini     PetscBool isseqsbaij;
449a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
4503301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
4513301b35fSStefano Zampini     if (isseqsbaij) {
452a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
453a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
454a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
455a64f4aa4SStefano Zampini     } else {
456a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
457a64f4aa4SStefano Zampini       A_BB = tA_BB;
458a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
459a64f4aa4SStefano Zampini       A_IB = tA_IB;
460a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
461a64f4aa4SStefano Zampini       A_BI = tA_BI;
462f11841e3SStefano Zampini     }
463a58a30b4SStefano Zampini   } else {
4645a95e1ceSStefano Zampini     A_II = NULL;
4655a95e1ceSStefano Zampini     A_IB = NULL;
4665a95e1ceSStefano Zampini     A_BI = NULL;
4675a95e1ceSStefano Zampini     A_BB = NULL;
468b1b3d7a2SStefano Zampini   }
4695a95e1ceSStefano Zampini   S_all = NULL;
470b1b3d7a2SStefano Zampini 
471b1b3d7a2SStefano Zampini   /* determine interior problems */
4723dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4733dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
474b1b3d7a2SStefano Zampini     PetscBT                touched;
475b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
476b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
477b1b3d7a2SStefano Zampini 
4781633d1f0SStefano Zampini     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
479b1b3d7a2SStefano Zampini     /* get sizes */
480b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
481b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
482b1b3d7a2SStefano Zampini 
483b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
484b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
485b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
486b1b3d7a2SStefano Zampini 
487b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
488b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
489b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
490b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
491b1b3d7a2SStefano Zampini     }
492b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
493b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
494b1b3d7a2SStefano Zampini 
495b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
496b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
497b1b3d7a2SStefano Zampini     n_prev_added = n_B;
498b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
499b1b3d7a2SStefano Zampini       PetscInt n_added;
500b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5016c4ed002SBarry 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);
502b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
503b1b3d7a2SStefano Zampini       n_prev_added = n_added;
504b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
505b1b3d7a2SStefano Zampini       if (!n_added) break;
506b1b3d7a2SStefano Zampini     }
507b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
508b1b3d7a2SStefano Zampini 
509883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
510a9b99552SStefano 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);
511b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
512a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
513883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
514df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
515b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
516b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
517a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
518b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
519b1b3d7a2SStefano Zampini 
520b1b3d7a2SStefano Zampini       /* II block */
5217dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
522b1b3d7a2SStefano Zampini     }
523b1b3d7a2SStefano Zampini   } else {
524b1b3d7a2SStefano Zampini     PetscInt n_I;
525b1b3d7a2SStefano Zampini 
526b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
527b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
528a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
529b1b3d7a2SStefano Zampini 
530b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
531df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
532b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
533b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
534b1b3d7a2SStefano Zampini 
535b1b3d7a2SStefano Zampini       /* II block is the same */
536b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
537b1b3d7a2SStefano Zampini       AE_II = A_II;
538b1b3d7a2SStefano Zampini     }
539b1b3d7a2SStefano Zampini   }
5405a95e1ceSStefano Zampini 
541883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
5425a95e1ceSStefano Zampini   max_subset_size = 0;
543883469d8SStefano Zampini   local_size = 0;
5445a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
5455a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5465a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
547883469d8SStefano Zampini     local_size += subset_size;
548883469d8SStefano Zampini   }
549883469d8SStefano Zampini 
550883469d8SStefano Zampini   /* Work arrays for local indices */
551883469d8SStefano Zampini   extra = 0;
552683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
553df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
554a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
555883469d8SStefano Zampini   }
556683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
557883469d8SStefano Zampini   if (extra) {
558883469d8SStefano Zampini     const PetscInt *idxs;
559a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
560883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
561a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
562883469d8SStefano Zampini   }
563883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
564dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
565dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
566883469d8SStefano Zampini 
567883469d8SStefano Zampini   /* Get local indices in local numbering */
568883469d8SStefano Zampini   local_size = 0;
5695a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
570883469d8SStefano Zampini     PetscInt j;
571883469d8SStefano Zampini     const    PetscInt *idxs;
572883469d8SStefano Zampini 
5735a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5745a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
575eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
576eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
577eb595f79SStefano Zampini     auxnum2[i] = subset_size;
578883469d8SStefano Zampini     /* subset indices in local numbering */
579883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5805a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
581883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
582883469d8SStefano Zampini     local_size += subset_size;
583883469d8SStefano Zampini   }
584883469d8SStefano Zampini 
585f4f7d9d6SStefano Zampini   /* allocate extra workspace needed only for GETRI or SYTRF */
58611955456SStefano Zampini   use_potr = use_sytr = PETSC_FALSE;
58711955456SStefano Zampini   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
588f4f7d9d6SStefano Zampini     use_potr = PETSC_TRUE;
58911955456SStefano Zampini   } else if (sub_schurs->is_symmetric) {
59011955456SStefano Zampini     use_sytr = PETSC_TRUE;
59111955456SStefano Zampini   }
59211955456SStefano Zampini   if (local_size && !use_potr) {
59359ac4de7SStefano Zampini     PetscScalar  lwork,dummyscalar = 0.;
59459ac4de7SStefano Zampini     PetscBLASInt dummyint = 0;
595d2627357SStefano Zampini 
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);
599f4f7d9d6SStefano Zampini     if (use_sytr) {
600f4f7d9d6SStefano Zampini       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
601f4f7d9d6SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
602f4f7d9d6SStefano Zampini     } else {
60359ac4de7SStefano Zampini       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
604d2627357SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
605f4f7d9d6SStefano Zampini     }
606f4f7d9d6SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
607d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
608d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
609056290a2SStefano Zampini   } else {
610056290a2SStefano Zampini     Bwork = NULL;
611056290a2SStefano Zampini     pivots = NULL;
612d2627357SStefano Zampini   }
613d2627357SStefano Zampini 
614d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
615dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
616dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
617dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
618dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
6196583bcc1SStefano Zampini   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
620dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
621dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
622dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6236c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
624dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
625e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
626d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
627d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
628d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
629d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6302972d61bSStefano Zampini 
6315a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6325a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6335a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6345a95e1ceSStefano Zampini 
6355a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6365a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
6376c4ed002SBarry 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);
6385a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
639b1b3d7a2SStefano Zampini   }
640b1b3d7a2SStefano Zampini 
64172b8c272SStefano Zampini   if (change) {
64272b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
64372b8c272SStefano Zampini     IS                     change_primal_B;
64472b8c272SStefano Zampini     IS                     change_primal_all;
64572b8c272SStefano Zampini 
646b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
647b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
648b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
64972b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
65072b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
65172b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
652b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
65372b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
65472b8c272SStefano Zampini     }
65572b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
65672b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
65772b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
65872b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
65972b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
660b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
66172b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
66272b8c272SStefano Zampini       Mat change_sub;
66372b8c272SStefano Zampini 
664b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
66572b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
66672b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
66772b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
6687dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
66972b8c272SStefano Zampini       } else {
67072b8c272SStefano Zampini         Mat change_subt;
6717dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
67272b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
67372b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
67472b8c272SStefano Zampini       }
67572b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
67672b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
677e62b6521Sstefano_zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
678e62b6521Sstefano_zampini       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
67972b8c272SStefano Zampini     }
68072b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
68172b8c272SStefano Zampini   }
68272b8c272SStefano Zampini 
6835a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
6845a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
6855a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
6865a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
6875a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
6885a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
689aa83b6aeSStefano Zampini   }
690b1b3d7a2SStefano Zampini 
6915a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
692be83ff47SStefano Zampini   F = NULL;
693*d943a642SStefano Zampini   if (!sub_schurs->schur_explicit) {
694*d943a642SStefano Zampini     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
695*d943a642SStefano Zampini        it is not efficient, unless the economic version of the scaling is used */
6965a95e1ceSStefano Zampini     Mat         S_Ej_expl;
6975a95e1ceSStefano Zampini     PetscScalar *work;
6985a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
6995a95e1ceSStefano Zampini     PetscBool   Sdense;
7005a95e1ceSStefano Zampini 
7015a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
7025a95e1ceSStefano Zampini     local_size = 0;
703b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7045a95e1ceSStefano Zampini       IS  is_subset_B;
7055a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
7065a95e1ceSStefano Zampini 
7075a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7085a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7095a95e1ceSStefano Zampini       /* EE block */
7107dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7115a95e1ceSStefano Zampini       /* IE block */
7127dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7135a95e1ceSStefano Zampini       /* EI block */
714*d943a642SStefano Zampini       if (sub_schurs->is_symmetric) {
7155a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
716*d943a642SStefano Zampini       } else if (sub_schurs->is_hermitian) {
717*d943a642SStefano Zampini         ierr = MatCreateHermitianTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7185a95e1ceSStefano Zampini       } else {
7197dae84e0SHong Zhang         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7205a95e1ceSStefano Zampini       }
721a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7225a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7235a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7245a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7255a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
726b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
727b1b3d7a2SStefano Zampini         KSP ksp;
728b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7295a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
730b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
731b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
732b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
733b1b3d7a2SStefano Zampini         KSPType   ksp_type;
734b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7355a95e1ceSStefano Zampini         PetscBool ispcnone;
736b1b3d7a2SStefano Zampini 
737b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7385a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
739b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
740b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
741b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
742b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7435a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7445a95e1ceSStefano Zampini         if (!ispcnone) {
7455a95e1ceSStefano Zampini           PCType pc_type;
746b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
747b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7485a95e1ceSStefano Zampini         } else {
7495a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7505a95e1ceSStefano Zampini         }
751b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
752b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
7533ca39a21SBarry Smith           MatSolverType solver = NULL;
754ea799195SBarry Smith           ierr = PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);CHKERRQ(ierr);
755b1b3d7a2SStefano Zampini           if (solver) {
7563ca39a21SBarry Smith             ierr = PCFactorSetMatSolverType(schurpc,solver);CHKERRQ(ierr);
757b1b3d7a2SStefano Zampini           }
758b1b3d7a2SStefano Zampini         }
759b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
760b1b3d7a2SStefano Zampini       }
7615a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7625a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
763*d943a642SStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
7645a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
7655a95e1ceSStefano Zampini       if (Sdense) {
7665a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
7675a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
768b1b3d7a2SStefano Zampini         }
7695a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
7706c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
7715a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
772a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
7735a95e1ceSStefano Zampini       local_size += subset_size;
7745a95e1ceSStefano Zampini     }
7755a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
776b1b3d7a2SStefano Zampini     /* free */
777b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
778b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
7795a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
780883469d8SStefano Zampini   } else {
7815cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
7829d54b7f4SStefano Zampini     Vec         Dall = NULL;
783ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
7845cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
7856dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
7863fc34f97SStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE;
7873fc34f97SStefano Zampini     PetscBool   schur_has_vertices,factor_workaround;
78811955456SStefano Zampini     PetscBool   use_cholesky;
789883469d8SStefano Zampini 
790683d3df6SStefano Zampini     /* get sizes */
79181ea8064SStefano Zampini     n_I = 0;
79281ea8064SStefano Zampini     if (is_I_layer) {
793a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
79481ea8064SStefano Zampini     }
795683d3df6SStefano Zampini     economic = PETSC_FALSE;
796683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
797683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
798683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
7999d54b7f4SStefano Zampini     size_active_schur = local_size;
8009d54b7f4SStefano Zampini 
801f17d2ae1SStefano Zampini     /* import scaling vector (wrong formulation if we have 3D edges) */
8029d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
8039d54b7f4SStefano Zampini       const PetscScalar *array;
8049d54b7f4SStefano Zampini       PetscScalar       *array2;
8059d54b7f4SStefano Zampini       const PetscInt    *idxs;
8069d54b7f4SStefano Zampini       PetscInt          i;
8079d54b7f4SStefano Zampini 
8089d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8099d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8109d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8119d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8129d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8139d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8149d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8159d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8169d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8179d54b7f4SStefano Zampini     }
818d62866d3SStefano Zampini 
819683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
8203fc34f97SStefano Zampini     factor_workaround = PETSC_FALSE;
8213fc34f97SStefano Zampini     schur_has_vertices = PETSC_FALSE;
822683d3df6SStefano Zampini     cum = n_I+size_active_schur;
823683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
824683d3df6SStefano Zampini       const PetscInt* idxs;
825683d3df6SStefano Zampini       PetscInt        n_dir;
826683d3df6SStefano Zampini 
827683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
828683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
829683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
830683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
831683d3df6SStefano Zampini       cum += n_dir;
8323fc34f97SStefano Zampini       factor_workaround = PETSC_TRUE;
833d62866d3SStefano Zampini     }
834683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
835367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
836683d3df6SStefano Zampini       PetscInt n_v;
837683d3df6SStefano Zampini 
838683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
839683d3df6SStefano Zampini       if (n_v) {
840683d3df6SStefano Zampini         const PetscInt* idxs;
841683d3df6SStefano Zampini 
842683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
843683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
844683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
845683d3df6SStefano Zampini         cum += n_v;
846683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
8473fc34f97SStefano Zampini         schur_has_vertices = PETSC_TRUE;
848683d3df6SStefano Zampini       }
849683d3df6SStefano Zampini     }
850683d3df6SStefano Zampini     size_schur = cum - n_I;
851683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
852683d3df6SStefano Zampini     if (cum == n) {
853683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
854683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
855683d3df6SStefano Zampini     } else {
8567dae84e0SHong Zhang       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
857683d3df6SStefano Zampini     }
858e62b6521Sstefano_zampini     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
859e62b6521Sstefano_zampini     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
860ca92afb2SStefano Zampini 
861ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
862367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
863ca92afb2SStefano Zampini     if (benign_n) {
864ca92afb2SStefano Zampini       Vec                    v;
865ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
866ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
8675cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
8685cbda25cSStefano Zampini       PetscInt               sizeA;
869ca92afb2SStefano Zampini 
870ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
871ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
872ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
873ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
874ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
8755cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
8765cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
877282d6408SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
8785cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
8795cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
880ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
881ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
882ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
8835cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
884ca92afb2SStefano Zampini         PetscScalar    *array;
885ca92afb2SStefano Zampini         const PetscInt *idxs;
886ca92afb2SStefano Zampini         PetscInt       j,nz;
887ca92afb2SStefano Zampini 
888ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
889ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
890ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8915cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
892ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8935cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
8945cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
8955cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
896ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
89722db5ddcSStefano Zampini         for (j=0;j<size_schur;j++) {
89822db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
89922db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
90022db5ddcSStefano Zampini #else
90122db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
90222db5ddcSStefano Zampini #endif
90322db5ddcSStefano Zampini         }
904ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
905ca92afb2SStefano Zampini       }
9065cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9075cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
908ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
909ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
910ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
911ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
912ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
913ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
914ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
915ca92afb2SStefano Zampini     }
91611955456SStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);CHKERRQ(ierr);
9173301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9183301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
919883469d8SStefano Zampini 
92011955456SStefano Zampini     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
92111955456SStefano Zampini     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
92211955456SStefano Zampini 
923683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
924d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
925d47842beSStefano Zampini 
926f4f7d9d6SStefano Zampini     if (n_I) {
9270aa714b2SStefano Zampini       IS is_schur;
9285a05ddb0SStefano Zampini 
92911955456SStefano Zampini       if (use_cholesky) {
93088113c35SStefano Zampini         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
931883469d8SStefano Zampini       } else {
93288113c35SStefano Zampini         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
933883469d8SStefano Zampini       }
9348c09ecd8Sstefano_zampini       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
9358c09ecd8Sstefano_zampini 
936883469d8SStefano Zampini       /* subsets ordered last */
9375a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9385a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9395a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
940883469d8SStefano Zampini 
941883469d8SStefano Zampini       /* factorization step */
94211955456SStefano Zampini       if (use_cholesky) {
943883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
944be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
945be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
946be83ff47SStefano Zampini #endif
947883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
948a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
949883469d8SStefano Zampini       } else {
950883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
951be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
952be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
953be83ff47SStefano Zampini #endif
954883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
955883469d8SStefano Zampini       }
956ec968b6dSstefano_zampini       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
957883469d8SStefano Zampini 
95811955456SStefano Zampini #if defined(SUB_SCHURS_DEBUG)
95911955456SStefano Zampini       {
96011955456SStefano Zampini         PetscViewer viewer;
96111955456SStefano Zampini         char        filename[256];
96211955456SStefano Zampini         Mat         S;
96311955456SStefano Zampini         IS          is;
96411955456SStefano Zampini 
96511955456SStefano Zampini         sprintf(filename,"sub_schurs_Schur_r%d.m",PetscGlobalRank);
96611955456SStefano Zampini         ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
96711955456SStefano Zampini         ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
96811955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr);
96911955456SStefano Zampini         ierr = MatView(A,viewer);CHKERRQ(ierr);
97011955456SStefano Zampini         ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr);
97111955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)S,"S");CHKERRQ(ierr);
97211955456SStefano Zampini         ierr = MatView(S,viewer);CHKERRQ(ierr);
97311955456SStefano Zampini         ierr = MatDestroy(&S);CHKERRQ(ierr);
97411955456SStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);CHKERRQ(ierr);
97511955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)is,"I");CHKERRQ(ierr);
97611955456SStefano Zampini         ierr = ISView(is,viewer);CHKERRQ(ierr);
97711955456SStefano Zampini         ierr = ISDestroy(&is);CHKERRQ(ierr);
97811955456SStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);CHKERRQ(ierr);
97911955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)is,"B");CHKERRQ(ierr);
98011955456SStefano Zampini         ierr = ISView(is,viewer);CHKERRQ(ierr);
98111955456SStefano Zampini         ierr = ISDestroy(&is);CHKERRQ(ierr);
98211955456SStefano Zampini         ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
98311955456SStefano Zampini       }
98411955456SStefano Zampini #endif
98511955456SStefano Zampini 
986883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
9877c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
988b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
989b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
990b3cb21ddSStefano Zampini 
991d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
992683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
993683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
99403dfb2d7SStefano Zampini       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
99503dfb2d7SStefano Zampini         reuse_solvers = factor_workaround = PETSC_FALSE;
99603dfb2d7SStefano Zampini 
997df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
998ca92afb2SStefano Zampini 
99972b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
1000ca92afb2SStefano Zampini       if (benign_n) {
10015cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1002ca92afb2SStefano Zampini         Vec         v;
10035cbda25cSStefano Zampini         PetscInt    sizeA;
1004ca92afb2SStefano Zampini 
1005ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
10065cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
10075cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1008ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1009ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1010ca92afb2SStefano Zampini #endif
1011ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1012ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1013ca92afb2SStefano Zampini #endif
10145cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10155cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1016ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
10175cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
101847484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
1019ca92afb2SStefano Zampini           const PetscInt *idxs;
1020ca92afb2SStefano Zampini           PetscInt       j,nz;
102147484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
1022ca92afb2SStefano Zampini 
10235cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
10245cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
10255cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1026ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1027ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1028282d6408SStefano Zampini           /* p0 dof (eliminated) is excluded from the sum */
1029282d6408SStefano Zampini           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1030ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10315cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1032ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
103347484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
103447484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
103547484b83SStefano Zampini           B_k = 1;
103647484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
103747484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
10385cbda25cSStefano Zampini           sum  = 1.;
103947484b83SStefano 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));
1040ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10415cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10425cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1043282d6408SStefano Zampini           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
10445cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10455cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1046ca92afb2SStefano Zampini         }
1047a7414863SStefano Zampini         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1048a7414863SStefano Zampini           PetscInt k,j;
1049a7414863SStefano Zampini           for (k=0;k<size_schur;k++) {
1050a7414863SStefano Zampini             for (j=k;j<size_schur;j++) {
1051a7414863SStefano Zampini               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1052a7414863SStefano Zampini             }
1053a7414863SStefano Zampini           }
1054a7414863SStefano Zampini         }
1055a7414863SStefano Zampini 
10565cbda25cSStefano Zampini         /* restore defaults */
10575cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
10585cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
10595cbda25cSStefano Zampini #endif
10605cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
10615cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
10625cbda25cSStefano Zampini #endif
10635cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10645cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1065ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1066ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1067ca92afb2SStefano Zampini       }
1068a3df083aSStefano Zampini       if (!reuse_solvers) {
1069a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1070a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1071a3df083aSStefano Zampini         }
1072a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
10735cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
10745cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1075a3df083aSStefano Zampini       }
1076df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1077be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1078683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1079166598c1SStefano Zampini       factor_workaround = PETSC_FALSE;
1080df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1081be83ff47SStefano Zampini     }
1082be83ff47SStefano Zampini 
1083be83ff47SStefano Zampini     if (reuse_solvers) {
1084a00504b5SStefano Zampini       Mat                A_II,Afake;
108553892102SStefano Zampini       Vec                vec1_B;
1086df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
10873462e049SStefano Zampini       PetscInt           n_R;
1088d5574798SStefano Zampini 
1089df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1090df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1091e28d306cSStefano Zampini       } else {
1092df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1093d62866d3SStefano Zampini       }
1094df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1095be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1096d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1097d5574798SStefano Zampini       msolv_ctx->F = F;
1098683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1099683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1100683d3df6SStefano Zampini       {
1101683d3df6SStefano Zampini         PetscScalar *array;
1102683d3df6SStefano Zampini         PetscInt    n;
1103683d3df6SStefano Zampini 
1104683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1105683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1106683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1107683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1108683d3df6SStefano Zampini       }
11093fc34f97SStefano Zampini       msolv_ctx->has_vertices = schur_has_vertices;
1110d62866d3SStefano Zampini 
1111d62866d3SStefano Zampini       /* interior solver */
1112683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1113d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1114d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1115d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1116df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1117df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1118d62866d3SStefano Zampini 
1119d62866d3SStefano Zampini       /* correction solver */
11203462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1121d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1122d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1123df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1124df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
112553892102SStefano Zampini 
112653892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
112753892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
112853892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
11293fc34f97SStefano Zampini       if (!schur_has_vertices) {
113053892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
113153892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
113253892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
113353892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1134683d3df6SStefano Zampini       } else {
1135683d3df6SStefano Zampini         IS              is_B_all;
1136683d3df6SStefano Zampini         const PetscInt* idxs;
1137683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1138683d3df6SStefano Zampini 
1139683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1140683d3df6SStefano Zampini         dual = size_schur - n_v;
1141683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1142683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1143683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1144683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1145683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1146683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1147683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1148683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1149683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1150683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1151683d3df6SStefano Zampini       }
11523462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
11533462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
11548d692618SStefano Zampini       ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11558d692618SStefano Zampini       ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11563462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
11573462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1158683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1159ca92afb2SStefano Zampini 
1160ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1161ca92afb2SStefano Zampini       if (benign_n) {
11625cbda25cSStefano Zampini         PetscScalar *array;
11635cbda25cSStefano Zampini 
1164ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1165ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1166ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
11675cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
1168282d6408SStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
11695cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11705cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
11715cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11725cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1173ca92afb2SStefano Zampini       }
1174ada6e2d7SStefano Zampini     } else {
1175ada6e2d7SStefano Zampini       if (sub_schurs->reuse_solver) {
1176ada6e2d7SStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1177ada6e2d7SStefano Zampini       }
1178ada6e2d7SStefano Zampini       ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1179d5574798SStefano Zampini     }
118008122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
118153892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
11825db18549SStefano Zampini 
1183be83ff47SStefano Zampini     /* Work arrays */
1184be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
11855a94c880SStefano Zampini 
118604270d10SStefano Zampini     /* matrices for deluxe scaling and adaptive selection */
118712d906b1SStefano Zampini     if (compute_Stilda) {
11885a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
1189b449dd42SStefano Zampini         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11905a94c880SStefano Zampini       }
11919d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1192b449dd42SStefano Zampini         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
11935a94c880SStefano Zampini       }
119412d906b1SStefano Zampini     }
1195d2627357SStefano Zampini 
1196be83ff47SStefano Zampini     /* S_Ej_all */
1197be83ff47SStefano Zampini     cum = cum2 = 0;
1198be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
119965d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
120065d8bf0aSStefano Zampini       PetscInt j;
120165d8bf0aSStefano Zampini 
12025a95e1ceSStefano Zampini       /* get S_E */
1203b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1204683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1205be83ff47SStefano Zampini         PetscInt k;
1206be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1207be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1208be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1209683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1210be83ff47SStefano Zampini           }
1211be83ff47SStefano Zampini         }
121206a4e24aSStefano Zampini       } else { /* just copy to workspace */
1213be83ff47SStefano Zampini         PetscInt k;
1214be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1215be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1216be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1217be83ff47SStefano Zampini           }
1218be83ff47SStefano Zampini         }
12199087bf02SStefano Zampini       }
12205a95e1ceSStefano Zampini       /* insert S_E values */
1221be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1222b7ab4a40SStefano Zampini       if (sub_schurs->change) {
12238760537fSStefano Zampini         Mat change_sub,SEj,T;
122472b8c272SStefano Zampini 
122572b8c272SStefano Zampini         /* change basis */
122672b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
122772b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12288760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
12298760537fSStefano Zampini           Mat T2;
12308760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
12318760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
12328760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
123372b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
12348760537fSStefano Zampini         } else {
12358760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
123672b8c272SStefano Zampini         }
123772b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
123872b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
1239afa0f562SStefano Zampini         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
124072b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
124172b8c272SStefano Zampini       }
12429d54b7f4SStefano Zampini       if (deluxe) {
12435a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1244683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1245862806e4SStefano Zampini         if (compute_Stilda) {
1246f4f7d9d6SStefano Zampini           PetscInt k;
1247f4f7d9d6SStefano Zampini 
124808122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
124908122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1250f4f7d9d6SStefano Zampini           if (use_potr) {
12515a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12522972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
12535a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
12542972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
12556c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
12566c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
12576c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
12586c3e6151SStefano Zampini               }
12596c3e6151SStefano Zampini             }
1260f4f7d9d6SStefano Zampini           } else if (use_sytr) {
1261f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1262f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1263f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,work,&B_N,pivots,Bwork,&B_ierr));
1264f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1265f4f7d9d6SStefano Zampini             for (k=0;k<subset_size;k++) {
1266f4f7d9d6SStefano Zampini               for (j=k;j<subset_size;j++) {
1267f4f7d9d6SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
1268f4f7d9d6SStefano Zampini               }
1269f4f7d9d6SStefano Zampini             }
1270d6462365SStefano Zampini           } else {
1271d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1272d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1273d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1274d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
12752972d61bSStefano Zampini           }
127608122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
12775a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12789087bf02SStefano Zampini         }
12799d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
12809d54b7f4SStefano Zampini         Mat         SEj;
12819d54b7f4SStefano Zampini         Vec         D;
12829d54b7f4SStefano Zampini         PetscScalar *array;
12839d54b7f4SStefano Zampini 
12849d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12859d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
12869d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
12879d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1288f17d2ae1SStefano Zampini         ierr = VecShift(D,-1.);CHKERRQ(ierr);
12899d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
12909d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
12919d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
12929d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12939d54b7f4SStefano Zampini       }
1294be83ff47SStefano Zampini       cum += subset_size;
1295be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1296be83ff47SStefano Zampini     }
1297be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
12984a6c6b0dSStefano Zampini 
1299df4d28bfSStefano Zampini     if (solver_S) {
13007c2f51b8SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
13014a6c6b0dSStefano Zampini     }
1302683d3df6SStefano Zampini 
1303683d3df6SStefano Zampini     schur_factor = NULL;
130445951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1305683d3df6SStefano Zampini 
13069d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
13074a6c6b0dSStefano Zampini         PetscInt j;
13084a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
13095a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13104a6c6b0dSStefano Zampini       } else {
1311683d3df6SStefano Zampini         Mat S_all_inv=NULL;
13123fc34f97SStefano Zampini         if (solver_S) {
1313683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1314683d3df6SStefano 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 */
13153fc34f97SStefano Zampini           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1316683d3df6SStefano Zampini             PetscScalar *data;
1317683d3df6SStefano Zampini             PetscInt     nd = 0;
13186dba178dSStefano Zampini 
1319f4f7d9d6SStefano Zampini             if (!use_potr) {
1320683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1321683d3df6SStefano Zampini             }
13227c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1323683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1324683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1325683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1326683d3df6SStefano Zampini             }
13273fc34f97SStefano Zampini 
13283fc34f97SStefano Zampini             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
13293fc34f97SStefano Zampini             if (schur_has_vertices) {
13303fc34f97SStefano Zampini               Mat          M;
13313fc34f97SStefano Zampini               PetscScalar *tdata;
13323fc34f97SStefano Zampini               PetscInt     nv = 0, news;
13333fc34f97SStefano Zampini 
13343fc34f97SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
13353fc34f97SStefano Zampini               news = size_active_schur + nv;
13363fc34f97SStefano Zampini               ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr);
1337683d3df6SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13383fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
13393fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));CHKERRQ(ierr);
13403fc34f97SStefano Zampini               }
13413fc34f97SStefano Zampini               for (i=0;i<nv;i++) {
13423fc34f97SStefano Zampini                 PetscInt k = i+size_active_schur;
13433fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));CHKERRQ(ierr);
13443fc34f97SStefano Zampini               }
13453fc34f97SStefano Zampini 
13463fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr);
13473fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
13483fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
13493fc34f97SStefano Zampini               /* save the factors */
13503fc34f97SStefano Zampini               cum = 0;
13513fc34f97SStefano Zampini               ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr);
13523fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13533fc34f97SStefano Zampini                 ierr = PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1354683d3df6SStefano Zampini                 cum += size_active_schur - i;
1355683d3df6SStefano Zampini               }
13563fc34f97SStefano Zampini               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
13573fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
13583fc34f97SStefano Zampini               /* move back just the active dofs to the Schur complement */
13593fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13603fc34f97SStefano Zampini                 ierr = PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));CHKERRQ(ierr);
13613fc34f97SStefano Zampini               }
13623fc34f97SStefano Zampini               ierr = PetscFree(tdata);CHKERRQ(ierr);
13633fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13643fc34f97SStefano Zampini             } else { /* we can factorize and invert just the activedofs */
13653fc34f97SStefano Zampini               Mat         M;
13665002105bSStefano Zampini               PetscScalar *aux;
13673fc34f97SStefano Zampini 
13685002105bSStefano Zampini               ierr = PetscMalloc1(nd,&aux);CHKERRQ(ierr);
13695002105bSStefano Zampini               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
13703fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr);
13713fc34f97SStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
13723fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
13733fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
13743fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
13753fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13765002105bSStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);CHKERRQ(ierr);
13775002105bSStefano Zampini               ierr = MatZeroEntries(M);CHKERRQ(ierr);
13785002105bSStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13795002105bSStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);CHKERRQ(ierr);
13805002105bSStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
13815002105bSStefano Zampini               ierr = MatZeroEntries(M);CHKERRQ(ierr);
13825002105bSStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13835002105bSStefano Zampini               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
13845002105bSStefano Zampini               ierr = PetscFree(aux);CHKERRQ(ierr);
13853fc34f97SStefano Zampini             }
1386683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
13873fc34f97SStefano Zampini           } else { /* use MatFactor calls to invert S */
13886dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
13897c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1390683d3df6SStefano Zampini           }
1391df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1392683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1393683d3df6SStefano Zampini           S_all_inv = S_all;
1394683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1395be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1396be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1397f4f7d9d6SStefano Zampini           if (use_potr) {
1398be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1399be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1400be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1401be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1402f4f7d9d6SStefano Zampini           } else if (use_sytr) {
1403f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1404f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1405f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1406f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1407d6462365SStefano Zampini           } else {
1408d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1409d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1410d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1411d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1412be83ff47SStefano Zampini           }
1413be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1414683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1415be83ff47SStefano Zampini         }
1416be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1417be83ff47SStefano Zampini         cum = cum2 = 0;
1418683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1419be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1420be83ff47SStefano Zampini           PetscInt j;
1421862806e4SStefano Zampini 
1422862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1423be83ff47SStefano Zampini           /* get (St^-1)_E */
142472b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
142506a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1426a0b0af32SStefano Zampini           if (S_lower_triangular) {
1427be83ff47SStefano Zampini             PetscInt k;
1428b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1429be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1430be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1431be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
14326c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1433be83ff47SStefano Zampini                 }
1434be83ff47SStefano Zampini               }
143572b8c272SStefano Zampini             } else {
143672b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
143772b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
143872b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
143972b8c272SStefano Zampini                 }
144072b8c272SStefano Zampini               }
144172b8c272SStefano Zampini             }
144272b8c272SStefano Zampini           } else {
1443be83ff47SStefano Zampini             PetscInt k;
1444be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1445be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1446be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1447be83ff47SStefano Zampini               }
1448be83ff47SStefano Zampini             }
1449be83ff47SStefano Zampini           }
1450b7ab4a40SStefano Zampini           if (sub_schurs->change) {
14518760537fSStefano Zampini             Mat change_sub,SEj,T;
145272b8c272SStefano Zampini 
145372b8c272SStefano Zampini             /* change basis */
145472b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
145572b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
14568760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
14578760537fSStefano Zampini               Mat T2;
14588760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
14598760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
146072b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
14618760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
14628760537fSStefano Zampini             } else {
14638760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
146472b8c272SStefano Zampini             }
14658760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
14668760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
146772b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1468afa0f562SStefano Zampini             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
146972b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
147072b8c272SStefano Zampini           }
1471be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
14725a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1473be83ff47SStefano Zampini           cum += subset_size;
1474be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1475883469d8SStefano Zampini         }
1476683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1477df4d28bfSStefano Zampini         if (solver_S) {
14783fc34f97SStefano Zampini           if (schur_has_vertices) {
14793fc34f97SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
14803fc34f97SStefano Zampini           } else {
14817c2f51b8SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
14825db18549SStefano Zampini           }
14833fc34f97SStefano Zampini         }
1484683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1485683d3df6SStefano Zampini       }
1486683d3df6SStefano Zampini 
14873fc34f97SStefano Zampini       /* move back factors if needed */
14883fc34f97SStefano Zampini       if (schur_has_vertices) {
1489683d3df6SStefano Zampini         Mat      S_tmp;
14903fc34f97SStefano Zampini         PetscInt nd = 0;
1491683d3df6SStefano Zampini 
14926542c05fSStefano Zampini         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
14937c2f51b8SStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1494f4f7d9d6SStefano Zampini         if (use_potr) {
1495683d3df6SStefano Zampini           PetscScalar *data;
1496683d3df6SStefano Zampini 
1497683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
14983fc34f97SStefano Zampini           ierr = PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1499683d3df6SStefano Zampini 
1500683d3df6SStefano Zampini           if (S_lower_triangular) {
1501683d3df6SStefano Zampini             cum = 0;
1502683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1503683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1504683d3df6SStefano Zampini               cum += size_active_schur-i;
1505683d3df6SStefano Zampini 	    }
1506683d3df6SStefano Zampini           } else {
1507683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1508683d3df6SStefano Zampini           }
1509683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1510683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1511683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1512683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1513683d3df6SStefano Zampini 	    }
1514683d3df6SStefano Zampini           }
15156dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1516683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1517683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
15186c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1519683d3df6SStefano Zampini           }
1520683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
15216542c05fSStefano Zampini         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
15223fc34f97SStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
15239087bf02SStefano Zampini       }
1524367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1525367aa537SStefano Zampini       PetscScalar *data;
1526367aa537SStefano Zampini       PetscInt    nd = 0;
1527367aa537SStefano Zampini 
1528367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1529367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1530367aa537SStefano Zampini       }
15317c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1532367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1533367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1534367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1535367aa537SStefano Zampini       }
1536367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1537367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
15386c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1539367aa537SStefano Zampini       }
1540367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
15413fc34f97SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
15424a6c6b0dSStefano Zampini     }
15434a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1544683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
15459d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
15464a6c6b0dSStefano Zampini   }
15475a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1548be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1549a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1550a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1551a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1552a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1553a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
15545db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15555db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15565a95e1ceSStefano Zampini   if (compute_Stilda) {
15575a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15585a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15599d54b7f4SStefano Zampini     if (deluxe) {
15605a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15615a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156208122e43SStefano Zampini     }
15639d54b7f4SStefano Zampini   }
15645db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
15655db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
15663927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
15673927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15685a95e1ceSStefano Zampini 
15695db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
15707dae84e0SHong Zhang   ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
157141fd5443SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
157241fd5443SStefano Zampini     ierr = MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
157341fd5443SStefano Zampini   } else {
157404270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
157541fd5443SStefano Zampini   }
157608122e43SStefano Zampini 
1577f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
15785a95e1ceSStefano Zampini   if (compute_Stilda) {
15795a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1580a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
158104270d10SStefano Zampini     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
158204270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
15839d54b7f4SStefano Zampini     if (deluxe) {
15845a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
158508122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
158604270d10SStefano Zampini       ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
158704270d10SStefano Zampini       ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
15889d54b7f4SStefano Zampini     } else {
15899d54b7f4SStefano Zampini       PetscScalar *array;
15909d54b7f4SStefano Zampini       PetscInt    cum;
15919d54b7f4SStefano Zampini 
15929d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15939d54b7f4SStefano Zampini       cum = 0;
15949d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
15959d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
15969d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
15979d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1598f4f7d9d6SStefano Zampini         if (use_potr) {
15999d54b7f4SStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
16009d54b7f4SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
16019d54b7f4SStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
16029d54b7f4SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1603f4f7d9d6SStefano Zampini         } else if (use_sytr) {
1604f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1605f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1606f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1607f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1608f4f7d9d6SStefano Zampini         } else {
1609f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1610f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1611f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1612f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1613f4f7d9d6SStefano Zampini         }
16149d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
16159d54b7f4SStefano Zampini         cum += subset_size*subset_size;
16169d54b7f4SStefano Zampini       }
16179d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
16189d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
16199d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
16209d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
16219d54b7f4SStefano Zampini     }
162208122e43SStefano Zampini   }
162304270d10SStefano Zampini   ierr = MatDestroySubMatrices(1,&submats);CHKERRQ(ierr);
162411955456SStefano Zampini #if defined(SUB_SCHURS_DEBUG)
162511955456SStefano Zampini   {
162611955456SStefano Zampini     PetscViewer viewer;
162711955456SStefano Zampini     char        filename[256];
162811955456SStefano Zampini     PetscInt    cum;
162911955456SStefano Zampini 
163011955456SStefano Zampini     sprintf(filename,"sub_schurs_mats_r%d.m",PetscGlobalRank);
163111955456SStefano Zampini     ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
163211955456SStefano Zampini     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
163311955456SStefano Zampini     if (sub_schurs->S_Ej_all) {
163411955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");CHKERRQ(ierr);
163511955456SStefano Zampini       ierr = MatView(sub_schurs->S_Ej_all,viewer);CHKERRQ(ierr);
163611955456SStefano Zampini     }
163711955456SStefano Zampini     if (sub_schurs->sum_S_Ej_all) {
163811955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");CHKERRQ(ierr);
163911955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_all,viewer);CHKERRQ(ierr);
164011955456SStefano Zampini     }
164111955456SStefano Zampini     if (sub_schurs->sum_S_Ej_inv_all) {
164211955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");CHKERRQ(ierr);
164311955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_inv_all,viewer);CHKERRQ(ierr);
164411955456SStefano Zampini     }
164511955456SStefano Zampini     if (sub_schurs->sum_S_Ej_tilda_all) {
164611955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");CHKERRQ(ierr);
164711955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_tilda_all,viewer);CHKERRQ(ierr);
164811955456SStefano Zampini     }
164911955456SStefano Zampini     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
165011955456SStefano Zampini       IS   is;
165111955456SStefano Zampini       char name[16];
165211955456SStefano Zampini 
165311955456SStefano Zampini       sprintf(name,"IE%d",i);
165411955456SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
165511955456SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);CHKERRQ(ierr);
165611955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)is,name);CHKERRQ(ierr);
165711955456SStefano Zampini       ierr = ISView(is,viewer);CHKERRQ(ierr);
165811955456SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
165911955456SStefano Zampini       cum += subset_size;
166011955456SStefano Zampini     }
166111955456SStefano Zampini     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
166211955456SStefano Zampini   }
166311955456SStefano Zampini #endif
16643202ece2SStefano Zampini 
16655a95e1ceSStefano Zampini   /* free workspace */
16665a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
166706a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1668a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
16693202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1670dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
16715a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1672b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1673b1b3d7a2SStefano Zampini }
1674b1b3d7a2SStefano Zampini 
167588113c35SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1676b1b3d7a2SStefano Zampini {
16779bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
16785a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
167988113c35SStefano Zampini   PetscBool       is_sorted,ispetsc;
1680b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1681b1b3d7a2SStefano Zampini 
1682b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1683b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
16846c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1685b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
16866c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1687b1b3d7a2SStefano Zampini 
1688b1b3d7a2SStefano Zampini   /* reset any previous data */
1689b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1690b1b3d7a2SStefano Zampini 
16915a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
16929bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1693b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
169408122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1695b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1696b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
16978b6046baSStefano Zampini     if (copycc) {
16988b6046baSStefano Zampini       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
16998b6046baSStefano Zampini     } else {
1700c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1701b1b3d7a2SStefano Zampini       all_cc[i] = faces[i];
1702b1b3d7a2SStefano Zampini     }
17038b6046baSStefano Zampini   }
1704b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
17058b6046baSStefano Zampini     if (copycc) {
17068b6046baSStefano Zampini       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
17078b6046baSStefano Zampini     } else {
1708c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1709b1b3d7a2SStefano Zampini       all_cc[n_faces+i] = edges[i];
17108b6046baSStefano Zampini     }
171108122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1712b1b3d7a2SStefano Zampini   }
1713c8272957SStefano Zampini   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1714c8272957SStefano Zampini   sub_schurs->is_vertices = vertices;
1715c8272957SStefano Zampini   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1716d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1717d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1718b1b3d7a2SStefano Zampini 
1719df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
172088113c35SStefano Zampini   ierr = PetscStrallocpy(prefix,&sub_schurs->prefix);CHKERRQ(ierr);
1721883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
172288113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);CHKERRQ(ierr);
172388113c35SStefano Zampini #elif defined(PETSC_HAVE_MKL_PARDISO)
172488113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);CHKERRQ(ierr);
172588113c35SStefano Zampini #else
172688113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);CHKERRQ(ierr);
1727df4d28bfSStefano Zampini #endif
172888113c35SStefano Zampini #if defined(PETSC_USE_COMPLEX)
172988113c35SStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
173088113c35SStefano Zampini #else
173188113c35SStefano Zampini   sub_schurs->is_hermitian = PETSC_TRUE;
1732883469d8SStefano Zampini #endif
173388113c35SStefano Zampini   sub_schurs->is_posdef    = PETSC_TRUE;
173411955456SStefano Zampini   sub_schurs->is_symmetric = PETSC_TRUE;
173588113c35SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
173688113c35SStefano Zampini   ierr = PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,64,NULL);CHKERRQ(ierr);
173711955456SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);CHKERRQ(ierr);
173888113c35SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
173988113c35SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
174088113c35SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
174188113c35SStefano Zampini   ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);CHKERRQ(ierr);
1742168c3995SStefano Zampini   sub_schurs->schur_explicit = (PetscBool)!ispetsc;
1743b1b3d7a2SStefano Zampini 
174411955456SStefano Zampini   /* for reals, symmetric and hermitian are synonims */
174511955456SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
174611955456SStefano Zampini   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
174711955456SStefano Zampini   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
174811955456SStefano Zampini #endif
174911955456SStefano Zampini 
1750b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1751b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1752b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1753b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
17545db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
17555db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
17565db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
17575db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
17585a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1759b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1760b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1761b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
176208122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1763b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1764b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1765b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1766b1b3d7a2SStefano Zampini }
1767b1b3d7a2SStefano Zampini 
176834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
176934a97f8cSStefano Zampini {
177034a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
177134a97f8cSStefano Zampini   PetscErrorCode  ierr;
177234a97f8cSStefano Zampini 
177334a97f8cSStefano Zampini   PetscFunctionBegin;
177434a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
17755ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
177634a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
177734a97f8cSStefano Zampini   PetscFunctionReturn(0);
177834a97f8cSStefano Zampini }
177934a97f8cSStefano Zampini 
178034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
178134a97f8cSStefano Zampini {
178234a97f8cSStefano Zampini   PetscInt       i;
178334a97f8cSStefano Zampini   PetscErrorCode ierr;
178434a97f8cSStefano Zampini 
178534a97f8cSStefano Zampini   PetscFunctionBegin;
1786aea80f77Sstefano_zampini   if (!sub_schurs) PetscFunctionReturn(0);
178788113c35SStefano Zampini   ierr = PetscFree(sub_schurs->prefix);CHKERRQ(ierr);
17881e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1789b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1790b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1791b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
17925db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
17935db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
179441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
179541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
179608122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1797a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
17985db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1799d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1800d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
180108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
180234a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1803b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
180434a97f8cSStefano Zampini   }
18055ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1806b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
18073dc780c3SStefano Zampini   }
1808df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1809df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1810d62866d3SStefano Zampini   }
1811df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
181272b8c272SStefano Zampini   if (sub_schurs->change) {
181372b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
181472b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1815b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
181672b8c272SStefano Zampini     }
181772b8c272SStefano Zampini   }
181872b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1819b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
182034a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
182134a97f8cSStefano Zampini   PetscFunctionReturn(0);
182234a97f8cSStefano Zampini }
182334a97f8cSStefano Zampini 
1824aea80f77Sstefano_zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1825aea80f77Sstefano_zampini {
1826aea80f77Sstefano_zampini   PetscErrorCode ierr;
1827aea80f77Sstefano_zampini 
1828aea80f77Sstefano_zampini   PetscFunctionBegin;
1829aea80f77Sstefano_zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1830aea80f77Sstefano_zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1831aea80f77Sstefano_zampini   PetscFunctionReturn(0);
1832aea80f77Sstefano_zampini }
1833aea80f77Sstefano_zampini 
18342a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
183534a97f8cSStefano Zampini {
183634a97f8cSStefano Zampini   PetscInt       i,j,n;
183734a97f8cSStefano Zampini   PetscErrorCode ierr;
183834a97f8cSStefano Zampini 
183934a97f8cSStefano Zampini   PetscFunctionBegin;
184034a97f8cSStefano Zampini   n = 0;
184134a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
184234a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
184334a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
184434a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
184534a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
184634a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
184734a97f8cSStefano Zampini         queue_tip[n] = dof;
184834a97f8cSStefano Zampini         n++;
184934a97f8cSStefano Zampini       }
185034a97f8cSStefano Zampini     }
185134a97f8cSStefano Zampini   }
185234a97f8cSStefano Zampini   *n_added = n;
185334a97f8cSStefano Zampini   PetscFunctionReturn(0);
185434a97f8cSStefano Zampini }
1855