xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 04270d10403b2d22cab08fba741a972f7ff8fba2)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
308122e43SStefano Zampini #include <petscblaslapack.h>
434a97f8cSStefano Zampini 
53fc34f97SStefano Zampini /* this is declared in dense.h */
63fc34f97SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
73202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
85ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
9df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
10df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
11d62866d3SStefano Zampini 
12ca92afb2SStefano Zampini /* if v2 is not present, correction is done in-place */
135cbda25cSStefano Zampini PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
14ca92afb2SStefano Zampini {
15ca92afb2SStefano Zampini   PetscScalar    *array;
16ca92afb2SStefano Zampini   PetscScalar    *array2;
17ca92afb2SStefano Zampini   PetscErrorCode ierr;
18ca92afb2SStefano Zampini 
19ca92afb2SStefano Zampini   PetscFunctionBegin;
20ca92afb2SStefano Zampini   if (!ctx->benign_n) PetscFunctionReturn(0);
215cbda25cSStefano Zampini   if (sol && full) {
225cbda25cSStefano Zampini     PetscInt n_I,size_schur;
235cbda25cSStefano Zampini 
245cbda25cSStefano Zampini     /* get sizes */
25282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
265cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
275cbda25cSStefano Zampini     n_I = n_I - size_schur;
285cbda25cSStefano Zampini     /* get schur sol from array */
295cbda25cSStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
305cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
315cbda25cSStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
325cbda25cSStefano Zampini     /* apply interior sol correction */
33282d6408SStefano Zampini     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
345cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
355cbda25cSStefano Zampini     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
365cbda25cSStefano Zampini   }
37ca92afb2SStefano Zampini   if (v2) {
38ca92afb2SStefano Zampini     PetscInt nl;
39ca92afb2SStefano Zampini 
40ca92afb2SStefano Zampini     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
41ca92afb2SStefano Zampini     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
42ca92afb2SStefano Zampini     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
43ca92afb2SStefano Zampini     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
44ca92afb2SStefano Zampini   } else {
45ca92afb2SStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
46ca92afb2SStefano Zampini     array2 = array;
47ca92afb2SStefano Zampini   }
48ca92afb2SStefano Zampini   if (!sol) { /* change rhs */
49ca92afb2SStefano Zampini     PetscInt n;
50ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
51ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
52ca92afb2SStefano Zampini       const PetscInt *cols;
53ca92afb2SStefano Zampini       PetscInt       nz,i;
54ca92afb2SStefano Zampini 
55ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
56ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
57ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
5822db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
5922db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
6022db5ddcSStefano Zampini #else
61ca92afb2SStefano Zampini       sum = -sum/nz;
6222db5ddcSStefano Zampini #endif
63ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
64ca92afb2SStefano Zampini       ctx->benign_save_vals[n] = array2[cols[nz-1]];
65ca92afb2SStefano Zampini       array2[cols[nz-1]] = sum;
66ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
67ca92afb2SStefano Zampini     }
68ca92afb2SStefano Zampini   } else {
69ca92afb2SStefano Zampini     PetscInt n;
70ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
71ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
72ca92afb2SStefano Zampini       const PetscInt *cols;
73ca92afb2SStefano Zampini       PetscInt       nz,i;
74ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
75ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
76ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
7722db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
7822db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
7922db5ddcSStefano Zampini #else
80ca92afb2SStefano Zampini       sum = -sum/nz;
8122db5ddcSStefano Zampini #endif
82ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
83ca92afb2SStefano Zampini       array2[cols[nz-1]] = ctx->benign_save_vals[n];
84ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
85ca92afb2SStefano Zampini     }
86ca92afb2SStefano Zampini   }
87ca92afb2SStefano Zampini   if (v2) {
88ca92afb2SStefano Zampini     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
8972b8c272SStefano Zampini     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
90ca92afb2SStefano Zampini   } else {
91ca92afb2SStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
92ca92afb2SStefano Zampini   }
935cbda25cSStefano Zampini   if (!sol && full) {
945cbda25cSStefano Zampini     Vec      usedv;
955cbda25cSStefano Zampini     PetscInt n_I,size_schur;
965cbda25cSStefano Zampini 
975cbda25cSStefano Zampini     /* get sizes */
98282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
995cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
1005cbda25cSStefano Zampini     n_I = n_I - size_schur;
1015cbda25cSStefano Zampini     /* compute schur rhs correction */
1025cbda25cSStefano Zampini     if (v2) {
1035cbda25cSStefano Zampini       usedv = v2;
1045cbda25cSStefano Zampini     } else {
1055cbda25cSStefano Zampini       usedv = v;
1065cbda25cSStefano Zampini     }
1075cbda25cSStefano Zampini     /* apply schur rhs correction */
1085cbda25cSStefano Zampini     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
1095cbda25cSStefano Zampini     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1105cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
1115cbda25cSStefano Zampini     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
112282d6408SStefano Zampini     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1135cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1145cbda25cSStefano Zampini   }
115ca92afb2SStefano Zampini   PetscFunctionReturn(0);
116ca92afb2SStefano Zampini }
117ca92afb2SStefano Zampini 
118df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
119d62866d3SStefano Zampini {
120df4d28bfSStefano Zampini   PCBDDCReuseSolvers ctx;
121683d3df6SStefano Zampini   PetscBool          copy = PETSC_FALSE;
122d62866d3SStefano Zampini   PetscErrorCode     ierr;
123d62866d3SStefano Zampini 
124d62866d3SStefano Zampini   PetscFunctionBegin;
125d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
126683d3df6SStefano Zampini   if (full) {
127d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
128d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
129d62866d3SStefano Zampini #endif
1305cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1315cbda25cSStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
1325cbda25cSStefano Zampini #endif
133683d3df6SStefano Zampini     copy = ctx->has_vertices;
134d4933d67SStefano Zampini   } else { /* interior solver */
1356dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
136683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
1376dba178dSStefano Zampini #endif
138d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
139d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
140d4933d67SStefano Zampini #endif
141683d3df6SStefano Zampini     copy = PETSC_TRUE;
142683d3df6SStefano Zampini   }
143683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
144683d3df6SStefano Zampini   if (copy) {
145ca92afb2SStefano Zampini     PetscInt    n;
146df4d28bfSStefano Zampini     PetscScalar *array,*array_solver;
147ca92afb2SStefano Zampini 
148ca92afb2SStefano Zampini     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
149683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
150df4d28bfSStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
151df4d28bfSStefano Zampini     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
152df4d28bfSStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
153683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
154683d3df6SStefano Zampini 
1555cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
156683d3df6SStefano Zampini     if (transpose) {
157683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
158683d3df6SStefano Zampini     } else {
159683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
160683d3df6SStefano Zampini     }
1615cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
162683d3df6SStefano Zampini 
163683d3df6SStefano Zampini     /* get back data to caller worskpace */
164df4d28bfSStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
165683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
166df4d28bfSStefano Zampini     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
167683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
168df4d28bfSStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
169683d3df6SStefano Zampini   } else {
170ca92afb2SStefano Zampini     if (ctx->benign_n) {
1715cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
172ca92afb2SStefano Zampini       if (transpose) {
173ca92afb2SStefano Zampini         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
174ca92afb2SStefano Zampini       } else {
175ca92afb2SStefano Zampini         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
176ca92afb2SStefano Zampini       }
1775cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
178ca92afb2SStefano Zampini     } else {
179e28d306cSStefano Zampini       if (transpose) {
180e28d306cSStefano Zampini         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
181e28d306cSStefano Zampini       } else {
1826816873aSStefano Zampini         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
183e28d306cSStefano Zampini       }
184683d3df6SStefano Zampini     }
185ca92afb2SStefano Zampini   }
1865cbda25cSStefano Zampini   /* restore defaults */
1875cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1885cbda25cSStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
1895cbda25cSStefano Zampini #endif
190d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
191d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
192d4933d67SStefano Zampini #endif
193d62866d3SStefano Zampini   PetscFunctionReturn(0);
194d62866d3SStefano Zampini }
195d62866d3SStefano Zampini 
196df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
197e28d306cSStefano Zampini {
198e28d306cSStefano Zampini   PetscErrorCode   ierr;
199e28d306cSStefano Zampini 
200e28d306cSStefano Zampini   PetscFunctionBegin;
201df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
202e28d306cSStefano Zampini   PetscFunctionReturn(0);
203e28d306cSStefano Zampini }
204e28d306cSStefano Zampini 
205df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
206e28d306cSStefano Zampini {
207e28d306cSStefano Zampini   PetscErrorCode   ierr;
208e28d306cSStefano Zampini 
209e28d306cSStefano Zampini   PetscFunctionBegin;
210df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
211683d3df6SStefano Zampini   PetscFunctionReturn(0);
212683d3df6SStefano Zampini }
213683d3df6SStefano Zampini 
214df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
215683d3df6SStefano Zampini {
216683d3df6SStefano Zampini   PetscErrorCode   ierr;
217683d3df6SStefano Zampini 
218683d3df6SStefano Zampini   PetscFunctionBegin;
219df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
220683d3df6SStefano Zampini   PetscFunctionReturn(0);
221683d3df6SStefano Zampini }
222683d3df6SStefano Zampini 
223df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
224683d3df6SStefano Zampini {
225683d3df6SStefano Zampini   PetscErrorCode   ierr;
226683d3df6SStefano Zampini 
227683d3df6SStefano Zampini   PetscFunctionBegin;
228df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
229e28d306cSStefano Zampini   PetscFunctionReturn(0);
230e28d306cSStefano Zampini }
231e28d306cSStefano Zampini 
232df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
233d62866d3SStefano Zampini {
234ca92afb2SStefano Zampini   PetscInt       i;
235d62866d3SStefano Zampini   PetscErrorCode ierr;
236d62866d3SStefano Zampini 
237d62866d3SStefano Zampini   PetscFunctionBegin;
238d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
239d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
240d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
241d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
242d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
24353892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
24453892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
245d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
24653892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
24753892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
248ca92afb2SStefano Zampini   for (i=0;i<reuse->benign_n;i++) {
249ca92afb2SStefano Zampini     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
250ca92afb2SStefano Zampini   }
251ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
252ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
2535cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
2545cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
2555cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
2565cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
257d62866d3SStefano Zampini   PetscFunctionReturn(0);
258d62866d3SStefano Zampini }
259d5574798SStefano Zampini 
2605ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
2613202ece2SStefano Zampini {
2623202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
2633202ece2SStefano Zampini   KSP            ksp;
2643202ece2SStefano Zampini   PC             pc;
2653202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
2663202ece2SStefano Zampini   PetscReal      fill = 2.0;
267f11841e3SStefano Zampini   PetscInt       n_I;
2683202ece2SStefano Zampini   PetscMPIInt    size;
2693202ece2SStefano Zampini   PetscErrorCode ierr;
2703202ece2SStefano Zampini 
2713202ece2SStefano Zampini   PetscFunctionBegin;
2723202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
2736c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
274f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
275f11841e3SStefano Zampini     PetscBool Sdense;
276f11841e3SStefano Zampini 
277f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
2786c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
279f11841e3SStefano Zampini   }
2803202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
2813202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
2823202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
2833202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
2843202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
2853202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
2863202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
2873202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
288f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
289f11841e3SStefano Zampini   if (n_I) {
2903202ece2SStefano Zampini     if (!Bdense) {
2913202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
2923202ece2SStefano Zampini     } else {
2933202ece2SStefano Zampini       Bd = B;
2943202ece2SStefano Zampini     }
2953202ece2SStefano Zampini 
2963202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
2973202ece2SStefano Zampini       Mat fact;
2983202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
2993202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
3003202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
3013202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
3023202ece2SStefano Zampini     } else {
30307b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
30407b1e237SStefano Zampini 
30507b1e237SStefano Zampini       if (ex) {
3063202ece2SStefano Zampini         Mat Ainvd;
3073202ece2SStefano Zampini 
3083202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
3093202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
3103202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
31107b1e237SStefano Zampini       } else {
31207b1e237SStefano Zampini         Vec         sol,rhs;
31307b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
31407b1e237SStefano Zampini         PetscInt    i,nrhs,n;
31507b1e237SStefano Zampini 
31607b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
31707b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
31807b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
31907b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
32007b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
32107b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
32207b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
32307b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
32407b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
32507b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
32607b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
32707b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
32807b1e237SStefano Zampini         }
32907b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
33007b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
33107b1e237SStefano Zampini       }
3323202ece2SStefano Zampini     }
3335ec10c6aSStefano Zampini     if (!Bdense & !issym) {
3343202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3353202ece2SStefano Zampini     }
3365ec10c6aSStefano Zampini 
3375ec10c6aSStefano Zampini     if (!issym) {
3383202ece2SStefano Zampini       if (!Cdense) {
3393202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
3403202ece2SStefano Zampini       } else {
3413202ece2SStefano Zampini         Cd = C;
3423202ece2SStefano Zampini       }
3435ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3443202ece2SStefano Zampini       if (!Cdense) {
3453202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
3463202ece2SStefano Zampini       }
3475ec10c6aSStefano Zampini     } else {
3485ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3495ec10c6aSStefano Zampini       if (!Bdense) {
3505ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3515ec10c6aSStefano Zampini       }
3525ec10c6aSStefano Zampini     }
3535ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
354f11841e3SStefano Zampini   }
3553202ece2SStefano Zampini 
3563202ece2SStefano Zampini   if (D) {
3573202ece2SStefano Zampini     Mat       Dd;
3583202ece2SStefano Zampini     PetscBool Ddense;
3593202ece2SStefano Zampini 
3603202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
3613202ece2SStefano Zampini     if (!Ddense) {
3623202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
3633202ece2SStefano Zampini     } else {
3643202ece2SStefano Zampini       Dd = D;
3653202ece2SStefano Zampini     }
366f11841e3SStefano Zampini     if (n_I) {
3673202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
368f11841e3SStefano Zampini     } else {
369f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
370f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
371f11841e3SStefano Zampini       } else {
372f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
373f11841e3SStefano Zampini       }
374f11841e3SStefano Zampini     }
3753202ece2SStefano Zampini     if (!Ddense) {
3763202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
3773202ece2SStefano Zampini     }
3783202ece2SStefano Zampini   } else {
3793202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
3803202ece2SStefano Zampini   }
3813202ece2SStefano Zampini   PetscFunctionReturn(0);
3823202ece2SStefano Zampini }
38334a97f8cSStefano Zampini 
38491af6908SStefano 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)
385b1b3d7a2SStefano Zampini {
386be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
387be83ff47SStefano Zampini   Mat                    S_all;
3885a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
3895db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
390b7ab4a40SStefano Zampini   IS                     is_I,is_I_layer;
391dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
392dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
393dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
3945a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
395683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
39608122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
39706a4b1faSStefano Zampini   PetscScalar            *Bwork;
3985a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
3995a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
4005a95e1ceSStefano Zampini   MPI_Comm               comm_n;
4019d54b7f4SStefano Zampini   PetscBool              deluxe = PETSC_TRUE, use_getr = PETSC_FALSE;
402b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
403b1b3d7a2SStefano Zampini 
404b1b3d7a2SStefano Zampini   PetscFunctionBegin;
405a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
406a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
407a64f4aa4SStefano Zampini 
4083301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_TRUE;
4093301b35fSStefano Zampini   sub_schurs->is_posdef    = PETSC_TRUE;
410e62b6521Sstefano_zampini   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
411e62b6521Sstefano_zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
412e62b6521Sstefano_zampini   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
413e62b6521Sstefano_zampini   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
414e62b6521Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
415e62b6521Sstefano_zampini 
416e62b6521Sstefano_zampini   /* convert matrix if needed */
417e62b6521Sstefano_zampini   if (Ain) {
418e62b6521Sstefano_zampini     PetscBool isseqaij;
419e62b6521Sstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
420a64f4aa4SStefano Zampini     if (isseqaij) {
421a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
422a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4233301b35fSStefano Zampini     } else {
424a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
425a64f4aa4SStefano Zampini     }
426a64f4aa4SStefano Zampini   }
4273301b35fSStefano Zampini 
428a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
429a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
430df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
431df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
432a64f4aa4SStefano Zampini   }
433a64f4aa4SStefano Zampini 
4345a95e1ceSStefano Zampini   /* preliminary checks */
435af25d912SStefano 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");
4365a95e1ceSStefano Zampini 
4375a95e1ceSStefano Zampini   /* restrict work on active processes */
4385a95e1ceSStefano Zampini   color = 0;
4395a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4405a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4415a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4425a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4435a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4445a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
4455a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
4465a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4475a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
4485a95e1ceSStefano Zampini     PetscFunctionReturn(0);
4495a95e1ceSStefano Zampini   }
45081ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
4515a95e1ceSStefano Zampini 
452b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
453df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
454a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
4553301b35fSStefano Zampini     PetscBool isseqsbaij;
456a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
4573301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
4583301b35fSStefano Zampini     if (isseqsbaij) {
459a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
460a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
461a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
462a64f4aa4SStefano Zampini     } else {
463a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
464a64f4aa4SStefano Zampini       A_BB = tA_BB;
465a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
466a64f4aa4SStefano Zampini       A_IB = tA_IB;
467a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
468a64f4aa4SStefano Zampini       A_BI = tA_BI;
469f11841e3SStefano Zampini     }
470a58a30b4SStefano Zampini   } else {
4715a95e1ceSStefano Zampini     A_II = NULL;
4725a95e1ceSStefano Zampini     A_IB = NULL;
4735a95e1ceSStefano Zampini     A_BI = NULL;
4745a95e1ceSStefano Zampini     A_BB = NULL;
475b1b3d7a2SStefano Zampini   }
4765a95e1ceSStefano Zampini   S_all = NULL;
477b1b3d7a2SStefano Zampini 
478b1b3d7a2SStefano Zampini   /* determine interior problems */
4793dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4803dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
481b1b3d7a2SStefano Zampini     PetscBT                touched;
482b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
483b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
484b1b3d7a2SStefano Zampini 
4851633d1f0SStefano Zampini     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
486b1b3d7a2SStefano Zampini     /* get sizes */
487b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
488b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
489b1b3d7a2SStefano Zampini 
490b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
491b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
492b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
493b1b3d7a2SStefano Zampini 
494b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
495b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
496b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
497b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
498b1b3d7a2SStefano Zampini     }
499b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
500b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
501b1b3d7a2SStefano Zampini 
502b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
503b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
504b1b3d7a2SStefano Zampini     n_prev_added = n_B;
505b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
506b1b3d7a2SStefano Zampini       PetscInt n_added;
507b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5086c4ed002SBarry 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);
509b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
510b1b3d7a2SStefano Zampini       n_prev_added = n_added;
511b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
512b1b3d7a2SStefano Zampini       if (!n_added) break;
513b1b3d7a2SStefano Zampini     }
514b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
515b1b3d7a2SStefano Zampini 
516883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
517a9b99552SStefano 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);
518b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
519a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
520883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
521df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
522b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
523b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
524a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
525b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
526b1b3d7a2SStefano Zampini 
527b1b3d7a2SStefano Zampini       /* II block */
5287dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
529b1b3d7a2SStefano Zampini     }
530b1b3d7a2SStefano Zampini   } else {
531b1b3d7a2SStefano Zampini     PetscInt n_I;
532b1b3d7a2SStefano Zampini 
533b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
534b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
535a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
536b1b3d7a2SStefano Zampini 
537b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
538df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
539b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
540b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
541b1b3d7a2SStefano Zampini 
542b1b3d7a2SStefano Zampini       /* II block is the same */
543b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
544b1b3d7a2SStefano Zampini       AE_II = A_II;
545b1b3d7a2SStefano Zampini     }
546b1b3d7a2SStefano Zampini   }
5475a95e1ceSStefano Zampini 
548883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
5495a95e1ceSStefano Zampini   max_subset_size = 0;
550883469d8SStefano Zampini   local_size = 0;
5515a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
5525a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5535a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
554883469d8SStefano Zampini     local_size += subset_size;
555883469d8SStefano Zampini   }
556883469d8SStefano Zampini 
557883469d8SStefano Zampini   /* Work arrays for local indices */
558883469d8SStefano Zampini   extra = 0;
559683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
560df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
561a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
562883469d8SStefano Zampini   }
563683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
564883469d8SStefano Zampini   if (extra) {
565883469d8SStefano Zampini     const PetscInt *idxs;
566a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
567883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
568a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
569883469d8SStefano Zampini   }
570883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
571dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
572dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
573883469d8SStefano Zampini 
574883469d8SStefano Zampini   /* Get local indices in local numbering */
575883469d8SStefano Zampini   local_size = 0;
5765a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
577883469d8SStefano Zampini     PetscInt j;
578883469d8SStefano Zampini     const    PetscInt *idxs;
579883469d8SStefano Zampini 
5805a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5815a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
582eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
583eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
584eb595f79SStefano Zampini     auxnum2[i] = subset_size;
585883469d8SStefano Zampini     /* subset indices in local numbering */
586883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
5875a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
588883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
589883469d8SStefano Zampini     local_size += subset_size;
590883469d8SStefano Zampini   }
591883469d8SStefano Zampini 
5925a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
59306a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
59459ac4de7SStefano Zampini     PetscScalar  lwork,dummyscalar = 0.;
59559ac4de7SStefano Zampini     PetscBLASInt dummyint = 0;
596d2627357SStefano Zampini 
59706a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
598d2627357SStefano Zampini     B_lwork = -1;
599d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
600d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
60159ac4de7SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
602d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
603d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
604d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
605d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
606056290a2SStefano Zampini   } else {
607056290a2SStefano Zampini     Bwork = NULL;
608056290a2SStefano Zampini     pivots = NULL;
609d2627357SStefano Zampini   }
610d2627357SStefano Zampini 
611d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
612dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
613dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
614dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
615dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
6166583bcc1SStefano Zampini   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
617dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
618dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
619dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6206c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
621dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
622e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
623d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
624d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
625d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
626d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6272972d61bSStefano Zampini 
6285a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6295a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6305a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6315a95e1ceSStefano Zampini 
6325a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6335a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
6346c4ed002SBarry 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);
6355a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
636b1b3d7a2SStefano Zampini   }
637b1b3d7a2SStefano Zampini 
63872b8c272SStefano Zampini   if (change) {
63972b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
64072b8c272SStefano Zampini     IS                     change_primal_B;
64172b8c272SStefano Zampini     IS                     change_primal_all;
64272b8c272SStefano Zampini 
643b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
644b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
645b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
64672b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
64772b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
64872b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
649b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
65072b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
65172b8c272SStefano Zampini     }
65272b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
65372b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
65472b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
65572b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
65672b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
657b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
65872b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
65972b8c272SStefano Zampini       Mat change_sub;
66072b8c272SStefano Zampini 
661b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
66272b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
66372b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
66472b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
6657dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
66672b8c272SStefano Zampini       } else {
66772b8c272SStefano Zampini         Mat change_subt;
6687dae84e0SHong Zhang         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
66972b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
67072b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
67172b8c272SStefano Zampini       }
67272b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
67372b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
674e62b6521Sstefano_zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
675e62b6521Sstefano_zampini       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
67672b8c272SStefano Zampini     }
67772b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
67872b8c272SStefano Zampini   }
67972b8c272SStefano Zampini 
6805a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
6815a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
6825a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
6835a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
6845a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
6855a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
686aa83b6aeSStefano Zampini   }
687b1b3d7a2SStefano Zampini 
6885a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
689be83ff47SStefano Zampini   F = NULL;
690df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) { /* this code branch is used when MatFactor with Schur complemnent support is not present; it is not very efficient, unless the economic version of the scaling is required */
6915a95e1ceSStefano Zampini     Mat         S_Ej_expl;
6925a95e1ceSStefano Zampini     PetscScalar *work;
6935a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
6945a95e1ceSStefano Zampini     PetscBool   Sdense;
6955a95e1ceSStefano Zampini 
6965a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
6975a95e1ceSStefano Zampini     local_size = 0;
698b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
6995a95e1ceSStefano Zampini       IS  is_subset_B;
7005a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
7015a95e1ceSStefano Zampini 
7025a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7035a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7045a95e1ceSStefano Zampini       /* EE block */
7057dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7065a95e1ceSStefano Zampini       /* IE block */
7077dae84e0SHong Zhang       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7085a95e1ceSStefano Zampini       /* EI block */
7095a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
7105a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7115a95e1ceSStefano Zampini       } else {
7127dae84e0SHong Zhang         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7135a95e1ceSStefano Zampini       }
714a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7155a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7165a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7175a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7185a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
719b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
720b1b3d7a2SStefano Zampini         KSP ksp;
721b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7225a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
723b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
724b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
725b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
726b1b3d7a2SStefano Zampini         KSPType   ksp_type;
727b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7285a95e1ceSStefano Zampini         PetscBool ispcnone;
729b1b3d7a2SStefano Zampini 
730b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7315a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
732b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
733b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
734b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
735b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7365a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7375a95e1ceSStefano Zampini         if (!ispcnone) {
7385a95e1ceSStefano Zampini           PCType pc_type;
739b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
740b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7415a95e1ceSStefano Zampini         } else {
7425a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7435a95e1ceSStefano Zampini         }
744b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
745b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
746b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
747b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
748b1b3d7a2SStefano Zampini           if (solver) {
749b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
750b1b3d7a2SStefano Zampini           }
751b1b3d7a2SStefano Zampini         }
752b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
753b1b3d7a2SStefano Zampini       }
7545a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7555a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
7565a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
7575a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
7585a95e1ceSStefano Zampini       if (Sdense) {
7595a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
7605a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
761b1b3d7a2SStefano Zampini         }
7625a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
7636c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
7645a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
765a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
7665a95e1ceSStefano Zampini       local_size += subset_size;
7675a95e1ceSStefano Zampini     }
7685a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
769b1b3d7a2SStefano Zampini     /* free */
770b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
771b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
7725a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
773883469d8SStefano Zampini   } else {
7745cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
7759d54b7f4SStefano Zampini     Vec         Dall = NULL;
776ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
7775cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
7786dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
7793fc34f97SStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE;
7803fc34f97SStefano Zampini     PetscBool   schur_has_vertices,factor_workaround;
781d4933d67SStefano Zampini     char        solver[256];
782883469d8SStefano Zampini 
783683d3df6SStefano Zampini     /* get sizes */
78481ea8064SStefano Zampini     n_I = 0;
78581ea8064SStefano Zampini     if (is_I_layer) {
786a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
78781ea8064SStefano Zampini     }
788683d3df6SStefano Zampini     economic = PETSC_FALSE;
789683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
790683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
791683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
7929d54b7f4SStefano Zampini     size_active_schur = local_size;
7939d54b7f4SStefano Zampini 
794f17d2ae1SStefano Zampini     /* import scaling vector (wrong formulation if we have 3D edges) */
7959d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
7969d54b7f4SStefano Zampini       const PetscScalar *array;
7979d54b7f4SStefano Zampini       PetscScalar       *array2;
7989d54b7f4SStefano Zampini       const PetscInt    *idxs;
7999d54b7f4SStefano Zampini       PetscInt          i;
8009d54b7f4SStefano Zampini 
8019d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8029d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8039d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8049d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8059d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8069d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8079d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8089d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8099d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8109d54b7f4SStefano Zampini     }
811d62866d3SStefano Zampini 
812683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
8133fc34f97SStefano Zampini     factor_workaround = PETSC_FALSE;
8143fc34f97SStefano Zampini     schur_has_vertices = PETSC_FALSE;
815683d3df6SStefano Zampini     cum = n_I+size_active_schur;
816683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
817683d3df6SStefano Zampini       const PetscInt* idxs;
818683d3df6SStefano Zampini       PetscInt        n_dir;
819683d3df6SStefano Zampini 
820683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
821683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
822683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
823683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
824683d3df6SStefano Zampini       cum += n_dir;
8253fc34f97SStefano Zampini       factor_workaround = PETSC_TRUE;
826d62866d3SStefano Zampini     }
827683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
828367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
829683d3df6SStefano Zampini       PetscInt n_v;
830683d3df6SStefano Zampini 
831683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
832683d3df6SStefano Zampini       if (n_v) {
833683d3df6SStefano Zampini         const PetscInt* idxs;
834683d3df6SStefano Zampini 
835683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
836683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
837683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
838683d3df6SStefano Zampini         cum += n_v;
839683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
8403fc34f97SStefano Zampini         schur_has_vertices = PETSC_TRUE;
841683d3df6SStefano Zampini       }
842683d3df6SStefano Zampini     }
843683d3df6SStefano Zampini     size_schur = cum - n_I;
844683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
845683d3df6SStefano Zampini     /* get working mat (TODO: factorize without actually permuting)  */
846683d3df6SStefano Zampini     if (cum == n) {
847683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
848683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
849683d3df6SStefano Zampini     } else {
8507dae84e0SHong Zhang       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
851683d3df6SStefano Zampini     }
852e62b6521Sstefano_zampini     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
853e62b6521Sstefano_zampini     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
854ca92afb2SStefano Zampini 
855ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
856367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
857ca92afb2SStefano Zampini     if (benign_n) {
858ca92afb2SStefano Zampini       Vec                    v;
859ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
860ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
8615cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
8625cbda25cSStefano Zampini       PetscInt               sizeA;
863ca92afb2SStefano Zampini 
864ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
865ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
866ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
867ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
868ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
8695cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
8705cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
871282d6408SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
8725cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
8735cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
874ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
875ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
876ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
8775cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
878ca92afb2SStefano Zampini         PetscScalar    *array;
879ca92afb2SStefano Zampini         const PetscInt *idxs;
880ca92afb2SStefano Zampini         PetscInt       j,nz;
881ca92afb2SStefano Zampini 
882ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
883ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
884ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8855cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
886ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8875cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
8885cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
8895cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
890ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
89122db5ddcSStefano Zampini         for (j=0;j<size_schur;j++) {
89222db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
89322db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
89422db5ddcSStefano Zampini #else
89522db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
89622db5ddcSStefano Zampini #endif
89722db5ddcSStefano Zampini         }
898ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
899ca92afb2SStefano Zampini       }
9005cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9015cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
902ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
903ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
904ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
905ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
906ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
907ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
908ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
909ca92afb2SStefano Zampini     }
9103301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
9113301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9123301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
913df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
914d4933d67SStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
915df4d28bfSStefano Zampini #else
916df4d28bfSStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr);
917df4d28bfSStefano Zampini #endif
918e62b6521Sstefano_zampini     ierr = PetscOptionsGetString(NULL,((PetscObject)A)->prefix,"-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
919883469d8SStefano Zampini 
920683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
921d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
922d47842beSStefano Zampini 
923683d3df6SStefano Zampini     if (n_I) { /* TODO add ordering from options */
9245a05ddb0SStefano Zampini       IS is_schur;
9255a05ddb0SStefano Zampini 
9269ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
927d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
928883469d8SStefano Zampini       } else {
929d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
930883469d8SStefano Zampini       }
9318c09ecd8Sstefano_zampini       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
9328c09ecd8Sstefano_zampini 
933883469d8SStefano Zampini       /* subsets ordered last */
9345a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9355a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9365a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
937883469d8SStefano Zampini 
938883469d8SStefano Zampini       /* factorization step */
9399ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
940883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
941be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
942be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
943be83ff47SStefano Zampini #endif
944883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
945a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
946883469d8SStefano Zampini       } else {
947883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
948be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
949be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
950be83ff47SStefano Zampini #endif
951883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
952883469d8SStefano Zampini       }
953ec968b6dSstefano_zampini       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
954883469d8SStefano Zampini 
955883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
9567c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
957b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
958b3cb21ddSStefano Zampini       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
959b3cb21ddSStefano Zampini 
960d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
961683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
962683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
963df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
964ca92afb2SStefano Zampini 
96572b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
966ca92afb2SStefano Zampini       if (benign_n) {
9675cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
968ca92afb2SStefano Zampini         Vec         v;
9695cbda25cSStefano Zampini         PetscInt    sizeA;
970ca92afb2SStefano Zampini 
971ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
9725cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
9735cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
974ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
975ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
976ca92afb2SStefano Zampini #endif
977ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
978ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
979ca92afb2SStefano Zampini #endif
9805cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9815cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
982ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
9835cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
98447484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
985ca92afb2SStefano Zampini           const PetscInt *idxs;
986ca92afb2SStefano Zampini           PetscInt       j,nz;
98747484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
988ca92afb2SStefano Zampini 
9895cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
9905cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
9915cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
992ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
993ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
994282d6408SStefano Zampini           /* p0 dof (eliminated) is excluded from the sum */
995282d6408SStefano Zampini           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
996ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9975cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
998ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
99947484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
100047484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
100147484b83SStefano Zampini           B_k = 1;
100247484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
100347484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
10045cbda25cSStefano Zampini           sum  = 1.;
100547484b83SStefano 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));
1006ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10075cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10085cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1009282d6408SStefano Zampini           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
10105cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10115cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1012ca92afb2SStefano Zampini         }
10135cbda25cSStefano Zampini   /* restore defaults */
10145cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
10155cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
10165cbda25cSStefano Zampini #endif
10175cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
10185cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
10195cbda25cSStefano Zampini #endif
10205cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10215cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1022ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1023ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1024ca92afb2SStefano Zampini       }
1025a3df083aSStefano Zampini       if (!reuse_solvers) {
1026a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1027a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1028a3df083aSStefano Zampini         }
1029a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
10305cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
10315cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1032a3df083aSStefano Zampini       }
1033df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1034be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1035683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1036166598c1SStefano Zampini       factor_workaround = PETSC_FALSE;
1037df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1038be83ff47SStefano Zampini     }
1039be83ff47SStefano Zampini 
1040be83ff47SStefano Zampini     if (reuse_solvers) {
1041a00504b5SStefano Zampini       Mat                A_II,Afake;
104253892102SStefano Zampini       Vec                vec1_B;
1043df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
10443462e049SStefano Zampini       PetscInt           n_R;
1045d5574798SStefano Zampini 
1046df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1047df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1048e28d306cSStefano Zampini       } else {
1049df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1050d62866d3SStefano Zampini       }
1051df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1052be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1053d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1054d5574798SStefano Zampini       msolv_ctx->F = F;
1055683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1056683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1057683d3df6SStefano Zampini       {
1058683d3df6SStefano Zampini         PetscScalar *array;
1059683d3df6SStefano Zampini         PetscInt    n;
1060683d3df6SStefano Zampini 
1061683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1062683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1063683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1064683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1065683d3df6SStefano Zampini       }
10663fc34f97SStefano Zampini       msolv_ctx->has_vertices = schur_has_vertices;
1067d62866d3SStefano Zampini 
1068d62866d3SStefano Zampini       /* interior solver */
1069683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1070d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1071d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1072d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1073df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1074df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1075d62866d3SStefano Zampini 
1076d62866d3SStefano Zampini       /* correction solver */
10773462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1078d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1079d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1080df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1081df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
108253892102SStefano Zampini 
108353892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
108453892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
108553892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
10863fc34f97SStefano Zampini       if (!schur_has_vertices) {
108753892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
108853892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
108953892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
109053892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1091683d3df6SStefano Zampini       } else {
1092683d3df6SStefano Zampini         IS              is_B_all;
1093683d3df6SStefano Zampini         const PetscInt* idxs;
1094683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1095683d3df6SStefano Zampini 
1096683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1097683d3df6SStefano Zampini         dual = size_schur - n_v;
1098683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1099683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1100683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1101683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1102683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1103683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1104683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1105683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1106683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1107683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1108683d3df6SStefano Zampini       }
11093462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
11103462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
11113462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
11123462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1113683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1114ca92afb2SStefano Zampini 
1115ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1116ca92afb2SStefano Zampini       if (benign_n) {
11175cbda25cSStefano Zampini         PetscScalar *array;
11185cbda25cSStefano Zampini 
1119ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1120ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1121ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
11225cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
1123282d6408SStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
11245cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11255cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
11265cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11275cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1128ca92afb2SStefano Zampini       }
1129d5574798SStefano Zampini     }
113008122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
113153892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
11325db18549SStefano Zampini 
1133be83ff47SStefano Zampini     /* Work arrays */
1134be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
11355a94c880SStefano Zampini 
1136*04270d10SStefano Zampini     /* matrices for deluxe scaling and adaptive selection */
1137*04270d10SStefano Zampini     if (!sub_schurs->sum_S_Ej_all) {
1138*04270d10SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1139*04270d10SStefano Zampini       ierr = MatSetSizes(sub_schurs->sum_S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1140*04270d10SStefano Zampini       ierr = MatSetType(sub_schurs->sum_S_Ej_all,MATAIJ);CHKERRQ(ierr);
1141*04270d10SStefano Zampini       ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_all,0,nnz);CHKERRQ(ierr);
1142*04270d10SStefano Zampini     }
114312d906b1SStefano Zampini     if (compute_Stilda) {
11445a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
1145*04270d10SStefano Zampini         ierr = MatDuplicate(sub_schurs->sum_S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11465a94c880SStefano Zampini       }
11479d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1148*04270d10SStefano Zampini         ierr = MatDuplicate(sub_schurs->sum_S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
11495a94c880SStefano Zampini       }
115012d906b1SStefano Zampini     }
1151d2627357SStefano Zampini 
1152be83ff47SStefano Zampini     /* S_Ej_all */
1153be83ff47SStefano Zampini     cum = cum2 = 0;
1154be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
115565d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
115665d8bf0aSStefano Zampini       PetscInt j;
115765d8bf0aSStefano Zampini 
11585a95e1ceSStefano Zampini       /* get S_E */
1159b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1160683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1161be83ff47SStefano Zampini         PetscInt k;
1162be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1163be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1164be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1165683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1166be83ff47SStefano Zampini           }
1167be83ff47SStefano Zampini         }
116806a4e24aSStefano Zampini       } else { /* just copy to workspace */
1169be83ff47SStefano Zampini         PetscInt k;
1170be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1171be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1172be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1173be83ff47SStefano Zampini           }
1174be83ff47SStefano Zampini         }
11759087bf02SStefano Zampini       }
11765a95e1ceSStefano Zampini       /* insert S_E values */
1177be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1178b7ab4a40SStefano Zampini       if (sub_schurs->change) {
11798760537fSStefano Zampini         Mat            change_sub,SEj,T;
118072b8c272SStefano Zampini 
118172b8c272SStefano Zampini         /* change basis */
118272b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
118372b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
11848760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
11858760537fSStefano Zampini           Mat T2;
11868760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
11878760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
11888760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
118972b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
11908760537fSStefano Zampini         } else {
11918760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
119272b8c272SStefano Zampini         }
119372b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
119472b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
1195afa0f562SStefano Zampini         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
119672b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
119772b8c272SStefano Zampini       }
11989d54b7f4SStefano Zampini       if (deluxe) {
11995a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1200683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1201862806e4SStefano Zampini         if (compute_Stilda) {
120208122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
120308122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
120406a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
12056c3e6151SStefano Zampini             PetscInt k;
12065a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12072972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
12085a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
12092972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
12106c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
12116c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
12126c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
12136c3e6151SStefano Zampini               }
12146c3e6151SStefano Zampini             }
1215d6462365SStefano Zampini           } else {
1216d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1217d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1218d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1219d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
12202972d61bSStefano Zampini           }
122108122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
12225a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12239087bf02SStefano Zampini         }
12249d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
12259d54b7f4SStefano Zampini         Mat         SEj;
12269d54b7f4SStefano Zampini         Vec         D;
12279d54b7f4SStefano Zampini         PetscScalar *array;
12289d54b7f4SStefano Zampini 
12299d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12309d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
12319d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
12329d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1233f17d2ae1SStefano Zampini         ierr = VecShift(D,-1.);CHKERRQ(ierr);
12349d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
12359d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
12369d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
12379d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12389d54b7f4SStefano Zampini       }
1239be83ff47SStefano Zampini       cum += subset_size;
1240be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1241be83ff47SStefano Zampini     }
1242be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
12434a6c6b0dSStefano Zampini 
1244df4d28bfSStefano Zampini     if (solver_S) {
12457c2f51b8SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
12464a6c6b0dSStefano Zampini     }
1247683d3df6SStefano Zampini 
1248683d3df6SStefano Zampini     schur_factor = NULL;
124945951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1250683d3df6SStefano Zampini 
12519d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
12524a6c6b0dSStefano Zampini         PetscInt j;
12534a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
12545a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12554a6c6b0dSStefano Zampini       } else {
1256683d3df6SStefano Zampini         Mat S_all_inv=NULL;
12573fc34f97SStefano Zampini         if (solver_S) {
1258683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1259683d3df6SStefano 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 */
12603fc34f97SStefano Zampini           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1261683d3df6SStefano Zampini             PetscScalar *data;
1262683d3df6SStefano Zampini             PetscInt     nd = 0;
12636dba178dSStefano Zampini 
1264683d3df6SStefano Zampini             if (!sub_schurs->is_posdef) {
1265683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1266683d3df6SStefano Zampini             }
12677c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1268683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1269683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1270683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1271683d3df6SStefano Zampini             }
12723fc34f97SStefano Zampini 
12733fc34f97SStefano Zampini             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
12743fc34f97SStefano Zampini             if (schur_has_vertices) {
12753fc34f97SStefano Zampini               Mat          M;
12763fc34f97SStefano Zampini               PetscScalar *tdata;
12773fc34f97SStefano Zampini               PetscInt     nv = 0, news;
12783fc34f97SStefano Zampini 
12793fc34f97SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
12803fc34f97SStefano Zampini               news = size_active_schur + nv;
12813fc34f97SStefano Zampini               ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr);
1282683d3df6SStefano Zampini               for (i=0;i<size_active_schur;i++) {
12833fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
12843fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));CHKERRQ(ierr);
12853fc34f97SStefano Zampini               }
12863fc34f97SStefano Zampini               for (i=0;i<nv;i++) {
12873fc34f97SStefano Zampini                 PetscInt k = i+size_active_schur;
12883fc34f97SStefano Zampini                 ierr = PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));CHKERRQ(ierr);
12893fc34f97SStefano Zampini               }
12903fc34f97SStefano Zampini 
12913fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr);
12923fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
12933fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
12943fc34f97SStefano Zampini               /* save the factors */
12953fc34f97SStefano Zampini               cum = 0;
12963fc34f97SStefano Zampini               ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr);
12973fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
12983fc34f97SStefano Zampini                 ierr = PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1299683d3df6SStefano Zampini                 cum += size_active_schur - i;
1300683d3df6SStefano Zampini               }
13013fc34f97SStefano Zampini               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
13023fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
13033fc34f97SStefano Zampini               /* move back just the active dofs to the Schur complement */
13043fc34f97SStefano Zampini               for (i=0;i<size_active_schur;i++) {
13053fc34f97SStefano Zampini                 ierr = PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));CHKERRQ(ierr);
13063fc34f97SStefano Zampini               }
13073fc34f97SStefano Zampini               ierr = PetscFree(tdata);CHKERRQ(ierr);
13083fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13093fc34f97SStefano Zampini             } else { /* we can factorize and invert just the activedofs */
13103fc34f97SStefano Zampini               Mat M;
13113fc34f97SStefano Zampini 
13123fc34f97SStefano Zampini               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr);
13133fc34f97SStefano Zampini               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
13143fc34f97SStefano Zampini               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
13153fc34f97SStefano Zampini               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
13163fc34f97SStefano Zampini               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
13173fc34f97SStefano Zampini               ierr = MatDestroy(&M);CHKERRQ(ierr);
13183fc34f97SStefano Zampini               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
13193fc34f97SStefano Zampini             }
1320683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
13213fc34f97SStefano Zampini           } else { /* use MatFactor calls to invert S */
13226dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
13237c2f51b8SStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1324683d3df6SStefano Zampini           }
1325df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1326683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1327683d3df6SStefano Zampini           S_all_inv = S_all;
1328683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1329be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1330be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
133106a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1332be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1333be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1334be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1335be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1336d6462365SStefano Zampini           } else {
1337d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1338d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1339d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1340d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1341be83ff47SStefano Zampini           }
1342be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1343683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1344be83ff47SStefano Zampini         }
1345be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1346be83ff47SStefano Zampini         cum = cum2 = 0;
1347683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1348be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1349be83ff47SStefano Zampini           PetscInt j;
1350862806e4SStefano Zampini 
1351862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1352be83ff47SStefano Zampini           /* get (St^-1)_E */
135372b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
135406a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1355a0b0af32SStefano Zampini           if (S_lower_triangular) {
1356be83ff47SStefano Zampini             PetscInt k;
1357b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1358be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1359be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1360be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
13616c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1362be83ff47SStefano Zampini                 }
1363be83ff47SStefano Zampini               }
136472b8c272SStefano Zampini             } else {
136572b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
136672b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
136772b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
136872b8c272SStefano Zampini                 }
136972b8c272SStefano Zampini               }
137072b8c272SStefano Zampini             }
137172b8c272SStefano Zampini           } else {
1372be83ff47SStefano Zampini             PetscInt k;
1373be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1374be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1375be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1376be83ff47SStefano Zampini               }
1377be83ff47SStefano Zampini             }
1378be83ff47SStefano Zampini           }
1379b7ab4a40SStefano Zampini           if (sub_schurs->change) {
13808760537fSStefano Zampini             Mat change_sub,SEj,T;
138172b8c272SStefano Zampini 
138272b8c272SStefano Zampini             /* change basis */
138372b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
138472b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
13858760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
13868760537fSStefano Zampini               Mat T2;
13878760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
13888760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
138972b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
13908760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
13918760537fSStefano Zampini             } else {
13928760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
139372b8c272SStefano Zampini             }
13948760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
13958760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
139672b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1397afa0f562SStefano Zampini             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
139872b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
139972b8c272SStefano Zampini           }
1400be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
14015a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1402be83ff47SStefano Zampini           cum += subset_size;
1403be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1404883469d8SStefano Zampini         }
1405683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1406df4d28bfSStefano Zampini         if (solver_S) {
14073fc34f97SStefano Zampini           if (schur_has_vertices) {
14083fc34f97SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
14093fc34f97SStefano Zampini           } else {
14107c2f51b8SStefano Zampini             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
14115db18549SStefano Zampini           }
14123fc34f97SStefano Zampini         }
1413683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1414683d3df6SStefano Zampini       }
1415683d3df6SStefano Zampini 
14163fc34f97SStefano Zampini       /* move back factors if needed */
14173fc34f97SStefano Zampini       if (schur_has_vertices) {
1418683d3df6SStefano Zampini         Mat      S_tmp;
14193fc34f97SStefano Zampini         PetscInt nd = 0;
1420683d3df6SStefano Zampini 
1421df4d28bfSStefano Zampini         if (!solver_S) {
1422683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1423683d3df6SStefano Zampini         }
14247c2f51b8SStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1425683d3df6SStefano Zampini         if (sub_schurs->is_posdef) {
1426683d3df6SStefano Zampini           PetscScalar *data;
1427683d3df6SStefano Zampini 
1428683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
14293fc34f97SStefano Zampini           ierr = PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1430683d3df6SStefano Zampini 
1431683d3df6SStefano Zampini           if (S_lower_triangular) {
1432683d3df6SStefano Zampini             cum = 0;
1433683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1434683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1435683d3df6SStefano Zampini               cum += size_active_schur-i;
1436683d3df6SStefano Zampini 	    }
1437683d3df6SStefano Zampini           } else {
1438683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1439683d3df6SStefano Zampini           }
1440683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1441683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1442683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1443683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1444683d3df6SStefano Zampini 	    }
1445683d3df6SStefano Zampini           }
14466dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1447683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1448683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
14496c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1450683d3df6SStefano Zampini           }
1451683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1452683d3df6SStefano Zampini         } else {
1453683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1454683d3df6SStefano Zampini         }
14553fc34f97SStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
14569087bf02SStefano Zampini       }
1457367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1458367aa537SStefano Zampini       PetscScalar *data;
1459367aa537SStefano Zampini       PetscInt    nd = 0;
1460367aa537SStefano Zampini 
1461367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1462367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1463367aa537SStefano Zampini       }
14647c2f51b8SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1465367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1466367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1467367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1468367aa537SStefano Zampini       }
1469367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1470367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
14716c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1472367aa537SStefano Zampini       }
1473367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
14743fc34f97SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
14754a6c6b0dSStefano Zampini     }
14764a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1477683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
14789d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
14794a6c6b0dSStefano Zampini   }
14805a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1481be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1482a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1483a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1484a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1485a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1486a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
14875db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14885db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14895a95e1ceSStefano Zampini   if (compute_Stilda) {
14905a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14915a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14929d54b7f4SStefano Zampini     if (deluxe) {
14935a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14945a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
149508122e43SStefano Zampini     }
14969d54b7f4SStefano Zampini   }
1497a1337663SStefano Zampini 
14985db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
14995db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
15003927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
15013927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15025a95e1ceSStefano Zampini 
15035db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
15047dae84e0SHong Zhang   ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
1505*04270d10SStefano Zampini   ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
150608122e43SStefano Zampini 
1507f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
15085a95e1ceSStefano Zampini   if (compute_Stilda) {
15095a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1510a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1511*04270d10SStefano Zampini     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1512*04270d10SStefano Zampini     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
15139d54b7f4SStefano Zampini     if (deluxe) {
15145a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
151508122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1516*04270d10SStefano Zampini       ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1517*04270d10SStefano Zampini       ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
15189d54b7f4SStefano Zampini     } else {
15199d54b7f4SStefano Zampini       PetscScalar *array;
15209d54b7f4SStefano Zampini       PetscInt    cum;
15219d54b7f4SStefano Zampini 
15229d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15239d54b7f4SStefano Zampini       cum = 0;
15249d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
15259d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
15269d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
15279d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
15289d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
15299d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
15309d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
15319d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
15329d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
15339d54b7f4SStefano Zampini         cum += subset_size*subset_size;
15349d54b7f4SStefano Zampini       }
15359d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15369d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
15379d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
15389d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
15399d54b7f4SStefano Zampini     }
154008122e43SStefano Zampini   }
1541*04270d10SStefano Zampini   ierr = MatDestroySubMatrices(1,&submats);CHKERRQ(ierr);
15423202ece2SStefano Zampini 
15435a95e1ceSStefano Zampini   /* free workspace */
15445a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
154506a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1546a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
15473202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1548dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
15495a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1550b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1551b1b3d7a2SStefano Zampini }
1552b1b3d7a2SStefano Zampini 
15538b6046baSStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1554b1b3d7a2SStefano Zampini {
15559bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
15565a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1557b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1558b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1559b1b3d7a2SStefano Zampini 
1560b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1561b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
15626c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1563b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
15646c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1565b1b3d7a2SStefano Zampini 
1566b1b3d7a2SStefano Zampini   /* reset any previous data */
1567b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1568b1b3d7a2SStefano Zampini 
15695a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
15709bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1571b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
157208122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1573b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1574b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
15758b6046baSStefano Zampini     if (copycc) {
15768b6046baSStefano Zampini       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
15778b6046baSStefano Zampini     } else {
1578c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1579b1b3d7a2SStefano Zampini       all_cc[i] = faces[i];
1580b1b3d7a2SStefano Zampini     }
15818b6046baSStefano Zampini   }
1582b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
15838b6046baSStefano Zampini     if (copycc) {
15848b6046baSStefano Zampini       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
15858b6046baSStefano Zampini     } else {
1586c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1587b1b3d7a2SStefano Zampini       all_cc[n_faces+i] = edges[i];
15888b6046baSStefano Zampini     }
158908122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1590b1b3d7a2SStefano Zampini   }
1591c8272957SStefano Zampini   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1592c8272957SStefano Zampini   sub_schurs->is_vertices = vertices;
1593c8272957SStefano Zampini   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1594d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1595d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1596b1b3d7a2SStefano Zampini 
1597df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
1598df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_FALSE;
1599883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1600df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1601df4d28bfSStefano Zampini #endif
1602df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1603df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1604883469d8SStefano Zampini #endif
1605b1b3d7a2SStefano Zampini 
1606b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1607b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1608b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1609b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
16105db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
16115db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
16125db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
16135db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
16145a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1615b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1616df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1617b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1618b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1619b96c3477SStefano Zampini     }
16209bb4a8caSStefano Zampini   }
1621b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1622b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
162308122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1624b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1625b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1626b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1627b1b3d7a2SStefano Zampini }
1628b1b3d7a2SStefano Zampini 
162934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
163034a97f8cSStefano Zampini {
163134a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
163234a97f8cSStefano Zampini   PetscErrorCode  ierr;
163334a97f8cSStefano Zampini 
163434a97f8cSStefano Zampini   PetscFunctionBegin;
163534a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
16365ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
163734a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
163834a97f8cSStefano Zampini   PetscFunctionReturn(0);
163934a97f8cSStefano Zampini }
164034a97f8cSStefano Zampini 
164134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
164234a97f8cSStefano Zampini {
164334a97f8cSStefano Zampini   PetscInt       i;
164434a97f8cSStefano Zampini   PetscErrorCode ierr;
164534a97f8cSStefano Zampini 
164634a97f8cSStefano Zampini   PetscFunctionBegin;
1647aea80f77Sstefano_zampini   if (!sub_schurs) PetscFunctionReturn(0);
16481e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1649b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1650b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1651b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
16525db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
16535db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
165441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
165541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
165608122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1657a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
16585db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1659d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1660d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
166108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
166234a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1663b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
166434a97f8cSStefano Zampini   }
16655ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1666b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
16673dc780c3SStefano Zampini   }
1668df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1669df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1670d62866d3SStefano Zampini   }
1671df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
167272b8c272SStefano Zampini   if (sub_schurs->change) {
167372b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
167472b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1675b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
167672b8c272SStefano Zampini     }
167772b8c272SStefano Zampini   }
167872b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1679b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
168034a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
168134a97f8cSStefano Zampini   PetscFunctionReturn(0);
168234a97f8cSStefano Zampini }
168334a97f8cSStefano Zampini 
1684aea80f77Sstefano_zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1685aea80f77Sstefano_zampini {
1686aea80f77Sstefano_zampini   PetscErrorCode ierr;
1687aea80f77Sstefano_zampini 
1688aea80f77Sstefano_zampini   PetscFunctionBegin;
1689aea80f77Sstefano_zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1690aea80f77Sstefano_zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1691aea80f77Sstefano_zampini   PetscFunctionReturn(0);
1692aea80f77Sstefano_zampini }
1693aea80f77Sstefano_zampini 
16942a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
169534a97f8cSStefano Zampini {
169634a97f8cSStefano Zampini   PetscInt       i,j,n;
169734a97f8cSStefano Zampini   PetscErrorCode ierr;
169834a97f8cSStefano Zampini 
169934a97f8cSStefano Zampini   PetscFunctionBegin;
170034a97f8cSStefano Zampini   n = 0;
170134a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
170234a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
170334a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
170434a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
170534a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
170634a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
170734a97f8cSStefano Zampini         queue_tip[n] = dof;
170834a97f8cSStefano Zampini         n++;
170934a97f8cSStefano Zampini       }
171034a97f8cSStefano Zampini     }
171134a97f8cSStefano Zampini   }
171234a97f8cSStefano Zampini   *n_added = n;
171334a97f8cSStefano Zampini   PetscFunctionReturn(0);
171434a97f8cSStefano Zampini }
1715