xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 7f9db97b55e3f4b4485c94e04a216bedeaa8ef29)
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;
402*7f9db97bSStefano Zampini   PetscBool              print_schurs = PETSC_FALSE;
403b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
404b1b3d7a2SStefano Zampini 
405b1b3d7a2SStefano Zampini   PetscFunctionBegin;
406a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
407a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
408e62b6521Sstefano_zampini   /* convert matrix if needed */
409e62b6521Sstefano_zampini   if (Ain) {
410e62b6521Sstefano_zampini     PetscBool isseqaij;
411e62b6521Sstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
412a64f4aa4SStefano Zampini     if (isseqaij) {
413a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
414a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4153301b35fSStefano Zampini     } else {
416a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
417a64f4aa4SStefano Zampini     }
418a64f4aa4SStefano Zampini   }
4193301b35fSStefano Zampini 
420a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
421a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
422df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
423df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
424a64f4aa4SStefano Zampini   }
425a64f4aa4SStefano Zampini 
4265a95e1ceSStefano Zampini   /* preliminary checks */
427af25d912SStefano 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");
4285a95e1ceSStefano Zampini 
42988113c35SStefano Zampini   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
43088113c35SStefano Zampini 
431*7f9db97bSStefano Zampini   /* debug prints */
432*7f9db97bSStefano Zampini   if (sub_schurs->debug) {
433*7f9db97bSStefano Zampini     PetscMPIInt size,rank;
434*7f9db97bSStefano Zampini     PetscInt    nr,*print_schurs_ranks;
435*7f9db97bSStefano Zampini     PetscBool   flg;
436*7f9db97bSStefano Zampini 
437*7f9db97bSStefano Zampini     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);CHKERRQ(ierr);
438*7f9db97bSStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
439*7f9db97bSStefano Zampini     nr   = size;
440*7f9db97bSStefano Zampini     ierr = PetscMalloc1(nr,&print_schurs_ranks);CHKERRQ(ierr);
441*7f9db97bSStefano Zampini     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
442*7f9db97bSStefano Zampini     ierr = PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);CHKERRQ(ierr);
443*7f9db97bSStefano Zampini     if (!flg) print_schurs = PETSC_TRUE;
444*7f9db97bSStefano Zampini     else {
445*7f9db97bSStefano Zampini       for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; }
446*7f9db97bSStefano Zampini     }
447*7f9db97bSStefano Zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
448*7f9db97bSStefano Zampini     ierr = PetscFree(print_schurs_ranks);CHKERRQ(ierr);
449*7f9db97bSStefano Zampini   }
450*7f9db97bSStefano Zampini 
451*7f9db97bSStefano Zampini 
4525a95e1ceSStefano Zampini   /* restrict work on active processes */
4535a95e1ceSStefano Zampini   color = 0;
4545a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4555a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4565a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4575a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4585a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4595a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
4605a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
4615a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4625a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
4635a95e1ceSStefano Zampini     PetscFunctionReturn(0);
4645a95e1ceSStefano Zampini   }
46581ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
4665a95e1ceSStefano Zampini 
467b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
468df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
469a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
4703301b35fSStefano Zampini     PetscBool isseqsbaij;
471a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
4723301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
4733301b35fSStefano Zampini     if (isseqsbaij) {
474a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
475a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
476a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
477a64f4aa4SStefano Zampini     } else {
478a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
479a64f4aa4SStefano Zampini       A_BB = tA_BB;
480a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
481a64f4aa4SStefano Zampini       A_IB = tA_IB;
482a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
483a64f4aa4SStefano Zampini       A_BI = tA_BI;
484f11841e3SStefano Zampini     }
485a58a30b4SStefano Zampini   } else {
4865a95e1ceSStefano Zampini     A_II = NULL;
4875a95e1ceSStefano Zampini     A_IB = NULL;
4885a95e1ceSStefano Zampini     A_BI = NULL;
4895a95e1ceSStefano Zampini     A_BB = NULL;
490b1b3d7a2SStefano Zampini   }
4915a95e1ceSStefano Zampini   S_all = NULL;
492b1b3d7a2SStefano Zampini 
493b1b3d7a2SStefano Zampini   /* determine interior problems */
4943dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4953dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
496b1b3d7a2SStefano Zampini     PetscBT                touched;
497b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
498b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
499b1b3d7a2SStefano Zampini 
5001633d1f0SStefano Zampini     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
501b1b3d7a2SStefano Zampini     /* get sizes */
502b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
503b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
504b1b3d7a2SStefano Zampini 
505b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
506b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
507b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
508b1b3d7a2SStefano Zampini 
509b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
510b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
511b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
512b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
513b1b3d7a2SStefano Zampini     }
514b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
515b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
516b1b3d7a2SStefano Zampini 
517b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
518b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
519b1b3d7a2SStefano Zampini     n_prev_added = n_B;
520b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
521b1b3d7a2SStefano Zampini       PetscInt n_added;
522b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5236c4ed002SBarry 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);
524b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
525b1b3d7a2SStefano Zampini       n_prev_added = n_added;
526b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
527b1b3d7a2SStefano Zampini       if (!n_added) break;
528b1b3d7a2SStefano Zampini     }
529b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
530b1b3d7a2SStefano Zampini 
531883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
532a9b99552SStefano 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);
533b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
534a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
535883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
536df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
537b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
538b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
539a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
540b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
541b1b3d7a2SStefano Zampini 
542b1b3d7a2SStefano Zampini       /* II block */
5437dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
544b1b3d7a2SStefano Zampini     }
545b1b3d7a2SStefano Zampini   } else {
546b1b3d7a2SStefano Zampini     PetscInt n_I;
547b1b3d7a2SStefano Zampini 
548b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
549b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
550a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
551b1b3d7a2SStefano Zampini 
552b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
553df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
554b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
555b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
556b1b3d7a2SStefano Zampini 
557b1b3d7a2SStefano Zampini       /* II block is the same */
558b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
559b1b3d7a2SStefano Zampini       AE_II = A_II;
560b1b3d7a2SStefano Zampini     }
561b1b3d7a2SStefano Zampini   }
5625a95e1ceSStefano Zampini 
563883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
5645a95e1ceSStefano Zampini   max_subset_size = 0;
565883469d8SStefano Zampini   local_size = 0;
5665a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
5675a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5685a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
569883469d8SStefano Zampini     local_size += subset_size;
570883469d8SStefano Zampini   }
571883469d8SStefano Zampini 
572883469d8SStefano Zampini   /* Work arrays for local indices */
573883469d8SStefano Zampini   extra = 0;
574683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
575df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
576a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
577883469d8SStefano Zampini   }
578683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
579883469d8SStefano Zampini   if (extra) {
580883469d8SStefano Zampini     const PetscInt *idxs;
581a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
582883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
583a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
584883469d8SStefano Zampini   }
585883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
586dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
587dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
588883469d8SStefano Zampini 
589883469d8SStefano Zampini   /* Get local indices in local numbering */
590883469d8SStefano Zampini   local_size = 0;
5915a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
592883469d8SStefano Zampini     PetscInt j;
593883469d8SStefano Zampini     const    PetscInt *idxs;
594883469d8SStefano Zampini 
5955a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5965a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
597eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
598eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
599eb595f79SStefano Zampini     auxnum2[i] = subset_size;
600883469d8SStefano Zampini     /* subset indices in local numbering */
601883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
6025a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
603883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
604883469d8SStefano Zampini     local_size += subset_size;
605883469d8SStefano Zampini   }
606883469d8SStefano Zampini 
607f4f7d9d6SStefano Zampini   /* allocate extra workspace needed only for GETRI or SYTRF */
60811955456SStefano Zampini   use_potr = use_sytr = PETSC_FALSE;
60911955456SStefano Zampini   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
610f4f7d9d6SStefano Zampini     use_potr = PETSC_TRUE;
61111955456SStefano Zampini   } else if (sub_schurs->is_symmetric) {
61211955456SStefano Zampini     use_sytr = PETSC_TRUE;
61311955456SStefano Zampini   }
61411955456SStefano Zampini   if (local_size && !use_potr) {
61559ac4de7SStefano Zampini     PetscScalar  lwork,dummyscalar = 0.;
61659ac4de7SStefano Zampini     PetscBLASInt dummyint = 0;
617d2627357SStefano Zampini 
618d2627357SStefano Zampini     B_lwork = -1;
619d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
620d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
621f4f7d9d6SStefano Zampini     if (use_sytr) {
622f4f7d9d6SStefano Zampini       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
623f4f7d9d6SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
624f4f7d9d6SStefano Zampini     } else {
62559ac4de7SStefano Zampini       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
626d2627357SStefano Zampini       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
627f4f7d9d6SStefano Zampini     }
628f4f7d9d6SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
629d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
630d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
631056290a2SStefano Zampini   } else {
632056290a2SStefano Zampini     Bwork = NULL;
633056290a2SStefano Zampini     pivots = NULL;
634d2627357SStefano Zampini   }
635d2627357SStefano Zampini 
636d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
637dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
638dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
639dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
640dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
6416583bcc1SStefano Zampini   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
642dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
643dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
644dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6456c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
646dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
647e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
648d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
649d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
650d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
651d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6522972d61bSStefano Zampini 
6535a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6545a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6555a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6565a95e1ceSStefano Zampini 
6575a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6585a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
6596c4ed002SBarry 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);
6605a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
661b1b3d7a2SStefano Zampini   }
662b1b3d7a2SStefano Zampini 
66372b8c272SStefano Zampini   if (change) {
66472b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
66572b8c272SStefano Zampini     IS                     change_primal_B;
66672b8c272SStefano Zampini     IS                     change_primal_all;
66772b8c272SStefano Zampini 
668b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
669b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
670b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
67172b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
67272b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
67372b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
674b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
67572b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
67672b8c272SStefano Zampini     }
67772b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
67872b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
67972b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
68072b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
68172b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
682b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
68372b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
68472b8c272SStefano Zampini       Mat change_sub;
68572b8c272SStefano Zampini 
686b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
68772b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
68872b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
68972b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
6907dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
69172b8c272SStefano Zampini       } else {
69272b8c272SStefano Zampini         Mat change_subt;
6937dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
69472b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
69572b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
69672b8c272SStefano Zampini       }
69772b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
69872b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
699e62b6521Sstefano_zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
700e62b6521Sstefano_zampini       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
70172b8c272SStefano Zampini     }
70272b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
70372b8c272SStefano Zampini   }
70472b8c272SStefano Zampini 
7055a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
7065a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
7075a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
7085a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
7095a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
7105a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
711aa83b6aeSStefano Zampini   }
712b1b3d7a2SStefano Zampini 
7135a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
714be83ff47SStefano Zampini   F = NULL;
715d943a642SStefano Zampini   if (!sub_schurs->schur_explicit) {
716d943a642SStefano Zampini     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
717d943a642SStefano Zampini        it is not efficient, unless the economic version of the scaling is used */
7185a95e1ceSStefano Zampini     Mat         S_Ej_expl;
7195a95e1ceSStefano Zampini     PetscScalar *work;
7205a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
7215a95e1ceSStefano Zampini     PetscBool   Sdense;
7225a95e1ceSStefano Zampini 
7235a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
7245a95e1ceSStefano Zampini     local_size = 0;
725b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7265a95e1ceSStefano Zampini       IS  is_subset_B;
7275a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
7285a95e1ceSStefano Zampini 
7295a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7305a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7315a95e1ceSStefano Zampini       /* EE block */
7327dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7335a95e1ceSStefano Zampini       /* IE block */
7347dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7355a95e1ceSStefano Zampini       /* EI block */
736d943a642SStefano Zampini       if (sub_schurs->is_symmetric) {
7375a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
738d943a642SStefano Zampini       } else if (sub_schurs->is_hermitian) {
739d943a642SStefano Zampini         ierr = MatCreateHermitianTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7405a95e1ceSStefano Zampini       } else {
7417dae84e0SHong Zhang         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7425a95e1ceSStefano Zampini       }
743a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7445a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7455a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7465a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7475a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
748b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
749b1b3d7a2SStefano Zampini         KSP ksp;
750b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7515a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
752b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
753b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
754b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
755b1b3d7a2SStefano Zampini         KSPType   ksp_type;
756b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7575a95e1ceSStefano Zampini         PetscBool ispcnone;
758b1b3d7a2SStefano Zampini 
759b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7605a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
761b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
762b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
763b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
764b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7655a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7665a95e1ceSStefano Zampini         if (!ispcnone) {
7675a95e1ceSStefano Zampini           PCType pc_type;
768b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
769b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7705a95e1ceSStefano Zampini         } else {
7715a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7725a95e1ceSStefano Zampini         }
773b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
774b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
7753ca39a21SBarry Smith           MatSolverType solver = NULL;
776ea799195SBarry Smith           ierr = PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);CHKERRQ(ierr);
777b1b3d7a2SStefano Zampini           if (solver) {
7783ca39a21SBarry Smith             ierr = PCFactorSetMatSolverType(schurpc,solver);CHKERRQ(ierr);
779b1b3d7a2SStefano Zampini           }
780b1b3d7a2SStefano Zampini         }
781b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
782b1b3d7a2SStefano Zampini       }
7835a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7845a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
785d943a642SStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
7865a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
7875a95e1ceSStefano Zampini       if (Sdense) {
7885a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
7895a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
790b1b3d7a2SStefano Zampini         }
7915a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
7926c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
7935a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
794a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
7955a95e1ceSStefano Zampini       local_size += subset_size;
7965a95e1ceSStefano Zampini     }
7975a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
798b1b3d7a2SStefano Zampini     /* free */
799b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
800b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
8015a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
802883469d8SStefano Zampini   } else {
8035cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
8049d54b7f4SStefano Zampini     Vec         Dall = NULL;
805ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
8065cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
8076dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
8083fc34f97SStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE;
8093fc34f97SStefano Zampini     PetscBool   schur_has_vertices,factor_workaround;
81011955456SStefano Zampini     PetscBool   use_cholesky;
811883469d8SStefano Zampini 
812683d3df6SStefano Zampini     /* get sizes */
81381ea8064SStefano Zampini     n_I = 0;
81481ea8064SStefano Zampini     if (is_I_layer) {
815a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
81681ea8064SStefano Zampini     }
817683d3df6SStefano Zampini     economic = PETSC_FALSE;
818683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
819683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
820683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
8219d54b7f4SStefano Zampini     size_active_schur = local_size;
8229d54b7f4SStefano Zampini 
823f17d2ae1SStefano Zampini     /* import scaling vector (wrong formulation if we have 3D edges) */
8249d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
8259d54b7f4SStefano Zampini       const PetscScalar *array;
8269d54b7f4SStefano Zampini       PetscScalar       *array2;
8279d54b7f4SStefano Zampini       const PetscInt    *idxs;
8289d54b7f4SStefano Zampini       PetscInt          i;
8299d54b7f4SStefano Zampini 
8309d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8319d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8329d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8339d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8349d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8359d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8369d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8379d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8389d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8399d54b7f4SStefano Zampini     }
840d62866d3SStefano Zampini 
841683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
8423fc34f97SStefano Zampini     factor_workaround = PETSC_FALSE;
8433fc34f97SStefano Zampini     schur_has_vertices = PETSC_FALSE;
844683d3df6SStefano Zampini     cum = n_I+size_active_schur;
845683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
846683d3df6SStefano Zampini       const PetscInt* idxs;
847683d3df6SStefano Zampini       PetscInt        n_dir;
848683d3df6SStefano Zampini 
849683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
850683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
851683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
852683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
853683d3df6SStefano Zampini       cum += n_dir;
8543fc34f97SStefano Zampini       factor_workaround = PETSC_TRUE;
855d62866d3SStefano Zampini     }
856683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
857367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
858683d3df6SStefano Zampini       PetscInt n_v;
859683d3df6SStefano Zampini 
860683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
861683d3df6SStefano Zampini       if (n_v) {
862683d3df6SStefano Zampini         const PetscInt* idxs;
863683d3df6SStefano Zampini 
864683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
865683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
866683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
867683d3df6SStefano Zampini         cum += n_v;
868683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
8693fc34f97SStefano Zampini         schur_has_vertices = PETSC_TRUE;
870683d3df6SStefano Zampini       }
871683d3df6SStefano Zampini     }
872683d3df6SStefano Zampini     size_schur = cum - n_I;
873683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
874683d3df6SStefano Zampini     if (cum == n) {
875683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
876683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
877683d3df6SStefano Zampini     } else {
8787dae84e0SHong Zhang       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
879683d3df6SStefano Zampini     }
880e62b6521Sstefano_zampini     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
881e62b6521Sstefano_zampini     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
882ca92afb2SStefano Zampini 
883ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
884367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
885ca92afb2SStefano Zampini     if (benign_n) {
886ca92afb2SStefano Zampini       Vec                    v;
887ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
888ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
8895cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
8905cbda25cSStefano Zampini       PetscInt               sizeA;
891ca92afb2SStefano Zampini 
892ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
893ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
894ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
895ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
896ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
8975cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
8985cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
899282d6408SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
9005cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9015cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
902ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
903ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
904ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
9055cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
906ca92afb2SStefano Zampini         PetscScalar    *array;
907ca92afb2SStefano Zampini         const PetscInt *idxs;
908ca92afb2SStefano Zampini         PetscInt       j,nz;
909ca92afb2SStefano Zampini 
910ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
911ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
912ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9135cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
914ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9155cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
9165cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
9175cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
918ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
91922db5ddcSStefano Zampini         for (j=0;j<size_schur;j++) {
92022db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
92122db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
92222db5ddcSStefano Zampini #else
92322db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
92422db5ddcSStefano Zampini #endif
92522db5ddcSStefano Zampini         }
926ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
927ca92afb2SStefano Zampini       }
9285cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9295cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
930ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
931ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
932ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
933ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
934ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
935ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
936ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
937ca92afb2SStefano Zampini     }
93811955456SStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);CHKERRQ(ierr);
9393301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9403301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
941883469d8SStefano Zampini 
94211955456SStefano Zampini     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
94311955456SStefano Zampini     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
94411955456SStefano Zampini 
945683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
946d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
947d47842beSStefano Zampini 
948f4f7d9d6SStefano Zampini     if (n_I) {
9490aa714b2SStefano Zampini       IS is_schur;
9505a05ddb0SStefano Zampini 
95111955456SStefano Zampini       if (use_cholesky) {
95288113c35SStefano Zampini         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
953883469d8SStefano Zampini       } else {
95488113c35SStefano Zampini         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
955883469d8SStefano Zampini       }
9568c09ecd8Sstefano_zampini       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
9578c09ecd8Sstefano_zampini 
958883469d8SStefano Zampini       /* subsets ordered last */
9595a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9605a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9615a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
962883469d8SStefano Zampini 
963883469d8SStefano Zampini       /* factorization step */
96411955456SStefano Zampini       if (use_cholesky) {
965883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
966be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
967be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
968be83ff47SStefano Zampini #endif
969883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
970a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
971883469d8SStefano Zampini       } else {
972883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
973be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
974be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
975be83ff47SStefano Zampini #endif
976883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
977883469d8SStefano Zampini       }
978ec968b6dSstefano_zampini       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
979883469d8SStefano Zampini 
980*7f9db97bSStefano Zampini       if (print_schurs) {
98111955456SStefano Zampini         PetscViewer viewer;
98211955456SStefano Zampini         char        filename[256];
98311955456SStefano Zampini         Mat         S;
98411955456SStefano Zampini         IS          is;
98511955456SStefano Zampini 
986*7f9db97bSStefano Zampini         ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
98711955456SStefano Zampini         ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
98811955456SStefano Zampini         ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
98911955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr);
99011955456SStefano Zampini         ierr = MatView(A,viewer);CHKERRQ(ierr);
99111955456SStefano Zampini         ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr);
99211955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)S,"S");CHKERRQ(ierr);
99311955456SStefano Zampini         ierr = MatView(S,viewer);CHKERRQ(ierr);
99411955456SStefano Zampini         ierr = MatDestroy(&S);CHKERRQ(ierr);
99511955456SStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);CHKERRQ(ierr);
99611955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)is,"I");CHKERRQ(ierr);
99711955456SStefano Zampini         ierr = ISView(is,viewer);CHKERRQ(ierr);
99811955456SStefano Zampini         ierr = ISDestroy(&is);CHKERRQ(ierr);
99911955456SStefano Zampini         ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);CHKERRQ(ierr);
100011955456SStefano Zampini         ierr = PetscObjectSetName((PetscObject)is,"B");CHKERRQ(ierr);
100111955456SStefano Zampini         ierr = ISView(is,viewer);CHKERRQ(ierr);
100211955456SStefano Zampini         ierr = ISDestroy(&is);CHKERRQ(ierr);
100311955456SStefano Zampini         ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
100411955456SStefano Zampini       }
100511955456SStefano Zampini 
1006883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
10077c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1008b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
1009b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
1010b3cb21ddSStefano Zampini 
1011d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
1012683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1013683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
101403dfb2d7SStefano Zampini       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
101503dfb2d7SStefano Zampini         reuse_solvers = factor_workaround = PETSC_FALSE;
101603dfb2d7SStefano Zampini 
1017df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
1018ca92afb2SStefano Zampini 
101972b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
1020ca92afb2SStefano Zampini       if (benign_n) {
10215cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1022ca92afb2SStefano Zampini         Vec         v;
10235cbda25cSStefano Zampini         PetscInt    sizeA;
1024ca92afb2SStefano Zampini 
1025ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
10265cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
10275cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1028ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1029ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1030ca92afb2SStefano Zampini #endif
1031ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1032ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1033ca92afb2SStefano Zampini #endif
10345cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10355cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1036ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
10375cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
103847484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
1039ca92afb2SStefano Zampini           const PetscInt *idxs;
1040ca92afb2SStefano Zampini           PetscInt       j,nz;
104147484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
1042ca92afb2SStefano Zampini 
10435cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
10445cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
10455cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1046ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1047ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1048282d6408SStefano Zampini           /* p0 dof (eliminated) is excluded from the sum */
1049282d6408SStefano Zampini           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1050ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10515cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1052ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
105347484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
105447484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
105547484b83SStefano Zampini           B_k = 1;
105647484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
105747484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
10585cbda25cSStefano Zampini           sum  = 1.;
105947484b83SStefano 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));
1060ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10615cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10625cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1063282d6408SStefano Zampini           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
10645cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10655cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1066ca92afb2SStefano Zampini         }
1067a7414863SStefano Zampini         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1068a7414863SStefano Zampini           PetscInt k,j;
1069a7414863SStefano Zampini           for (k=0;k<size_schur;k++) {
1070a7414863SStefano Zampini             for (j=k;j<size_schur;j++) {
1071a7414863SStefano Zampini               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1072a7414863SStefano Zampini             }
1073a7414863SStefano Zampini           }
1074a7414863SStefano Zampini         }
1075a7414863SStefano Zampini 
10765cbda25cSStefano Zampini         /* restore defaults */
10775cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
10785cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
10795cbda25cSStefano Zampini #endif
10805cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
10815cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
10825cbda25cSStefano Zampini #endif
10835cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10845cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1085ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1086ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1087ca92afb2SStefano Zampini       }
1088a3df083aSStefano Zampini       if (!reuse_solvers) {
1089a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1090a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1091a3df083aSStefano Zampini         }
1092a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
10935cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
10945cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1095a3df083aSStefano Zampini       }
1096df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1097be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1098683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1099166598c1SStefano Zampini       factor_workaround = PETSC_FALSE;
1100df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1101be83ff47SStefano Zampini     }
1102be83ff47SStefano Zampini 
1103be83ff47SStefano Zampini     if (reuse_solvers) {
1104a00504b5SStefano Zampini       Mat                A_II,Afake;
110553892102SStefano Zampini       Vec                vec1_B;
1106df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
11073462e049SStefano Zampini       PetscInt           n_R;
1108d5574798SStefano Zampini 
1109df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1110df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1111e28d306cSStefano Zampini       } else {
1112df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1113d62866d3SStefano Zampini       }
1114df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1115be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1116d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1117d5574798SStefano Zampini       msolv_ctx->F = F;
1118683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1119683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1120683d3df6SStefano Zampini       {
1121683d3df6SStefano Zampini         PetscScalar *array;
1122683d3df6SStefano Zampini         PetscInt    n;
1123683d3df6SStefano Zampini 
1124683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1125683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1126683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1127683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1128683d3df6SStefano Zampini       }
11293fc34f97SStefano Zampini       msolv_ctx->has_vertices = schur_has_vertices;
1130d62866d3SStefano Zampini 
1131d62866d3SStefano Zampini       /* interior solver */
1132683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1133d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1134d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1135d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1136df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1137df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1138d62866d3SStefano Zampini 
1139d62866d3SStefano Zampini       /* correction solver */
11403462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1141d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1142d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1143df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1144df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
114553892102SStefano Zampini 
114653892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
114753892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
114853892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
11493fc34f97SStefano Zampini       if (!schur_has_vertices) {
115053892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
115153892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
115253892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
115353892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1154683d3df6SStefano Zampini       } else {
1155683d3df6SStefano Zampini         IS              is_B_all;
1156683d3df6SStefano Zampini         const PetscInt* idxs;
1157683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1158683d3df6SStefano Zampini 
1159683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1160683d3df6SStefano Zampini         dual = size_schur - n_v;
1161683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1162683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1163683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1164683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1165683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1166683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1167683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1168683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1169683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1170683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1171683d3df6SStefano Zampini       }
11723462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
11733462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
11748d692618SStefano Zampini       ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11758d692618SStefano Zampini       ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11763462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
11773462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1178683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1179ca92afb2SStefano Zampini 
1180ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1181ca92afb2SStefano Zampini       if (benign_n) {
11825cbda25cSStefano Zampini         PetscScalar *array;
11835cbda25cSStefano Zampini 
1184ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1185ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1186ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
11875cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
1188282d6408SStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
11895cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11905cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
11915cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11925cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1193ca92afb2SStefano Zampini       }
1194ada6e2d7SStefano Zampini     } else {
1195ada6e2d7SStefano Zampini       if (sub_schurs->reuse_solver) {
1196ada6e2d7SStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1197ada6e2d7SStefano Zampini       }
1198ada6e2d7SStefano Zampini       ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1199d5574798SStefano Zampini     }
120008122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
120153892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
12025db18549SStefano Zampini 
1203be83ff47SStefano Zampini     /* Work arrays */
1204be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
12055a94c880SStefano Zampini 
120604270d10SStefano Zampini     /* matrices for deluxe scaling and adaptive selection */
120712d906b1SStefano Zampini     if (compute_Stilda) {
12085a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
1209b449dd42SStefano Zampini         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
12105a94c880SStefano Zampini       }
12119d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1212b449dd42SStefano Zampini         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
12135a94c880SStefano Zampini       }
121412d906b1SStefano Zampini     }
1215d2627357SStefano Zampini 
1216be83ff47SStefano Zampini     /* S_Ej_all */
1217be83ff47SStefano Zampini     cum = cum2 = 0;
1218be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
121965d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
122065d8bf0aSStefano Zampini       PetscInt j;
122165d8bf0aSStefano Zampini 
12225a95e1ceSStefano Zampini       /* get S_E */
1223b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1224683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1225be83ff47SStefano Zampini         PetscInt k;
1226be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1227be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1228be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1229683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1230be83ff47SStefano Zampini           }
1231be83ff47SStefano Zampini         }
123206a4e24aSStefano Zampini       } else { /* just copy to workspace */
1233be83ff47SStefano Zampini         PetscInt k;
1234be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1235be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1236be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1237be83ff47SStefano Zampini           }
1238be83ff47SStefano Zampini         }
12399087bf02SStefano Zampini       }
12405a95e1ceSStefano Zampini       /* insert S_E values */
1241be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1242b7ab4a40SStefano Zampini       if (sub_schurs->change) {
12438760537fSStefano Zampini         Mat change_sub,SEj,T;
124472b8c272SStefano Zampini 
124572b8c272SStefano Zampini         /* change basis */
124672b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
124772b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12488760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
12498760537fSStefano Zampini           Mat T2;
12508760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
12518760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
12528760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
125372b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
12548760537fSStefano Zampini         } else {
12558760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
125672b8c272SStefano Zampini         }
125772b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
125872b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
1259afa0f562SStefano Zampini         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
126072b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
126172b8c272SStefano Zampini       }
12629d54b7f4SStefano Zampini       if (deluxe) {
12635a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1264683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1265862806e4SStefano Zampini         if (compute_Stilda) {
1266f4f7d9d6SStefano Zampini           PetscInt k;
1267f4f7d9d6SStefano Zampini 
126808122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
126908122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1270f4f7d9d6SStefano Zampini           if (use_potr) {
12715a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12722972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
12735a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
12742972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
12756c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
12766c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
12776c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
12786c3e6151SStefano Zampini               }
12796c3e6151SStefano Zampini             }
1280f4f7d9d6SStefano Zampini           } else if (use_sytr) {
1281f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1282f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1283f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,work,&B_N,pivots,Bwork,&B_ierr));
1284f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1285f4f7d9d6SStefano Zampini             for (k=0;k<subset_size;k++) {
1286f4f7d9d6SStefano Zampini               for (j=k;j<subset_size;j++) {
1287f4f7d9d6SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
1288f4f7d9d6SStefano Zampini               }
1289f4f7d9d6SStefano Zampini             }
1290d6462365SStefano Zampini           } else {
1291d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1292d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1293d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1294d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
12952972d61bSStefano Zampini           }
129608122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
12975a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12989087bf02SStefano Zampini         }
12999d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
13009d54b7f4SStefano Zampini         Mat         SEj;
13019d54b7f4SStefano Zampini         Vec         D;
13029d54b7f4SStefano Zampini         PetscScalar *array;
13039d54b7f4SStefano Zampini 
13049d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
13059d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
13069d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
13079d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1308f17d2ae1SStefano Zampini         ierr = VecShift(D,-1.);CHKERRQ(ierr);
13099d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
13109d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
13119d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
13129d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13139d54b7f4SStefano Zampini       }
1314be83ff47SStefano Zampini       cum += subset_size;
1315be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1316be83ff47SStefano Zampini     }
1317be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
13184a6c6b0dSStefano Zampini 
1319df4d28bfSStefano Zampini     if (solver_S) {
13207c2f51b8SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
13214a6c6b0dSStefano Zampini     }
1322683d3df6SStefano Zampini 
1323683d3df6SStefano Zampini     schur_factor = NULL;
132445951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1325683d3df6SStefano Zampini 
13269d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
13274a6c6b0dSStefano Zampini         PetscInt j;
13284a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
13295a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13304a6c6b0dSStefano Zampini       } else {
1331683d3df6SStefano Zampini         Mat S_all_inv=NULL;
13323fc34f97SStefano Zampini         if (solver_S) {
1333683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1334683d3df6SStefano 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 */
13353fc34f97SStefano Zampini           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1336683d3df6SStefano Zampini             PetscScalar *data;
1337683d3df6SStefano Zampini             PetscInt     nd = 0;
13386dba178dSStefano Zampini 
1339f4f7d9d6SStefano Zampini             if (!use_potr) {
1340683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1341683d3df6SStefano Zampini             }
13427c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1343683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1344683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1345683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1346683d3df6SStefano Zampini             }
13473fc34f97SStefano Zampini 
13483fc34f97SStefano Zampini             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
13493fc34f97SStefano Zampini             if (schur_has_vertices) {
13503fc34f97SStefano Zampini               Mat          M;
13513fc34f97SStefano Zampini               PetscScalar *tdata;
13523fc34f97SStefano Zampini               PetscInt     nv = 0, news;
13533fc34f97SStefano Zampini 
13543fc34f97SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
13553fc34f97SStefano Zampini               news = size_active_schur + nv;
13563fc34f97SStefano Zampini               ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr);
1357683d3df6SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13583fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
13593fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));CHKERRQ(ierr);
13603fc34f97SStefano Zampini               }
13613fc34f97SStefano Zampini               for (i=0;i<nv;i++) {
13623fc34f97SStefano Zampini                 PetscInt k = i+size_active_schur;
13633fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));CHKERRQ(ierr);
13643fc34f97SStefano Zampini               }
13653fc34f97SStefano Zampini 
13663fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr);
13673fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
13683fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
13693fc34f97SStefano Zampini               /* save the factors */
13703fc34f97SStefano Zampini               cum = 0;
13713fc34f97SStefano Zampini               ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr);
13723fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13733fc34f97SStefano Zampini                 ierr = PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1374683d3df6SStefano Zampini                 cum += size_active_schur - i;
1375683d3df6SStefano Zampini               }
13763fc34f97SStefano Zampini               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
13773fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
13783fc34f97SStefano Zampini               /* move back just the active dofs to the Schur complement */
13793fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13803fc34f97SStefano Zampini                 ierr = PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));CHKERRQ(ierr);
13813fc34f97SStefano Zampini               }
13823fc34f97SStefano Zampini               ierr = PetscFree(tdata);CHKERRQ(ierr);
13833fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13843fc34f97SStefano Zampini             } else { /* we can factorize and invert just the activedofs */
13853fc34f97SStefano Zampini               Mat         M;
13865002105bSStefano Zampini               PetscScalar *aux;
13873fc34f97SStefano Zampini 
13885002105bSStefano Zampini               ierr = PetscMalloc1(nd,&aux);CHKERRQ(ierr);
13895002105bSStefano Zampini               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
13903fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr);
13913fc34f97SStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
13923fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
13933fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
13943fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
13953fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13965002105bSStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);CHKERRQ(ierr);
13975002105bSStefano Zampini               ierr = MatZeroEntries(M);CHKERRQ(ierr);
13985002105bSStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13995002105bSStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);CHKERRQ(ierr);
14005002105bSStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
14015002105bSStefano Zampini               ierr = MatZeroEntries(M);CHKERRQ(ierr);
14025002105bSStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
14035002105bSStefano Zampini               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
14045002105bSStefano Zampini               ierr = PetscFree(aux);CHKERRQ(ierr);
14053fc34f97SStefano Zampini             }
1406683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
14073fc34f97SStefano Zampini           } else { /* use MatFactor calls to invert S */
14086dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
14097c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1410683d3df6SStefano Zampini           }
1411df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1412683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1413683d3df6SStefano Zampini           S_all_inv = S_all;
1414683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1415be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1416be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1417f4f7d9d6SStefano Zampini           if (use_potr) {
1418be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1419be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1420be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1421be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1422f4f7d9d6SStefano Zampini           } else if (use_sytr) {
1423f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1424f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1425f4f7d9d6SStefano Zampini             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1426f4f7d9d6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1427d6462365SStefano Zampini           } else {
1428d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1429d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1430d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1431d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1432be83ff47SStefano Zampini           }
1433be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1434683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1435be83ff47SStefano Zampini         }
1436be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1437be83ff47SStefano Zampini         cum = cum2 = 0;
1438683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1439be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1440be83ff47SStefano Zampini           PetscInt j;
1441862806e4SStefano Zampini 
1442862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1443be83ff47SStefano Zampini           /* get (St^-1)_E */
144472b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
144506a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1446a0b0af32SStefano Zampini           if (S_lower_triangular) {
1447be83ff47SStefano Zampini             PetscInt k;
1448b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1449be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1450be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1451be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
14526c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1453be83ff47SStefano Zampini                 }
1454be83ff47SStefano Zampini               }
145572b8c272SStefano Zampini             } else {
145672b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
145772b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
145872b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
145972b8c272SStefano Zampini                 }
146072b8c272SStefano Zampini               }
146172b8c272SStefano Zampini             }
146272b8c272SStefano Zampini           } else {
1463be83ff47SStefano Zampini             PetscInt k;
1464be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1465be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1466be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1467be83ff47SStefano Zampini               }
1468be83ff47SStefano Zampini             }
1469be83ff47SStefano Zampini           }
1470b7ab4a40SStefano Zampini           if (sub_schurs->change) {
14718760537fSStefano Zampini             Mat change_sub,SEj,T;
147272b8c272SStefano Zampini 
147372b8c272SStefano Zampini             /* change basis */
147472b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
147572b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
14768760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
14778760537fSStefano Zampini               Mat T2;
14788760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
14798760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
148072b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
14818760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
14828760537fSStefano Zampini             } else {
14838760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
148472b8c272SStefano Zampini             }
14858760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
14868760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
148772b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1488afa0f562SStefano Zampini             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
148972b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
149072b8c272SStefano Zampini           }
1491be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
14925a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1493be83ff47SStefano Zampini           cum += subset_size;
1494be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1495883469d8SStefano Zampini         }
1496683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1497df4d28bfSStefano Zampini         if (solver_S) {
14983fc34f97SStefano Zampini           if (schur_has_vertices) {
14993fc34f97SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
15003fc34f97SStefano Zampini           } else {
15017c2f51b8SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
15025db18549SStefano Zampini           }
15033fc34f97SStefano Zampini         }
1504683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1505683d3df6SStefano Zampini       }
1506683d3df6SStefano Zampini 
15073fc34f97SStefano Zampini       /* move back factors if needed */
15083fc34f97SStefano Zampini       if (schur_has_vertices) {
1509683d3df6SStefano Zampini         Mat      S_tmp;
15103fc34f97SStefano Zampini         PetscInt nd = 0;
1511683d3df6SStefano Zampini 
15126542c05fSStefano Zampini         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
15137c2f51b8SStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1514f4f7d9d6SStefano Zampini         if (use_potr) {
1515683d3df6SStefano Zampini           PetscScalar *data;
1516683d3df6SStefano Zampini 
1517683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
15183fc34f97SStefano Zampini           ierr = PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1519683d3df6SStefano Zampini 
1520683d3df6SStefano Zampini           if (S_lower_triangular) {
1521683d3df6SStefano Zampini             cum = 0;
1522683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1523683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1524683d3df6SStefano Zampini               cum += size_active_schur-i;
1525683d3df6SStefano Zampini 	    }
1526683d3df6SStefano Zampini           } else {
1527683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1528683d3df6SStefano Zampini           }
1529683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1530683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1531683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1532683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1533683d3df6SStefano Zampini 	    }
1534683d3df6SStefano Zampini           }
15356dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1536683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1537683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
15386c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1539683d3df6SStefano Zampini           }
1540683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
15416542c05fSStefano Zampini         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
15423fc34f97SStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
15439087bf02SStefano Zampini       }
1544367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1545367aa537SStefano Zampini       PetscScalar *data;
1546367aa537SStefano Zampini       PetscInt    nd = 0;
1547367aa537SStefano Zampini 
1548367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1549367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1550367aa537SStefano Zampini       }
15517c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1552367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1553367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1554367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1555367aa537SStefano Zampini       }
1556367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1557367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
15586c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1559367aa537SStefano Zampini       }
1560367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
15613fc34f97SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
15624a6c6b0dSStefano Zampini     }
15634a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1564683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
15659d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
15664a6c6b0dSStefano Zampini   }
15675a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1568be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1569a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1570a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1571a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1572a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1573a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
15745db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15755db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15765a95e1ceSStefano Zampini   if (compute_Stilda) {
15775a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15785a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15799d54b7f4SStefano Zampini     if (deluxe) {
15805a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15815a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158208122e43SStefano Zampini     }
15839d54b7f4SStefano Zampini   }
15845db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
15855db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
15863927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
15873927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15885a95e1ceSStefano Zampini 
15895db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
15907dae84e0SHong Zhang   ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
159141fd5443SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
159241fd5443SStefano Zampini     ierr = MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
159341fd5443SStefano Zampini   } else {
159404270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
159541fd5443SStefano Zampini   }
159608122e43SStefano Zampini 
1597f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
15985a95e1ceSStefano Zampini   if (compute_Stilda) {
15995a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1600a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
160104270d10SStefano Zampini     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
160204270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
16039d54b7f4SStefano Zampini     if (deluxe) {
16045a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
160508122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
160604270d10SStefano Zampini       ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
160704270d10SStefano Zampini       ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
16089d54b7f4SStefano Zampini     } else {
16099d54b7f4SStefano Zampini       PetscScalar *array;
16109d54b7f4SStefano Zampini       PetscInt    cum;
16119d54b7f4SStefano Zampini 
16129d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
16139d54b7f4SStefano Zampini       cum = 0;
16149d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
16159d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
16169d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
16179d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1618f4f7d9d6SStefano Zampini         if (use_potr) {
16199d54b7f4SStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
16209d54b7f4SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
16219d54b7f4SStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
16229d54b7f4SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1623f4f7d9d6SStefano Zampini         } else if (use_sytr) {
1624f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1625f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1626f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1627f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1628f4f7d9d6SStefano Zampini         } else {
1629f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1630f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1631f4f7d9d6SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1632f4f7d9d6SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1633f4f7d9d6SStefano Zampini         }
16349d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
16359d54b7f4SStefano Zampini         cum += subset_size*subset_size;
16369d54b7f4SStefano Zampini       }
16379d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
16389d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
16399d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
16409d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
16419d54b7f4SStefano Zampini     }
164208122e43SStefano Zampini   }
164304270d10SStefano Zampini   ierr = MatDestroySubMatrices(1,&submats);CHKERRQ(ierr);
1644*7f9db97bSStefano Zampini   if (print_schurs) {
164511955456SStefano Zampini     PetscViewer viewer;
164611955456SStefano Zampini     char        filename[256];
164711955456SStefano Zampini     PetscInt    cum;
164811955456SStefano Zampini 
1649*7f9db97bSStefano Zampini     ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_mats_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
165011955456SStefano Zampini     ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
165111955456SStefano Zampini     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
165211955456SStefano Zampini     if (sub_schurs->S_Ej_all) {
165311955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");CHKERRQ(ierr);
165411955456SStefano Zampini       ierr = MatView(sub_schurs->S_Ej_all,viewer);CHKERRQ(ierr);
165511955456SStefano Zampini     }
165611955456SStefano Zampini     if (sub_schurs->sum_S_Ej_all) {
165711955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");CHKERRQ(ierr);
165811955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_all,viewer);CHKERRQ(ierr);
165911955456SStefano Zampini     }
166011955456SStefano Zampini     if (sub_schurs->sum_S_Ej_inv_all) {
166111955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");CHKERRQ(ierr);
166211955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_inv_all,viewer);CHKERRQ(ierr);
166311955456SStefano Zampini     }
166411955456SStefano Zampini     if (sub_schurs->sum_S_Ej_tilda_all) {
166511955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");CHKERRQ(ierr);
166611955456SStefano Zampini       ierr = MatView(sub_schurs->sum_S_Ej_tilda_all,viewer);CHKERRQ(ierr);
166711955456SStefano Zampini     }
166811955456SStefano Zampini     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
166911955456SStefano Zampini       IS   is;
167011955456SStefano Zampini       char name[16];
167111955456SStefano Zampini 
1672*7f9db97bSStefano Zampini       ierr = PetscSNPrintf(name,sizeof(name),"IE%D",i);CHKERRQ(ierr);
167311955456SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
167411955456SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);CHKERRQ(ierr);
167511955456SStefano Zampini       ierr = PetscObjectSetName((PetscObject)is,name);CHKERRQ(ierr);
167611955456SStefano Zampini       ierr = ISView(is,viewer);CHKERRQ(ierr);
167711955456SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
167811955456SStefano Zampini       cum += subset_size;
167911955456SStefano Zampini     }
168011955456SStefano Zampini     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
168111955456SStefano Zampini   }
16823202ece2SStefano Zampini 
16835a95e1ceSStefano Zampini   /* free workspace */
16845a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
168506a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1686a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
16873202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1688dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
16895a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1690b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1691b1b3d7a2SStefano Zampini }
1692b1b3d7a2SStefano Zampini 
169388113c35SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1694b1b3d7a2SStefano Zampini {
16959bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
16965a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
169788113c35SStefano Zampini   PetscBool       is_sorted,ispetsc;
1698b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1699b1b3d7a2SStefano Zampini 
1700b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1701b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
17026c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1703b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
17046c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1705b1b3d7a2SStefano Zampini 
1706b1b3d7a2SStefano Zampini   /* reset any previous data */
1707b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1708b1b3d7a2SStefano Zampini 
17095a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
17109bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1711b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
171208122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1713b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1714b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
17158b6046baSStefano Zampini     if (copycc) {
17168b6046baSStefano Zampini       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
17178b6046baSStefano Zampini     } else {
1718c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1719b1b3d7a2SStefano Zampini       all_cc[i] = faces[i];
1720b1b3d7a2SStefano Zampini     }
17218b6046baSStefano Zampini   }
1722b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
17238b6046baSStefano Zampini     if (copycc) {
17248b6046baSStefano Zampini       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
17258b6046baSStefano Zampini     } else {
1726c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1727b1b3d7a2SStefano Zampini       all_cc[n_faces+i] = edges[i];
17288b6046baSStefano Zampini     }
172908122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1730b1b3d7a2SStefano Zampini   }
1731c8272957SStefano Zampini   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1732c8272957SStefano Zampini   sub_schurs->is_vertices = vertices;
1733c8272957SStefano Zampini   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1734d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1735d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1736b1b3d7a2SStefano Zampini 
1737df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
173888113c35SStefano Zampini   ierr = PetscStrallocpy(prefix,&sub_schurs->prefix);CHKERRQ(ierr);
1739883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
174088113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);CHKERRQ(ierr);
174188113c35SStefano Zampini #elif defined(PETSC_HAVE_MKL_PARDISO)
174288113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);CHKERRQ(ierr);
174388113c35SStefano Zampini #else
174488113c35SStefano Zampini   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);CHKERRQ(ierr);
1745df4d28bfSStefano Zampini #endif
174688113c35SStefano Zampini #if defined(PETSC_USE_COMPLEX)
174788113c35SStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
174888113c35SStefano Zampini #else
174988113c35SStefano Zampini   sub_schurs->is_hermitian = PETSC_TRUE;
1750883469d8SStefano Zampini #endif
175188113c35SStefano Zampini   sub_schurs->is_posdef    = PETSC_TRUE;
175211955456SStefano Zampini   sub_schurs->is_symmetric = PETSC_TRUE;
1753*7f9db97bSStefano Zampini   sub_schurs->debug        = PETSC_FALSE;
175488113c35SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
175588113c35SStefano 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);
175611955456SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);CHKERRQ(ierr);
175788113c35SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
175888113c35SStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
1759*7f9db97bSStefano Zampini   ierr = PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);CHKERRQ(ierr);
176088113c35SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
176188113c35SStefano Zampini   ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);CHKERRQ(ierr);
1762168c3995SStefano Zampini   sub_schurs->schur_explicit = (PetscBool)!ispetsc;
1763b1b3d7a2SStefano Zampini 
176411955456SStefano Zampini   /* for reals, symmetric and hermitian are synonims */
176511955456SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
176611955456SStefano Zampini   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
176711955456SStefano Zampini   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
176811955456SStefano Zampini #endif
176911955456SStefano Zampini 
1770b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1771b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1772b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1773b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
17745db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
17755db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
17765db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
17775db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
17785a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1779b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1780b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1781b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
178208122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1783b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1784b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1785b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1786b1b3d7a2SStefano Zampini }
1787b1b3d7a2SStefano Zampini 
178834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
178934a97f8cSStefano Zampini {
179034a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
179134a97f8cSStefano Zampini   PetscErrorCode  ierr;
179234a97f8cSStefano Zampini 
179334a97f8cSStefano Zampini   PetscFunctionBegin;
179434a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
17955ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
179634a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
179734a97f8cSStefano Zampini   PetscFunctionReturn(0);
179834a97f8cSStefano Zampini }
179934a97f8cSStefano Zampini 
180034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
180134a97f8cSStefano Zampini {
180234a97f8cSStefano Zampini   PetscInt       i;
180334a97f8cSStefano Zampini   PetscErrorCode ierr;
180434a97f8cSStefano Zampini 
180534a97f8cSStefano Zampini   PetscFunctionBegin;
1806aea80f77Sstefano_zampini   if (!sub_schurs) PetscFunctionReturn(0);
180788113c35SStefano Zampini   ierr = PetscFree(sub_schurs->prefix);CHKERRQ(ierr);
18081e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1809b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1810b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1811b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
18125db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
18135db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
181441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
181541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
181608122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1817a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
18185db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1819d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1820d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
182108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
182234a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1823b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
182434a97f8cSStefano Zampini   }
18255ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1826b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
18273dc780c3SStefano Zampini   }
1828df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1829df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1830d62866d3SStefano Zampini   }
1831df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
183272b8c272SStefano Zampini   if (sub_schurs->change) {
183372b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
183472b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1835b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
183672b8c272SStefano Zampini     }
183772b8c272SStefano Zampini   }
183872b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1839b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
184034a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
184134a97f8cSStefano Zampini   PetscFunctionReturn(0);
184234a97f8cSStefano Zampini }
184334a97f8cSStefano Zampini 
1844aea80f77Sstefano_zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1845aea80f77Sstefano_zampini {
1846aea80f77Sstefano_zampini   PetscErrorCode ierr;
1847aea80f77Sstefano_zampini 
1848aea80f77Sstefano_zampini   PetscFunctionBegin;
1849aea80f77Sstefano_zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1850aea80f77Sstefano_zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1851aea80f77Sstefano_zampini   PetscFunctionReturn(0);
1852aea80f77Sstefano_zampini }
1853aea80f77Sstefano_zampini 
18542a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
185534a97f8cSStefano Zampini {
185634a97f8cSStefano Zampini   PetscInt       i,j,n;
185734a97f8cSStefano Zampini   PetscErrorCode ierr;
185834a97f8cSStefano Zampini 
185934a97f8cSStefano Zampini   PetscFunctionBegin;
186034a97f8cSStefano Zampini   n = 0;
186134a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
186234a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
186334a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
186434a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
186534a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
186634a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
186734a97f8cSStefano Zampini         queue_tip[n] = dof;
186834a97f8cSStefano Zampini         n++;
186934a97f8cSStefano Zampini       }
187034a97f8cSStefano Zampini     }
187134a97f8cSStefano Zampini   }
187234a97f8cSStefano Zampini   *n_added = n;
187334a97f8cSStefano Zampini   PetscFunctionReturn(0);
187434a97f8cSStefano Zampini }
1875