xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 47484b83076eeb5e6243fac95a9ce27ab7ed5f37)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
308122e43SStefano Zampini #include <petscblaslapack.h>
434a97f8cSStefano Zampini 
53202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
65ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
7df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
8df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
9d62866d3SStefano Zampini 
10ca92afb2SStefano Zampini /* if v2 is not present, correction is done in-place */
11ca92afb2SStefano Zampini #undef __FUNCT__
125cbda25cSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolversBenignAdapt"
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 */
255cbda25cSStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,NULL,&size_schur);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 */
335cbda25cSStefano Zampini     ierr = MatMult(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]];
58ca92afb2SStefano Zampini       sum = -sum/nz;
59ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
60ca92afb2SStefano Zampini       ctx->benign_save_vals[n] = array2[cols[nz-1]];
61ca92afb2SStefano Zampini       array2[cols[nz-1]] = sum;
62ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
63ca92afb2SStefano Zampini     }
64ca92afb2SStefano Zampini   } else {
65ca92afb2SStefano Zampini     PetscInt n;
66ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
67ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
68ca92afb2SStefano Zampini       const PetscInt *cols;
69ca92afb2SStefano Zampini       PetscInt       nz,i;
70ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
71ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
72ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
73ca92afb2SStefano Zampini       sum = -sum/nz;
74ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
75ca92afb2SStefano Zampini       array2[cols[nz-1]] = ctx->benign_save_vals[n];
76ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
77ca92afb2SStefano Zampini     }
78ca92afb2SStefano Zampini   }
79ca92afb2SStefano Zampini   if (v2) {
80ca92afb2SStefano Zampini     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
8172b8c272SStefano Zampini     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
82ca92afb2SStefano Zampini   } else {
83ca92afb2SStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
84ca92afb2SStefano Zampini   }
855cbda25cSStefano Zampini   if (!sol && full) {
865cbda25cSStefano Zampini     Vec      usedv;
875cbda25cSStefano Zampini     PetscInt n_I,size_schur;
885cbda25cSStefano Zampini 
895cbda25cSStefano Zampini     /* get sizes */
905cbda25cSStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,NULL,&size_schur);CHKERRQ(ierr);
915cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
925cbda25cSStefano Zampini     n_I = n_I - size_schur;
935cbda25cSStefano Zampini     /* compute schur rhs correction */
945cbda25cSStefano Zampini     if (v2) {
955cbda25cSStefano Zampini       usedv = v2;
965cbda25cSStefano Zampini     } else {
975cbda25cSStefano Zampini       usedv = v;
985cbda25cSStefano Zampini     }
995cbda25cSStefano Zampini     /* apply schur rhs correction */
1005cbda25cSStefano Zampini     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
1015cbda25cSStefano Zampini     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1025cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
1035cbda25cSStefano Zampini     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1045cbda25cSStefano Zampini     ierr = MatMultTransposeAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1055cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1065cbda25cSStefano Zampini   }
107ca92afb2SStefano Zampini   PetscFunctionReturn(0);
108ca92afb2SStefano Zampini }
109ca92afb2SStefano Zampini 
110d62866d3SStefano Zampini #undef __FUNCT__
111df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_Solve_Private"
112df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
113d62866d3SStefano Zampini {
114df4d28bfSStefano Zampini   PCBDDCReuseSolvers ctx;
115683d3df6SStefano Zampini   PetscBool          copy = PETSC_FALSE;
116d62866d3SStefano Zampini   PetscErrorCode     ierr;
117d62866d3SStefano Zampini 
118d62866d3SStefano Zampini   PetscFunctionBegin;
119d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
120683d3df6SStefano Zampini   if (full) {
121d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
122d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
123d62866d3SStefano Zampini #endif
1245cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1255cbda25cSStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
1265cbda25cSStefano Zampini #endif
127683d3df6SStefano Zampini     copy = ctx->has_vertices;
128d4933d67SStefano Zampini   } else { /* interior solver */
1296dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
130683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
1316dba178dSStefano Zampini #endif
132d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
133d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
134d4933d67SStefano Zampini #endif
135683d3df6SStefano Zampini     copy = PETSC_TRUE;
136683d3df6SStefano Zampini   }
137683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
138683d3df6SStefano Zampini   if (copy) {
139ca92afb2SStefano Zampini     PetscInt    n;
140df4d28bfSStefano Zampini     PetscScalar *array,*array_solver;
141ca92afb2SStefano Zampini 
142ca92afb2SStefano Zampini     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
143683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
144df4d28bfSStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
145df4d28bfSStefano Zampini     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
146df4d28bfSStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
147683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
148683d3df6SStefano Zampini 
1495cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
150683d3df6SStefano Zampini     if (transpose) {
151683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
152683d3df6SStefano Zampini     } else {
153683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
154683d3df6SStefano Zampini     }
1555cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
156683d3df6SStefano Zampini 
157683d3df6SStefano Zampini     /* get back data to caller worskpace */
158df4d28bfSStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
159683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
160df4d28bfSStefano Zampini     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
161683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
162df4d28bfSStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
163683d3df6SStefano Zampini   } else {
164ca92afb2SStefano Zampini     if (ctx->benign_n) {
1655cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
166ca92afb2SStefano Zampini       if (transpose) {
167ca92afb2SStefano Zampini         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
168ca92afb2SStefano Zampini       } else {
169ca92afb2SStefano Zampini         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
170ca92afb2SStefano Zampini       }
1715cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
172ca92afb2SStefano Zampini     } else {
173e28d306cSStefano Zampini       if (transpose) {
174e28d306cSStefano Zampini         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
175e28d306cSStefano Zampini       } else {
1766816873aSStefano Zampini         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
177e28d306cSStefano Zampini       }
178683d3df6SStefano Zampini     }
179ca92afb2SStefano Zampini   }
1805cbda25cSStefano Zampini   /* restore defaults */
1815cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1825cbda25cSStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
1835cbda25cSStefano Zampini #endif
184d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
185d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
186d4933d67SStefano Zampini #endif
187d62866d3SStefano Zampini   PetscFunctionReturn(0);
188d62866d3SStefano Zampini }
189d62866d3SStefano Zampini 
190d62866d3SStefano Zampini #undef __FUNCT__
191df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_Correction"
192df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
193e28d306cSStefano Zampini {
194e28d306cSStefano Zampini   PetscErrorCode   ierr;
195e28d306cSStefano Zampini 
196e28d306cSStefano Zampini   PetscFunctionBegin;
197df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
198e28d306cSStefano Zampini   PetscFunctionReturn(0);
199e28d306cSStefano Zampini }
200e28d306cSStefano Zampini 
201e28d306cSStefano Zampini #undef __FUNCT__
202df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_CorrectionTranspose"
203df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
204e28d306cSStefano Zampini {
205e28d306cSStefano Zampini   PetscErrorCode   ierr;
206e28d306cSStefano Zampini 
207e28d306cSStefano Zampini   PetscFunctionBegin;
208df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
209683d3df6SStefano Zampini   PetscFunctionReturn(0);
210683d3df6SStefano Zampini }
211683d3df6SStefano Zampini 
212683d3df6SStefano Zampini #undef __FUNCT__
213df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_Interior"
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 
223683d3df6SStefano Zampini #undef __FUNCT__
224df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_InteriorTranspose"
225df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
226683d3df6SStefano Zampini {
227683d3df6SStefano Zampini   PetscErrorCode   ierr;
228683d3df6SStefano Zampini 
229683d3df6SStefano Zampini   PetscFunctionBegin;
230df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
231e28d306cSStefano Zampini   PetscFunctionReturn(0);
232e28d306cSStefano Zampini }
233e28d306cSStefano Zampini 
234e28d306cSStefano Zampini #undef __FUNCT__
235df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolversReset"
236df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
237d62866d3SStefano Zampini {
238ca92afb2SStefano Zampini   PetscInt       i;
239d62866d3SStefano Zampini   PetscErrorCode ierr;
240d62866d3SStefano Zampini 
241d62866d3SStefano Zampini   PetscFunctionBegin;
242d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
243d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
244d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
245d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
246d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
24753892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
24853892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
249d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
25053892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
25153892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
252ca92afb2SStefano Zampini   for (i=0;i<reuse->benign_n;i++) {
253ca92afb2SStefano Zampini     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
254ca92afb2SStefano Zampini   }
255ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
256ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
2575cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
2585cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
2595cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
2605cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
261d62866d3SStefano Zampini   PetscFunctionReturn(0);
262d62866d3SStefano Zampini }
263d5574798SStefano Zampini 
264d5574798SStefano Zampini #undef __FUNCT__
2653202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
2665ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
2673202ece2SStefano Zampini {
2683202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
2693202ece2SStefano Zampini   KSP            ksp;
2703202ece2SStefano Zampini   PC             pc;
2713202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
2723202ece2SStefano Zampini   PetscReal      fill = 2.0;
273f11841e3SStefano Zampini   PetscInt       n_I;
2743202ece2SStefano Zampini   PetscMPIInt    size;
2753202ece2SStefano Zampini   PetscErrorCode ierr;
2763202ece2SStefano Zampini 
2773202ece2SStefano Zampini   PetscFunctionBegin;
2783202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
2796c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
280f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
281f11841e3SStefano Zampini     PetscBool Sdense;
282f11841e3SStefano Zampini 
283f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
2846c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
285f11841e3SStefano Zampini   }
2863202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
2873202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
2883202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
2893202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
2903202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
2913202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
2923202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
2933202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
294f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
295f11841e3SStefano Zampini   if (n_I) {
2963202ece2SStefano Zampini     if (!Bdense) {
2973202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
2983202ece2SStefano Zampini     } else {
2993202ece2SStefano Zampini       Bd = B;
3003202ece2SStefano Zampini     }
3013202ece2SStefano Zampini 
3023202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
3033202ece2SStefano Zampini       Mat fact;
3043202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
3053202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
3063202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
3073202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
3083202ece2SStefano Zampini     } else {
30907b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
31007b1e237SStefano Zampini 
31107b1e237SStefano Zampini       if (ex) {
3123202ece2SStefano Zampini         Mat Ainvd;
3133202ece2SStefano Zampini 
3143202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
3153202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
3163202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
31707b1e237SStefano Zampini       } else {
31807b1e237SStefano Zampini         Vec         sol,rhs;
31907b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
32007b1e237SStefano Zampini         PetscInt    i,nrhs,n;
32107b1e237SStefano Zampini 
32207b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
32307b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
32407b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
32507b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
32607b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
32707b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
32807b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
32907b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
33007b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
33107b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
33207b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
33307b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
33407b1e237SStefano Zampini         }
33507b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
33607b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
33707b1e237SStefano Zampini       }
3383202ece2SStefano Zampini     }
3395ec10c6aSStefano Zampini     if (!Bdense & !issym) {
3403202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3413202ece2SStefano Zampini     }
3425ec10c6aSStefano Zampini 
3435ec10c6aSStefano Zampini     if (!issym) {
3443202ece2SStefano Zampini       if (!Cdense) {
3453202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
3463202ece2SStefano Zampini       } else {
3473202ece2SStefano Zampini         Cd = C;
3483202ece2SStefano Zampini       }
3495ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3503202ece2SStefano Zampini       if (!Cdense) {
3513202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
3523202ece2SStefano Zampini       }
3535ec10c6aSStefano Zampini     } else {
3545ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3555ec10c6aSStefano Zampini       if (!Bdense) {
3565ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3575ec10c6aSStefano Zampini       }
3585ec10c6aSStefano Zampini     }
3595ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
360f11841e3SStefano Zampini   }
3613202ece2SStefano Zampini 
3623202ece2SStefano Zampini   if (D) {
3633202ece2SStefano Zampini     Mat       Dd;
3643202ece2SStefano Zampini     PetscBool Ddense;
3653202ece2SStefano Zampini 
3663202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
3673202ece2SStefano Zampini     if (!Ddense) {
3683202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
3693202ece2SStefano Zampini     } else {
3703202ece2SStefano Zampini       Dd = D;
3713202ece2SStefano Zampini     }
372f11841e3SStefano Zampini     if (n_I) {
3733202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
374f11841e3SStefano Zampini     } else {
375f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
376f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
377f11841e3SStefano Zampini       } else {
378f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
379f11841e3SStefano Zampini       }
380f11841e3SStefano Zampini     }
3813202ece2SStefano Zampini     if (!Ddense) {
3823202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
3833202ece2SStefano Zampini     }
3843202ece2SStefano Zampini   } else {
3853202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
3863202ece2SStefano Zampini   }
3873202ece2SStefano Zampini   PetscFunctionReturn(0);
3883202ece2SStefano Zampini }
38934a97f8cSStefano Zampini 
39034a97f8cSStefano Zampini #undef __FUNCT__
3911580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
39291af6908SStefano 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)
393b1b3d7a2SStefano Zampini {
394be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
395be83ff47SStefano Zampini   Mat                    S_all;
3965a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
3975db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
398b7ab4a40SStefano Zampini   IS                     is_I,is_I_layer;
399dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
400dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
401dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
4025a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
403683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
40408122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
40506a4b1faSStefano Zampini   PetscScalar            *Bwork;
4065a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
4075a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
4085a95e1ceSStefano Zampini   MPI_Comm               comm_n;
4099d54b7f4SStefano Zampini   PetscBool              deluxe = PETSC_TRUE, use_getr = PETSC_FALSE;
410b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
411b1b3d7a2SStefano Zampini 
412b1b3d7a2SStefano Zampini   PetscFunctionBegin;
413a64f4aa4SStefano Zampini   /* update info in sub_schurs */
414a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
415a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
4163301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
4173301b35fSStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
418a64f4aa4SStefano Zampini   if (Ain) {
419a64f4aa4SStefano Zampini     PetscBool isseqaij;
4203301b35fSStefano Zampini     /* determine if we are dealing with hermitian positive definite problems */
4213301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
4223301b35fSStefano Zampini     if (Ain->symmetric_set) {
4233301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->symmetric;
4243301b35fSStefano Zampini     }
4253301b35fSStefano Zampini #else
4263301b35fSStefano Zampini     if (Ain->hermitian_set) {
4273301b35fSStefano Zampini       sub_schurs->is_hermitian = Ain->hermitian;
4283301b35fSStefano Zampini     }
4293301b35fSStefano Zampini #endif
4303301b35fSStefano Zampini     if (Ain->spd_set) {
4313301b35fSStefano Zampini       sub_schurs->is_posdef = Ain->spd;
4323301b35fSStefano Zampini     }
433a64f4aa4SStefano Zampini 
4343301b35fSStefano Zampini     /* check */
435a64f4aa4SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
4361dac5d76SStefano Zampini     if (compute_Stilda && (!Ain->hermitian_set || !Ain->spd_set)) { /* these are lazy checks, maybe I should throw an error if they are not set */
4373301b35fSStefano Zampini       PetscInt lsize;
4383301b35fSStefano Zampini 
4393301b35fSStefano Zampini       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
4403301b35fSStefano Zampini       if (lsize) {
4413301b35fSStefano Zampini         PetscScalar val;
4423301b35fSStefano Zampini         PetscReal   norm;
4433301b35fSStefano Zampini         Vec         vec1,vec2,vec3;
4443301b35fSStefano Zampini 
4453301b35fSStefano Zampini         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
4463301b35fSStefano Zampini         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
4473301b35fSStefano Zampini         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
4483301b35fSStefano Zampini         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
4493301b35fSStefano Zampini #if !defined(PETSC_USE_COMPLEX)
4503301b35fSStefano Zampini         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
4513301b35fSStefano Zampini #else
4523301b35fSStefano Zampini         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
4533301b35fSStefano Zampini #endif
4543301b35fSStefano Zampini         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
4553301b35fSStefano Zampini         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
4561dac5d76SStefano Zampini         if (!Ain->hermitian_set) {
457862806e4SStefano Zampini           if (norm > PetscSqrtReal(PETSC_SMALL)) {
4583301b35fSStefano Zampini             sub_schurs->is_hermitian = PETSC_FALSE;
4593301b35fSStefano Zampini           } else {
4603301b35fSStefano Zampini             sub_schurs->is_hermitian = PETSC_TRUE;
4613301b35fSStefano Zampini           }
4621dac5d76SStefano Zampini         }
4631dac5d76SStefano Zampini         if (!Ain->spd_set && !benign_trick) {
4643301b35fSStefano Zampini           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
4653301b35fSStefano Zampini           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
46606a4e24aSStefano Zampini         }
4673301b35fSStefano Zampini         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
4683301b35fSStefano Zampini         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
4693301b35fSStefano Zampini         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
4703301b35fSStefano Zampini       } else {
4713301b35fSStefano Zampini         sub_schurs->is_hermitian = PETSC_TRUE;
4723301b35fSStefano Zampini         sub_schurs->is_posdef = PETSC_TRUE;
4733301b35fSStefano Zampini       }
4743301b35fSStefano Zampini     }
475a64f4aa4SStefano Zampini     if (isseqaij) {
476a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
477a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4783301b35fSStefano Zampini     } else {
479a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
480a64f4aa4SStefano Zampini     }
481a64f4aa4SStefano Zampini   }
4823301b35fSStefano Zampini 
483a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
484a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
485df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
486df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
487a64f4aa4SStefano Zampini   }
488a64f4aa4SStefano Zampini 
4895a95e1ceSStefano Zampini   /* preliminary checks */
490af25d912SStefano 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");
4915a95e1ceSStefano Zampini 
4925a95e1ceSStefano Zampini   /* restrict work on active processes */
4935a95e1ceSStefano Zampini   color = 0;
4945a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4955a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4965a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4975a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4985a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4995a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
5005a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
5015a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
5025a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
5035a95e1ceSStefano Zampini     PetscFunctionReturn(0);
5045a95e1ceSStefano Zampini   }
50581ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
5065a95e1ceSStefano Zampini 
507b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
508df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
509a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
5103301b35fSStefano Zampini     PetscBool isseqsbaij;
511a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
5123301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
5133301b35fSStefano Zampini     if (isseqsbaij) {
514a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
515a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
516a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
517a64f4aa4SStefano Zampini     } else {
518a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
519a64f4aa4SStefano Zampini       A_BB = tA_BB;
520a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
521a64f4aa4SStefano Zampini       A_IB = tA_IB;
522a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
523a64f4aa4SStefano Zampini       A_BI = tA_BI;
524f11841e3SStefano Zampini     }
525a58a30b4SStefano Zampini   } else {
5265a95e1ceSStefano Zampini     A_II = NULL;
5275a95e1ceSStefano Zampini     A_IB = NULL;
5285a95e1ceSStefano Zampini     A_BI = NULL;
5295a95e1ceSStefano Zampini     A_BB = NULL;
530b1b3d7a2SStefano Zampini   }
5315a95e1ceSStefano Zampini   S_all = NULL;
532b1b3d7a2SStefano Zampini 
533b1b3d7a2SStefano Zampini   /* determine interior problems */
5343dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
5353dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
536b1b3d7a2SStefano Zampini     PetscBT                touched;
537b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
538b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
539b1b3d7a2SStefano Zampini 
5406c4ed002SBarry Smith     if (!xadj || !adjncy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
541b1b3d7a2SStefano Zampini     /* get sizes */
542b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
543b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
544b1b3d7a2SStefano Zampini 
545b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
546b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
547b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
548b1b3d7a2SStefano Zampini 
549b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
550b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
551b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
552b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
553b1b3d7a2SStefano Zampini     }
554b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
555b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
556b1b3d7a2SStefano Zampini 
557b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
558b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
559b1b3d7a2SStefano Zampini     n_prev_added = n_B;
560b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
561b1b3d7a2SStefano Zampini       PetscInt n_added;
562b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5636c4ed002SBarry 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);
564b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
565b1b3d7a2SStefano Zampini       n_prev_added = n_added;
566b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
567b1b3d7a2SStefano Zampini       if (!n_added) break;
568b1b3d7a2SStefano Zampini     }
569b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
570b1b3d7a2SStefano Zampini 
571883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
572a9b99552SStefano 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);
573b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
574a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
575883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
576df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
577b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
578b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
579a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
580b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
581b1b3d7a2SStefano Zampini 
582b1b3d7a2SStefano Zampini       /* II block */
583b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
584b1b3d7a2SStefano Zampini     }
585b1b3d7a2SStefano Zampini   } else {
586b1b3d7a2SStefano Zampini     PetscInt n_I;
587b1b3d7a2SStefano Zampini 
588b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
589b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
590a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
591b1b3d7a2SStefano Zampini 
592b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
593df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
594b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
595b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
596b1b3d7a2SStefano Zampini 
597b1b3d7a2SStefano Zampini       /* II block is the same */
598b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
599b1b3d7a2SStefano Zampini       AE_II = A_II;
600b1b3d7a2SStefano Zampini     }
601b1b3d7a2SStefano Zampini   }
6025a95e1ceSStefano Zampini 
603883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
6045a95e1ceSStefano Zampini   max_subset_size = 0;
605883469d8SStefano Zampini   local_size = 0;
6065a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6075a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6085a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
609883469d8SStefano Zampini     local_size += subset_size;
610883469d8SStefano Zampini   }
611883469d8SStefano Zampini 
612883469d8SStefano Zampini   /* Work arrays for local indices */
613883469d8SStefano Zampini   extra = 0;
614683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
615df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
616a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
617883469d8SStefano Zampini   }
618683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
619883469d8SStefano Zampini   if (extra) {
620883469d8SStefano Zampini     const PetscInt *idxs;
621a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
622883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
623a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
624883469d8SStefano Zampini   }
625883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
626dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
627dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
628883469d8SStefano Zampini 
629883469d8SStefano Zampini   /* Get local indices in local numbering */
630883469d8SStefano Zampini   local_size = 0;
6315a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
632883469d8SStefano Zampini     PetscInt j;
633883469d8SStefano Zampini     const    PetscInt *idxs;
634883469d8SStefano Zampini 
6355a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
6365a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
637eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
638eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
639eb595f79SStefano Zampini     auxnum2[i] = subset_size;
640883469d8SStefano Zampini     /* subset indices in local numbering */
641883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
6425a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
643883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
644883469d8SStefano Zampini     local_size += subset_size;
645883469d8SStefano Zampini   }
646883469d8SStefano Zampini 
6475a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
648d2627357SStefano Zampini   Bwork = NULL;
649d2627357SStefano Zampini   pivots = NULL;
65006a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
651d2627357SStefano Zampini     PetscScalar lwork;
652d2627357SStefano Zampini 
65306a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
654d2627357SStefano Zampini     B_lwork = -1;
655d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
656d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
657d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
658d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
659d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
660d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
661d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
662d2627357SStefano Zampini   }
663d2627357SStefano Zampini 
664d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
665dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
666dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
667dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
668dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
669dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
670dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
671dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
672dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6736c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
674dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
675e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
676d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
677d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
678d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
679d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6802972d61bSStefano Zampini 
6815a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6825a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6835a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6845a95e1ceSStefano Zampini 
6855a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6865a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
6876c4ed002SBarry 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);
6885a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
689b1b3d7a2SStefano Zampini   }
690b1b3d7a2SStefano Zampini 
69172b8c272SStefano Zampini   if (change) {
69272b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
69372b8c272SStefano Zampini     IS                     change_primal_B;
69472b8c272SStefano Zampini     IS                     change_primal_all;
69572b8c272SStefano Zampini 
696b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
697b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
698b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
69972b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
70072b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
70172b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
702b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
70372b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
70472b8c272SStefano Zampini     }
70572b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
70672b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
70772b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
70872b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
70972b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
710b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
71172b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
71272b8c272SStefano Zampini       Mat change_sub;
71372b8c272SStefano Zampini 
714b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
71572b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
71672b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
71772b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
71872b8c272SStefano Zampini         ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
71972b8c272SStefano Zampini       } else {
72072b8c272SStefano Zampini         Mat change_subt;
72172b8c272SStefano Zampini         ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
72272b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
72372b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
72472b8c272SStefano Zampini       }
72572b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
72672b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
72772b8c272SStefano Zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
72872b8c272SStefano Zampini     }
72972b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
73072b8c272SStefano Zampini   }
73172b8c272SStefano Zampini 
7325a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
7335a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
7345a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
7355a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
7365a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
7375a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
738aa83b6aeSStefano Zampini   }
739b1b3d7a2SStefano Zampini 
7405a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
741be83ff47SStefano Zampini   F = NULL;
742df4d28bfSStefano 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 */
7435a95e1ceSStefano Zampini     Mat         S_Ej_expl;
7445a95e1ceSStefano Zampini     PetscScalar *work;
7455a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
7465a95e1ceSStefano Zampini     PetscBool   Sdense;
7475a95e1ceSStefano Zampini 
7485a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
7495a95e1ceSStefano Zampini     local_size = 0;
750b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7515a95e1ceSStefano Zampini       IS  is_subset_B;
7525a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
7535a95e1ceSStefano Zampini 
7545a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7555a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7565a95e1ceSStefano Zampini       /* EE block */
7575a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7585a95e1ceSStefano Zampini       /* IE block */
7595a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7605a95e1ceSStefano Zampini       /* EI block */
7615a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
7625a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7635a95e1ceSStefano Zampini       } else {
7645a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7655a95e1ceSStefano Zampini       }
766a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7675a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7685a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7695a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7705a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
771b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
772b1b3d7a2SStefano Zampini         KSP ksp;
773b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7745a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
775b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
776b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
777b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
778b1b3d7a2SStefano Zampini         KSPType   ksp_type;
779b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7805a95e1ceSStefano Zampini         PetscBool ispcnone;
781b1b3d7a2SStefano Zampini 
782b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7835a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
784b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
785b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
786b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
787b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7885a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7895a95e1ceSStefano Zampini         if (!ispcnone) {
7905a95e1ceSStefano Zampini           PCType pc_type;
791b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
792b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7935a95e1ceSStefano Zampini         } else {
7945a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7955a95e1ceSStefano Zampini         }
796b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
797b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
798b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
799b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
800b1b3d7a2SStefano Zampini           if (solver) {
801b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
802b1b3d7a2SStefano Zampini           }
803b1b3d7a2SStefano Zampini         }
804b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
805b1b3d7a2SStefano Zampini       }
8065a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
8075a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
8085a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
8095a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
8105a95e1ceSStefano Zampini       if (Sdense) {
8115a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
8125a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
813b1b3d7a2SStefano Zampini         }
8145a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
8156c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
8165a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
817a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
8185a95e1ceSStefano Zampini       local_size += subset_size;
8195a95e1ceSStefano Zampini     }
8205a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
821b1b3d7a2SStefano Zampini     /* free */
822b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
823b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
8245a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
825883469d8SStefano Zampini   } else {
8265cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
8279d54b7f4SStefano Zampini     Vec         Dall = NULL;
828ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
8295cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
8306dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
831df4d28bfSStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
832d4933d67SStefano Zampini     char        solver[256];
833883469d8SStefano Zampini 
834683d3df6SStefano Zampini     /* get sizes */
83581ea8064SStefano Zampini     n_I = 0;
83681ea8064SStefano Zampini     if (is_I_layer) {
837a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
83881ea8064SStefano Zampini     }
839683d3df6SStefano Zampini     economic = PETSC_FALSE;
840683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
841683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
842683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
8439d54b7f4SStefano Zampini     size_active_schur = local_size;
8449d54b7f4SStefano Zampini 
8459d54b7f4SStefano Zampini     /* import scaling vector */
8469d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
8479d54b7f4SStefano Zampini       const PetscScalar *array;
8489d54b7f4SStefano Zampini       PetscScalar       *array2;
8499d54b7f4SStefano Zampini       const PetscInt    *idxs;
8509d54b7f4SStefano Zampini       PetscInt          i;
8519d54b7f4SStefano Zampini 
8529d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8539d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8549d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8559d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8569d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8579d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8589d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8599d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8609d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8619d54b7f4SStefano Zampini     }
862d62866d3SStefano Zampini 
863683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
864683d3df6SStefano Zampini     cum = n_I+size_active_schur;
865683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
866683d3df6SStefano Zampini       const PetscInt* idxs;
867683d3df6SStefano Zampini       PetscInt        n_dir;
868683d3df6SStefano Zampini 
869683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
870683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
871683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
872683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
873683d3df6SStefano Zampini       cum += n_dir;
874d62866d3SStefano Zampini     }
875683d3df6SStefano Zampini     factor_workaround = PETSC_FALSE;
876683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
877367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
878683d3df6SStefano Zampini       PetscInt n_v;
879683d3df6SStefano Zampini 
880683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
881683d3df6SStefano Zampini       if (n_v) {
882683d3df6SStefano Zampini         const PetscInt* idxs;
883683d3df6SStefano Zampini 
884683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
885683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
886683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
887683d3df6SStefano Zampini         cum += n_v;
888683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
889683d3df6SStefano Zampini       }
890683d3df6SStefano Zampini     }
891683d3df6SStefano Zampini     size_schur = cum - n_I;
892683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
893683d3df6SStefano Zampini     /* get working mat (TODO: factorize without actually permuting)  */
894683d3df6SStefano Zampini     if (cum == n) {
895683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
896683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
897683d3df6SStefano Zampini     } else {
8986816873aSStefano Zampini       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
899683d3df6SStefano Zampini     }
900a9b99552SStefano Zampini     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
901ca92afb2SStefano Zampini 
902ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
903367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
904ca92afb2SStefano Zampini     if (benign_n) {
905ca92afb2SStefano Zampini       Vec                    v;
906ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
907ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
9085cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
9095cbda25cSStefano Zampini       PetscInt               sizeA;
910ca92afb2SStefano Zampini 
911ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
912ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
913ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
914ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
915ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
9165cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
9175cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
9185cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,benign_n,size_schur,NULL,&cs_AIB_mat);CHKERRQ(ierr);
9195cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9205cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
921ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
922ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
923ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
9245cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
925ca92afb2SStefano Zampini         PetscScalar    *array;
926ca92afb2SStefano Zampini         const PetscInt *idxs;
927ca92afb2SStefano Zampini         PetscInt       j,nz;
928ca92afb2SStefano Zampini 
929ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
930ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
931ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9325cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
933ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
9345cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
9355cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
9365cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
937ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
9385cbda25cSStefano Zampini         for (j=0;j<size_schur;j++) cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
939ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
940ca92afb2SStefano Zampini       }
9415cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9425cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
943ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
944ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
945ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
946ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
947ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
948ca92afb2SStefano Zampini     }
9493301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
9503301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9513301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
952df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
953d4933d67SStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
954df4d28bfSStefano Zampini #else
955df4d28bfSStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr);
956df4d28bfSStefano Zampini #endif
9578760537fSStefano Zampini     ierr = PetscOptionsGetString(NULL,"sub_schurs_","-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
958883469d8SStefano Zampini 
959683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
960d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
961d47842beSStefano Zampini 
962683d3df6SStefano Zampini     if (n_I) { /* TODO add ordering from options */
9635a05ddb0SStefano Zampini       IS is_schur;
9645a05ddb0SStefano Zampini 
9659ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
966d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
967883469d8SStefano Zampini       } else {
968d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
969883469d8SStefano Zampini       }
970883469d8SStefano Zampini       /* subsets ordered last */
9715a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9725a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9735a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
974883469d8SStefano Zampini 
975883469d8SStefano Zampini       /* factorization step */
9769ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
977883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
978be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
979be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
980be83ff47SStefano Zampini #endif
981e8ade678SStefano Zampini         if (sub_schurs->is_posdef) {
982e8ade678SStefano Zampini           ierr = MatFactorSetSchurComplementSolverType(F,1);CHKERRQ(ierr);
983e8ade678SStefano Zampini         } else {
984e8ade678SStefano Zampini           ierr = MatFactorSetSchurComplementSolverType(F,2);CHKERRQ(ierr);
985e8ade678SStefano Zampini         }
986883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
987a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
988883469d8SStefano Zampini       } else {
989883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
990be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
991be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
992be83ff47SStefano Zampini #endif
993883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
994883469d8SStefano Zampini       }
995883469d8SStefano Zampini 
996883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
9975a05ddb0SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
998d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
999683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1000683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1001df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
1002ca92afb2SStefano Zampini 
100372b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
1004ca92afb2SStefano Zampini       if (benign_n) {
10055cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1006ca92afb2SStefano Zampini         Vec         v;
10075cbda25cSStefano Zampini         PetscInt    sizeA;
1008ca92afb2SStefano Zampini 
1009ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
10105cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
10115cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1012ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1013ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1014ca92afb2SStefano Zampini #endif
1015ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1016ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1017ca92afb2SStefano Zampini #endif
10185cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10195cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1020ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
10215cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
1022*47484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
1023ca92afb2SStefano Zampini           const PetscInt *idxs;
1024ca92afb2SStefano Zampini           PetscInt       j,nz;
1025*47484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
1026ca92afb2SStefano Zampini 
10275cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
10285cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
10295cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1030ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1031ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10325cbda25cSStefano Zampini           for (j=0;j<nz;j++) sum += AIIm1_data[idxs[j]+sizeA*i];
1033ca92afb2SStefano Zampini           sum -= 1.; /* p0 dof (eliminated) is excluded from the sum */
1034ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10355cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1036ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
1037*47484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1038*47484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
10395cbda25cSStefano Zampini           sum  = -sum;
1040*47484b83SStefano Zampini           B_k = 1;
1041*47484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
1042*47484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
10435cbda25cSStefano Zampini           sum  = 1.;
1044*47484b83SStefano 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));
1045ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10465cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10475cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10485cbda25cSStefano Zampini           AIIm1_data[idxs[nz-1]+sizeA*i] = 0.;
10495cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10505cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1051ca92afb2SStefano Zampini         }
10525cbda25cSStefano Zampini   /* restore defaults */
10535cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
10545cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
10555cbda25cSStefano Zampini #endif
10565cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
10575cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
10585cbda25cSStefano Zampini #endif
10595cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10605cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1061ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1062ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1063ca92afb2SStefano Zampini       }
1064a3df083aSStefano Zampini       if (!reuse_solvers) {
1065a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1066a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1067a3df083aSStefano Zampini         }
1068a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
10695cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
10705cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1071a3df083aSStefano Zampini       }
1072df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1073be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1074683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1075df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1076be83ff47SStefano Zampini     }
1077be83ff47SStefano Zampini 
1078be83ff47SStefano Zampini     if (reuse_solvers) {
1079a00504b5SStefano Zampini       Mat                A_II,Afake;
108053892102SStefano Zampini       Vec                vec1_B;
1081df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
10823462e049SStefano Zampini       PetscInt           n_R;
1083d5574798SStefano Zampini 
1084df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1085df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1086e28d306cSStefano Zampini       } else {
1087df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1088d62866d3SStefano Zampini       }
1089df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1090be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1091d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1092d5574798SStefano Zampini       msolv_ctx->F = F;
1093683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1094683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1095683d3df6SStefano Zampini       {
1096683d3df6SStefano Zampini         PetscScalar *array;
1097683d3df6SStefano Zampini         PetscInt    n;
1098683d3df6SStefano Zampini 
1099683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1100683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1101683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1102683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1103683d3df6SStefano Zampini       }
1104683d3df6SStefano Zampini       msolv_ctx->has_vertices = factor_workaround;
1105d62866d3SStefano Zampini 
1106d62866d3SStefano Zampini       /* interior solver */
1107683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1108d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1109d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1110d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1111df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1112df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1113d62866d3SStefano Zampini 
1114d62866d3SStefano Zampini       /* correction solver */
11153462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1116d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1117d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1118df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1119df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
112053892102SStefano Zampini 
112153892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
112253892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
112353892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1124683d3df6SStefano Zampini       if (!factor_workaround) {
112553892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
112653892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
112753892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
112853892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1129683d3df6SStefano Zampini       } else {
1130683d3df6SStefano Zampini         IS              is_B_all;
1131683d3df6SStefano Zampini         const PetscInt* idxs;
1132683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1133683d3df6SStefano Zampini 
1134683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1135683d3df6SStefano Zampini         dual = size_schur - n_v;
1136683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1137683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1138683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1139683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1140683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1141683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1142683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1143683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1144683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1145683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1146683d3df6SStefano Zampini       }
11473462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
11483462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
11493462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
11503462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1151683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1152ca92afb2SStefano Zampini 
1153ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1154ca92afb2SStefano Zampini       if (benign_n) {
11555cbda25cSStefano Zampini         PetscScalar *array;
11565cbda25cSStefano Zampini 
1157ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1158ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1159ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
11605cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
11615cbda25cSStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,NULL,&msolv_ctx->benign_corr_work);CHKERRQ(ierr);
11625cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11635cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
11645cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11655cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1166ca92afb2SStefano Zampini       }
1167d5574798SStefano Zampini     }
116808122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
116953892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
11705db18549SStefano Zampini 
1171be83ff47SStefano Zampini     /* Work arrays */
1172be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
11735a94c880SStefano Zampini 
11745a94c880SStefano Zampini     /* matrices for adaptive selection */
117512d906b1SStefano Zampini     if (compute_Stilda) {
11765a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
11775a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11785a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
11795a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
11805a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
11815a94c880SStefano Zampini       }
11829d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
11835a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
11845a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
11855a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
11865a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
11875a94c880SStefano Zampini       }
118812d906b1SStefano Zampini     }
1189d2627357SStefano Zampini 
1190be83ff47SStefano Zampini     /* S_Ej_all */
1191be83ff47SStefano Zampini     cum = cum2 = 0;
1192be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
119365d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
119465d8bf0aSStefano Zampini       PetscInt j;
119565d8bf0aSStefano Zampini 
11965a95e1ceSStefano Zampini       /* get S_E */
1197b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1198683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1199be83ff47SStefano Zampini         PetscInt k;
1200be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1201be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1202be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1203683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1204be83ff47SStefano Zampini           }
1205be83ff47SStefano Zampini         }
120606a4e24aSStefano Zampini       } else { /* just copy to workspace */
1207be83ff47SStefano Zampini         PetscInt k;
1208be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1209be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1210be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1211be83ff47SStefano Zampini           }
1212be83ff47SStefano Zampini         }
12139087bf02SStefano Zampini       }
12145a95e1ceSStefano Zampini       /* insert S_E values */
1215be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1216b7ab4a40SStefano Zampini       if (sub_schurs->change) {
12178760537fSStefano Zampini         Mat            change_sub,SEj,T;
12188760537fSStefano Zampini         const PetscInt *idxs;
12198760537fSStefano Zampini         PetscInt       nz;
122072b8c272SStefano Zampini 
122172b8c272SStefano Zampini         /* change basis */
122272b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
122372b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12248760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
12258760537fSStefano Zampini           Mat T2;
12268760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
12278760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
12288760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
122972b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
12308760537fSStefano Zampini         } else {
12318760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
123272b8c272SStefano Zampini         }
123372b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
123472b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
12358760537fSStefano Zampini         /* currently there's no support for ZeroRowsColumns with A dense */
1236b7ab4a40SStefano Zampini         ierr = ISGetIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
1237b7ab4a40SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->change_primal_sub[i],&nz);CHKERRQ(ierr);
12388760537fSStefano Zampini         for (j=0;j<nz;j++) {
12398760537fSStefano Zampini           PetscInt k;
12408760537fSStefano Zampini 
12418760537fSStefano Zampini           ierr = PetscMemzero(work+subset_size*idxs[j],subset_size*sizeof(PetscScalar));CHKERRQ(ierr);
12428760537fSStefano Zampini           for (k=0;k<subset_size;k++) work[idxs[j]+k*subset_size] = 0.;
12438760537fSStefano Zampini           work[idxs[j]*(subset_size+1)] = 1.0;
124472b8c272SStefano Zampini         }
1245b7ab4a40SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
124672b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
124772b8c272SStefano Zampini       }
12489d54b7f4SStefano Zampini       if (deluxe) {
12495a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1250683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1251862806e4SStefano Zampini         if (compute_Stilda) {
125208122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
125308122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
125406a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
12556c3e6151SStefano Zampini             PetscInt k;
12565a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12572972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
12585a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
12592972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
12606c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
12616c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
12626c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
12636c3e6151SStefano Zampini               }
12646c3e6151SStefano Zampini             }
1265d6462365SStefano Zampini           } else {
1266d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1267d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1268d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1269d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
12702972d61bSStefano Zampini           }
127108122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
12725a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12739087bf02SStefano Zampini         }
12749d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
12759d54b7f4SStefano Zampini         Mat         SEj;
12769d54b7f4SStefano Zampini         Vec         D;
12779d54b7f4SStefano Zampini         PetscScalar *array;
12789d54b7f4SStefano Zampini 
12799d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12809d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
12819d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
12829d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
12839d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
12849d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
12859d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
12869d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12879d54b7f4SStefano Zampini       }
1288be83ff47SStefano Zampini       cum += subset_size;
1289be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1290be83ff47SStefano Zampini     }
1291be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
12924a6c6b0dSStefano Zampini 
1293df4d28bfSStefano Zampini     if (solver_S) {
12945a05ddb0SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
12954a6c6b0dSStefano Zampini     }
1296683d3df6SStefano Zampini 
1297683d3df6SStefano Zampini     schur_factor = NULL;
129845951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1299683d3df6SStefano Zampini 
13009d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
13014a6c6b0dSStefano Zampini         PetscInt j;
13024a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
13035a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
13044a6c6b0dSStefano Zampini       } else {
1305683d3df6SStefano Zampini         Mat S_all_inv=NULL;
1306df4d28bfSStefano Zampini         if (solver_S) { /* use MatFactor calls to invert S */
130772b8c272SStefano Zampini             /* let MatFactor factor the Schur complement */
13086dba178dSStefano Zampini           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
1309683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1310683d3df6SStefano 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 */
1311683d3df6SStefano Zampini           if (factor_workaround) {
1312683d3df6SStefano Zampini             PetscScalar *data;
1313683d3df6SStefano Zampini             PetscInt nd = 0;
13146dba178dSStefano Zampini 
1315683d3df6SStefano Zampini             if (!sub_schurs->is_posdef) {
1316683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1317683d3df6SStefano Zampini             }
13186dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1319683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1320683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1321683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1322683d3df6SStefano Zampini             }
1323683d3df6SStefano Zampini             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
1324683d3df6SStefano Zampini             cum = 0;
1325683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1326683d3df6SStefano Zampini               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1327683d3df6SStefano Zampini               cum += size_active_schur-i;
1328683d3df6SStefano Zampini             }
1329683d3df6SStefano Zampini             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
13306dba178dSStefano Zampini             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1331683d3df6SStefano Zampini             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1332683d3df6SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1333683d3df6SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
1334683d3df6SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1335683d3df6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1336683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1337683d3df6SStefano Zampini           } else {
13386dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
13396dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1340683d3df6SStefano Zampini           }
1341df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1342683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1343683d3df6SStefano Zampini           S_all_inv = S_all;
1344683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1345be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1346be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
134706a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1348be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1349be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1350be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1351be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1352d6462365SStefano Zampini           } else {
1353d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1354d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1355d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1356d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1357be83ff47SStefano Zampini           }
1358be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1359683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1360be83ff47SStefano Zampini         }
1361be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1362be83ff47SStefano Zampini         cum = cum2 = 0;
1363683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1364be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1365be83ff47SStefano Zampini           PetscInt j;
1366862806e4SStefano Zampini 
1367862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1368be83ff47SStefano Zampini           /* get (St^-1)_E */
136972b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
137006a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1371a0b0af32SStefano Zampini           if (S_lower_triangular) {
1372be83ff47SStefano Zampini             PetscInt k;
1373b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1374be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1375be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1376be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
13776c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1378be83ff47SStefano Zampini                 }
1379be83ff47SStefano Zampini               }
138072b8c272SStefano Zampini             } else {
138172b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
138272b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
138372b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
138472b8c272SStefano Zampini                 }
138572b8c272SStefano Zampini               }
138672b8c272SStefano Zampini             }
138772b8c272SStefano Zampini           } else {
1388be83ff47SStefano Zampini             PetscInt k;
1389be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1390be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1391be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1392be83ff47SStefano Zampini               }
1393be83ff47SStefano Zampini             }
1394be83ff47SStefano Zampini           }
1395b7ab4a40SStefano Zampini           if (sub_schurs->change) {
13968760537fSStefano Zampini             Mat change_sub,SEj,T;
139772b8c272SStefano Zampini 
139872b8c272SStefano Zampini             /* change basis */
139972b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
140072b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
14018760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
14028760537fSStefano Zampini               Mat T2;
14038760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
14048760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
140572b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
14068760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
14078760537fSStefano Zampini             } else {
14088760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
140972b8c272SStefano Zampini             }
14108760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
14118760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
141272b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
141372b8c272SStefano Zampini             {
141472b8c272SStefano Zampini               const PetscInt *idxs;
141572b8c272SStefano Zampini               PetscInt n,j;
1416b7ab4a40SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->change_primal_sub[i],&n);CHKERRQ(ierr);
1417b7ab4a40SStefano Zampini               ierr = ISGetIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
141872b8c272SStefano Zampini               for (j=0;j<n;j++) {
141972b8c272SStefano Zampini                 PetscInt k;
142072b8c272SStefano Zampini                 for (k=0;k<subset_size;k++) {
142172b8c272SStefano Zampini                   work[k+idxs[j]*subset_size] = work[idxs[j]+k*subset_size] = 0.;
142272b8c272SStefano Zampini                 }
142372b8c272SStefano Zampini                 work[idxs[j] + idxs[j]*subset_size] = 1./PETSC_SMALL;
142472b8c272SStefano Zampini               }
1425b7ab4a40SStefano Zampini               ierr = ISRestoreIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
142672b8c272SStefano Zampini             }
142772b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
142872b8c272SStefano Zampini           }
1429be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
14305a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1431be83ff47SStefano Zampini           cum += subset_size;
1432be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1433883469d8SStefano Zampini         }
1434683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1435df4d28bfSStefano Zampini         if (solver_S) {
14366dba178dSStefano Zampini           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
14375db18549SStefano Zampini         }
1438683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1439683d3df6SStefano Zampini       }
1440683d3df6SStefano Zampini 
1441683d3df6SStefano Zampini       /* move back factors to Schur data into F */
1442683d3df6SStefano Zampini       if (factor_workaround) {
1443683d3df6SStefano Zampini         Mat S_tmp;
1444683d3df6SStefano Zampini         PetscInt nv,nd=0;
1445683d3df6SStefano Zampini 
1446df4d28bfSStefano Zampini         if (!solver_S) {
1447683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1448683d3df6SStefano Zampini         }
1449683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
14506dba178dSStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1451683d3df6SStefano Zampini         if (sub_schurs->is_posdef) {
1452683d3df6SStefano Zampini           PetscScalar *data;
1453683d3df6SStefano Zampini 
1454683d3df6SStefano Zampini           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1455683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1456683d3df6SStefano Zampini 
1457683d3df6SStefano Zampini           if (S_lower_triangular) {
1458683d3df6SStefano Zampini             cum = 0;
1459683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1460683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1461683d3df6SStefano Zampini               cum += size_active_schur-i;
1462683d3df6SStefano Zampini 	    }
1463683d3df6SStefano Zampini           } else {
1464683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1465683d3df6SStefano Zampini           }
1466683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1467683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1468683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1469683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1470683d3df6SStefano Zampini 	    }
1471683d3df6SStefano Zampini           }
14726dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1473683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1474683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
14756c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1476683d3df6SStefano Zampini           }
1477683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1478683d3df6SStefano Zampini         } else {
1479683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1480683d3df6SStefano Zampini         }
14816dba178dSStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
14829087bf02SStefano Zampini       }
1483367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1484367aa537SStefano Zampini       PetscScalar *data;
1485367aa537SStefano Zampini       PetscInt    nd = 0;
1486367aa537SStefano Zampini 
1487367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1488367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1489367aa537SStefano Zampini       }
1490367aa537SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
1491367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1492367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1493367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1494367aa537SStefano Zampini       }
1495367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1496367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
14976c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1498367aa537SStefano Zampini       }
1499367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
1500367aa537SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
15014a6c6b0dSStefano Zampini     }
15024a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1503683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
15049d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
15054a6c6b0dSStefano Zampini   }
15065a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1507be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1508a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1509a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1510a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1511a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1512a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
15135db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15145db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15155a95e1ceSStefano Zampini   if (compute_Stilda) {
15165a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15175a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15189d54b7f4SStefano Zampini     if (deluxe) {
15195a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15205a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
152108122e43SStefano Zampini     }
15229d54b7f4SStefano Zampini   }
1523a1337663SStefano Zampini 
15245db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
15255db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
15263927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
15273927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15285a95e1ceSStefano Zampini 
15295db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
15305a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
1531dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
15325a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
15335a94c880SStefano Zampini   } else {
15345a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
15355a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
1536dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
15375a94c880SStefano Zampini   }
153808122e43SStefano Zampini 
1539f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
15405a95e1ceSStefano Zampini   if (compute_Stilda) {
15415a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1542a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15435a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1544dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
15459d54b7f4SStefano Zampini     if (deluxe) {
15465a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
154708122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15485a94c880SStefano Zampini       submats[0] = sub_schurs->sum_S_Ej_inv_all;
1549dc456d91SStefano Zampini       ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
15509d54b7f4SStefano Zampini     } else {
15519d54b7f4SStefano Zampini       PetscScalar *array;
15529d54b7f4SStefano Zampini       PetscInt    cum;
15539d54b7f4SStefano Zampini 
15549d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15559d54b7f4SStefano Zampini       cum = 0;
15569d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
15579d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
15589d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
15599d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
15609d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
15619d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
15629d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
15639d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
15649d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
15659d54b7f4SStefano Zampini         cum += subset_size*subset_size;
15669d54b7f4SStefano Zampini       }
15679d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15689d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
15699d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
15709d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
15719d54b7f4SStefano Zampini     }
157208122e43SStefano Zampini   }
15733202ece2SStefano Zampini 
15745a95e1ceSStefano Zampini   /* free workspace */
15755a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
157606a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1577a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
15783202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1579dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
15805a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1581b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1582b1b3d7a2SStefano Zampini }
1583b1b3d7a2SStefano Zampini 
1584b1b3d7a2SStefano Zampini #undef __FUNCT__
1585b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
1586a64f4aa4SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1587b1b3d7a2SStefano Zampini {
15889bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
15895a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1590b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1591b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1592b1b3d7a2SStefano Zampini 
1593b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1594b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
15956c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1596b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
15976c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1598b1b3d7a2SStefano Zampini 
1599b1b3d7a2SStefano Zampini   /* reset any previous data */
1600b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1601b1b3d7a2SStefano Zampini 
16025a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
16039bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1604b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
160508122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1606b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1607b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
1608b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
1609b1b3d7a2SStefano Zampini   }
1610b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
1611b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
161208122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1613b1b3d7a2SStefano Zampini   }
1614b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
1615b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
1616d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1617d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1618b1b3d7a2SStefano Zampini 
1619df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
1620df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_FALSE;
1621883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1622df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1623df4d28bfSStefano Zampini #endif
1624df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1625df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1626883469d8SStefano Zampini #endif
1627b1b3d7a2SStefano Zampini 
1628b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1629b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1630b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1631b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
16325db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
16335db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
16345db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
16355db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
16365a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1637b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1638df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1639b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1640b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1641b96c3477SStefano Zampini     }
16429bb4a8caSStefano Zampini   }
1643d62866d3SStefano Zampini   sub_schurs->is_vertices = vertices;
1644b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1645b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
164608122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1647b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1648b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1649b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1650b1b3d7a2SStefano Zampini }
1651b1b3d7a2SStefano Zampini 
1652b1b3d7a2SStefano Zampini #undef __FUNCT__
165334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
165434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
165534a97f8cSStefano Zampini {
165634a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
165734a97f8cSStefano Zampini   PetscErrorCode  ierr;
165834a97f8cSStefano Zampini 
165934a97f8cSStefano Zampini   PetscFunctionBegin;
166034a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
16615ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
166234a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
166334a97f8cSStefano Zampini   PetscFunctionReturn(0);
166434a97f8cSStefano Zampini }
166534a97f8cSStefano Zampini 
166634a97f8cSStefano Zampini #undef __FUNCT__
166734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
166834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
166934a97f8cSStefano Zampini {
167034a97f8cSStefano Zampini   PetscInt       i;
167134a97f8cSStefano Zampini   PetscErrorCode ierr;
167234a97f8cSStefano Zampini 
167334a97f8cSStefano Zampini   PetscFunctionBegin;
16741e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1675b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1676b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1677b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
16785db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
16795db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
168041c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
168141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
168208122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1683a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
16845db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1685d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1686d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
168708122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
168834a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1689b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
169034a97f8cSStefano Zampini   }
16915ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1692b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
16933dc780c3SStefano Zampini   }
1694df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1695df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1696d62866d3SStefano Zampini   }
1697df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
169872b8c272SStefano Zampini   if (sub_schurs->change) {
169972b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
170072b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1701b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
170272b8c272SStefano Zampini     }
170372b8c272SStefano Zampini   }
170472b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1705b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
170634a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
170734a97f8cSStefano Zampini   PetscFunctionReturn(0);
170834a97f8cSStefano Zampini }
170934a97f8cSStefano Zampini 
171034a97f8cSStefano Zampini #undef __FUNCT__
171134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
17122a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
171334a97f8cSStefano Zampini {
171434a97f8cSStefano Zampini   PetscInt       i,j,n;
171534a97f8cSStefano Zampini   PetscErrorCode ierr;
171634a97f8cSStefano Zampini 
171734a97f8cSStefano Zampini   PetscFunctionBegin;
171834a97f8cSStefano Zampini   n = 0;
171934a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
172034a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
172134a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
172234a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
172334a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
172434a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
172534a97f8cSStefano Zampini         queue_tip[n] = dof;
172634a97f8cSStefano Zampini         n++;
172734a97f8cSStefano Zampini       }
172834a97f8cSStefano Zampini     }
172934a97f8cSStefano Zampini   }
173034a97f8cSStefano Zampini   *n_added = n;
173134a97f8cSStefano Zampini   PetscFunctionReturn(0);
173234a97f8cSStefano Zampini }
1733