xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision e62b652183a3ad7bfed7e58f1c5d8560af14abc7)
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 */
25282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
265cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
275cbda25cSStefano Zampini     n_I = n_I - size_schur;
285cbda25cSStefano Zampini     /* get schur sol from array */
295cbda25cSStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
305cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
315cbda25cSStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
325cbda25cSStefano Zampini     /* apply interior sol correction */
33282d6408SStefano Zampini     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
345cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
355cbda25cSStefano Zampini     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
365cbda25cSStefano Zampini   }
37ca92afb2SStefano Zampini   if (v2) {
38ca92afb2SStefano Zampini     PetscInt nl;
39ca92afb2SStefano Zampini 
40ca92afb2SStefano Zampini     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
41ca92afb2SStefano Zampini     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
42ca92afb2SStefano Zampini     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
43ca92afb2SStefano Zampini     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
44ca92afb2SStefano Zampini   } else {
45ca92afb2SStefano Zampini     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
46ca92afb2SStefano Zampini     array2 = array;
47ca92afb2SStefano Zampini   }
48ca92afb2SStefano Zampini   if (!sol) { /* change rhs */
49ca92afb2SStefano Zampini     PetscInt n;
50ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
51ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
52ca92afb2SStefano Zampini       const PetscInt *cols;
53ca92afb2SStefano Zampini       PetscInt       nz,i;
54ca92afb2SStefano Zampini 
55ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
56ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
57ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
5822db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
5922db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
6022db5ddcSStefano Zampini #else
61ca92afb2SStefano Zampini       sum = -sum/nz;
6222db5ddcSStefano Zampini #endif
63ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
64ca92afb2SStefano Zampini       ctx->benign_save_vals[n] = array2[cols[nz-1]];
65ca92afb2SStefano Zampini       array2[cols[nz-1]] = sum;
66ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
67ca92afb2SStefano Zampini     }
68ca92afb2SStefano Zampini   } else {
69ca92afb2SStefano Zampini     PetscInt n;
70ca92afb2SStefano Zampini     for (n=0;n<ctx->benign_n;n++) {
71ca92afb2SStefano Zampini       PetscScalar    sum = 0.;
72ca92afb2SStefano Zampini       const PetscInt *cols;
73ca92afb2SStefano Zampini       PetscInt       nz,i;
74ca92afb2SStefano Zampini       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
75ca92afb2SStefano Zampini       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
76ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) sum += array[cols[i]];
7722db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
7822db5ddcSStefano Zampini       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
7922db5ddcSStefano Zampini #else
80ca92afb2SStefano Zampini       sum = -sum/nz;
8122db5ddcSStefano Zampini #endif
82ca92afb2SStefano Zampini       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
83ca92afb2SStefano Zampini       array2[cols[nz-1]] = ctx->benign_save_vals[n];
84ca92afb2SStefano Zampini       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
85ca92afb2SStefano Zampini     }
86ca92afb2SStefano Zampini   }
87ca92afb2SStefano Zampini   if (v2) {
88ca92afb2SStefano Zampini     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
8972b8c272SStefano Zampini     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
90ca92afb2SStefano Zampini   } else {
91ca92afb2SStefano Zampini     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
92ca92afb2SStefano Zampini   }
935cbda25cSStefano Zampini   if (!sol && full) {
945cbda25cSStefano Zampini     Vec      usedv;
955cbda25cSStefano Zampini     PetscInt n_I,size_schur;
965cbda25cSStefano Zampini 
975cbda25cSStefano Zampini     /* get sizes */
98282d6408SStefano Zampini     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
995cbda25cSStefano Zampini     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
1005cbda25cSStefano Zampini     n_I = n_I - size_schur;
1015cbda25cSStefano Zampini     /* compute schur rhs correction */
1025cbda25cSStefano Zampini     if (v2) {
1035cbda25cSStefano Zampini       usedv = v2;
1045cbda25cSStefano Zampini     } else {
1055cbda25cSStefano Zampini       usedv = v;
1065cbda25cSStefano Zampini     }
1075cbda25cSStefano Zampini     /* apply schur rhs correction */
1085cbda25cSStefano Zampini     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
1095cbda25cSStefano Zampini     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
1105cbda25cSStefano Zampini     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
1115cbda25cSStefano Zampini     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
112282d6408SStefano Zampini     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1135cbda25cSStefano Zampini     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1145cbda25cSStefano Zampini   }
115ca92afb2SStefano Zampini   PetscFunctionReturn(0);
116ca92afb2SStefano Zampini }
117ca92afb2SStefano Zampini 
118d62866d3SStefano Zampini #undef __FUNCT__
119df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_Solve_Private"
120df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
121d62866d3SStefano Zampini {
122df4d28bfSStefano Zampini   PCBDDCReuseSolvers ctx;
123683d3df6SStefano Zampini   PetscBool          copy = PETSC_FALSE;
124d62866d3SStefano Zampini   PetscErrorCode     ierr;
125d62866d3SStefano Zampini 
126d62866d3SStefano Zampini   PetscFunctionBegin;
127d62866d3SStefano Zampini   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
128683d3df6SStefano Zampini   if (full) {
129d62866d3SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
130d62866d3SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
131d62866d3SStefano Zampini #endif
1325cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1335cbda25cSStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
1345cbda25cSStefano Zampini #endif
135683d3df6SStefano Zampini     copy = ctx->has_vertices;
136d4933d67SStefano Zampini   } else { /* interior solver */
1376dba178dSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
138683d3df6SStefano Zampini     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
1396dba178dSStefano Zampini #endif
140d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
141d4933d67SStefano Zampini     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
142d4933d67SStefano Zampini #endif
143683d3df6SStefano Zampini     copy = PETSC_TRUE;
144683d3df6SStefano Zampini   }
145683d3df6SStefano Zampini   /* copy rhs into factored matrix workspace */
146683d3df6SStefano Zampini   if (copy) {
147ca92afb2SStefano Zampini     PetscInt    n;
148df4d28bfSStefano Zampini     PetscScalar *array,*array_solver;
149ca92afb2SStefano Zampini 
150ca92afb2SStefano Zampini     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
151683d3df6SStefano Zampini     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
152df4d28bfSStefano Zampini     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
153df4d28bfSStefano Zampini     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
154df4d28bfSStefano Zampini     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
155683d3df6SStefano Zampini     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
156683d3df6SStefano Zampini 
1575cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
158683d3df6SStefano Zampini     if (transpose) {
159683d3df6SStefano Zampini       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
160683d3df6SStefano Zampini     } else {
161683d3df6SStefano Zampini       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
162683d3df6SStefano Zampini     }
1635cbda25cSStefano Zampini     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
164683d3df6SStefano Zampini 
165683d3df6SStefano Zampini     /* get back data to caller worskpace */
166df4d28bfSStefano Zampini     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
167683d3df6SStefano Zampini     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
168df4d28bfSStefano Zampini     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
169683d3df6SStefano Zampini     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
170df4d28bfSStefano Zampini     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
171683d3df6SStefano Zampini   } else {
172ca92afb2SStefano Zampini     if (ctx->benign_n) {
1735cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
174ca92afb2SStefano Zampini       if (transpose) {
175ca92afb2SStefano Zampini         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
176ca92afb2SStefano Zampini       } else {
177ca92afb2SStefano Zampini         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
178ca92afb2SStefano Zampini       }
1795cbda25cSStefano Zampini       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
180ca92afb2SStefano Zampini     } else {
181e28d306cSStefano Zampini       if (transpose) {
182e28d306cSStefano Zampini         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
183e28d306cSStefano Zampini       } else {
1846816873aSStefano Zampini         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
185e28d306cSStefano Zampini       }
186683d3df6SStefano Zampini     }
187ca92afb2SStefano Zampini   }
1885cbda25cSStefano Zampini   /* restore defaults */
1895cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1905cbda25cSStefano Zampini   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
1915cbda25cSStefano Zampini #endif
192d4933d67SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
193d4933d67SStefano Zampini   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
194d4933d67SStefano Zampini #endif
195d62866d3SStefano Zampini   PetscFunctionReturn(0);
196d62866d3SStefano Zampini }
197d62866d3SStefano Zampini 
198d62866d3SStefano Zampini #undef __FUNCT__
199df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_Correction"
200df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
201e28d306cSStefano Zampini {
202e28d306cSStefano Zampini   PetscErrorCode   ierr;
203e28d306cSStefano Zampini 
204e28d306cSStefano Zampini   PetscFunctionBegin;
205df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
206e28d306cSStefano Zampini   PetscFunctionReturn(0);
207e28d306cSStefano Zampini }
208e28d306cSStefano Zampini 
209e28d306cSStefano Zampini #undef __FUNCT__
210df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_CorrectionTranspose"
211df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
212e28d306cSStefano Zampini {
213e28d306cSStefano Zampini   PetscErrorCode   ierr;
214e28d306cSStefano Zampini 
215e28d306cSStefano Zampini   PetscFunctionBegin;
216df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
217683d3df6SStefano Zampini   PetscFunctionReturn(0);
218683d3df6SStefano Zampini }
219683d3df6SStefano Zampini 
220683d3df6SStefano Zampini #undef __FUNCT__
221df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_Interior"
222df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
223683d3df6SStefano Zampini {
224683d3df6SStefano Zampini   PetscErrorCode   ierr;
225683d3df6SStefano Zampini 
226683d3df6SStefano Zampini   PetscFunctionBegin;
227df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
228683d3df6SStefano Zampini   PetscFunctionReturn(0);
229683d3df6SStefano Zampini }
230683d3df6SStefano Zampini 
231683d3df6SStefano Zampini #undef __FUNCT__
232df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolvers_InteriorTranspose"
233df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
234683d3df6SStefano Zampini {
235683d3df6SStefano Zampini   PetscErrorCode   ierr;
236683d3df6SStefano Zampini 
237683d3df6SStefano Zampini   PetscFunctionBegin;
238df4d28bfSStefano Zampini   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
239e28d306cSStefano Zampini   PetscFunctionReturn(0);
240e28d306cSStefano Zampini }
241e28d306cSStefano Zampini 
242e28d306cSStefano Zampini #undef __FUNCT__
243df4d28bfSStefano Zampini #define __FUNCT__ "PCBDDCReuseSolversReset"
244df4d28bfSStefano Zampini static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
245d62866d3SStefano Zampini {
246ca92afb2SStefano Zampini   PetscInt       i;
247d62866d3SStefano Zampini   PetscErrorCode ierr;
248d62866d3SStefano Zampini 
249d62866d3SStefano Zampini   PetscFunctionBegin;
250d62866d3SStefano Zampini   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
251d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
252d62866d3SStefano Zampini   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
253d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
254d62866d3SStefano Zampini   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
25553892102SStefano Zampini   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
25653892102SStefano Zampini   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
257d62866d3SStefano Zampini   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
25853892102SStefano Zampini   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
25953892102SStefano Zampini   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
260ca92afb2SStefano Zampini   for (i=0;i<reuse->benign_n;i++) {
261ca92afb2SStefano Zampini     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
262ca92afb2SStefano Zampini   }
263ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
264ca92afb2SStefano Zampini   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
2655cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
2665cbda25cSStefano Zampini   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
2675cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
2685cbda25cSStefano Zampini   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
269d62866d3SStefano Zampini   PetscFunctionReturn(0);
270d62866d3SStefano Zampini }
271d5574798SStefano Zampini 
272d5574798SStefano Zampini #undef __FUNCT__
2733202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
2745ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
2753202ece2SStefano Zampini {
2763202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
2773202ece2SStefano Zampini   KSP            ksp;
2783202ece2SStefano Zampini   PC             pc;
2793202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
2803202ece2SStefano Zampini   PetscReal      fill = 2.0;
281f11841e3SStefano Zampini   PetscInt       n_I;
2823202ece2SStefano Zampini   PetscMPIInt    size;
2833202ece2SStefano Zampini   PetscErrorCode ierr;
2843202ece2SStefano Zampini 
2853202ece2SStefano Zampini   PetscFunctionBegin;
2863202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
2876c4ed002SBarry Smith   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
288f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
289f11841e3SStefano Zampini     PetscBool Sdense;
290f11841e3SStefano Zampini 
291f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
2926c4ed002SBarry Smith     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
293f11841e3SStefano Zampini   }
2943202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
2953202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
2963202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
2973202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
2983202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
2993202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
3003202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
3013202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
302f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
303f11841e3SStefano Zampini   if (n_I) {
3043202ece2SStefano Zampini     if (!Bdense) {
3053202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
3063202ece2SStefano Zampini     } else {
3073202ece2SStefano Zampini       Bd = B;
3083202ece2SStefano Zampini     }
3093202ece2SStefano Zampini 
3103202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
3113202ece2SStefano Zampini       Mat fact;
3123202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
3133202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
3143202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
3153202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
3163202ece2SStefano Zampini     } else {
31707b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
31807b1e237SStefano Zampini 
31907b1e237SStefano Zampini       if (ex) {
3203202ece2SStefano Zampini         Mat Ainvd;
3213202ece2SStefano Zampini 
3223202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
3233202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
3243202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
32507b1e237SStefano Zampini       } else {
32607b1e237SStefano Zampini         Vec         sol,rhs;
32707b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
32807b1e237SStefano Zampini         PetscInt    i,nrhs,n;
32907b1e237SStefano Zampini 
33007b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
33107b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
33207b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
33307b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
33407b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
33507b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
33607b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
33707b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
33807b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
33907b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
34007b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
34107b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
34207b1e237SStefano Zampini         }
34307b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
34407b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
34507b1e237SStefano Zampini       }
3463202ece2SStefano Zampini     }
3475ec10c6aSStefano Zampini     if (!Bdense & !issym) {
3483202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3493202ece2SStefano Zampini     }
3505ec10c6aSStefano Zampini 
3515ec10c6aSStefano Zampini     if (!issym) {
3523202ece2SStefano Zampini       if (!Cdense) {
3533202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
3543202ece2SStefano Zampini       } else {
3553202ece2SStefano Zampini         Cd = C;
3563202ece2SStefano Zampini       }
3575ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3583202ece2SStefano Zampini       if (!Cdense) {
3593202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
3603202ece2SStefano Zampini       }
3615ec10c6aSStefano Zampini     } else {
3625ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
3635ec10c6aSStefano Zampini       if (!Bdense) {
3645ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
3655ec10c6aSStefano Zampini       }
3665ec10c6aSStefano Zampini     }
3675ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
368f11841e3SStefano Zampini   }
3693202ece2SStefano Zampini 
3703202ece2SStefano Zampini   if (D) {
3713202ece2SStefano Zampini     Mat       Dd;
3723202ece2SStefano Zampini     PetscBool Ddense;
3733202ece2SStefano Zampini 
3743202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
3753202ece2SStefano Zampini     if (!Ddense) {
3763202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
3773202ece2SStefano Zampini     } else {
3783202ece2SStefano Zampini       Dd = D;
3793202ece2SStefano Zampini     }
380f11841e3SStefano Zampini     if (n_I) {
3813202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
382f11841e3SStefano Zampini     } else {
383f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
384f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
385f11841e3SStefano Zampini       } else {
386f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
387f11841e3SStefano Zampini       }
388f11841e3SStefano Zampini     }
3893202ece2SStefano Zampini     if (!Ddense) {
3903202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
3913202ece2SStefano Zampini     }
3923202ece2SStefano Zampini   } else {
3933202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
3943202ece2SStefano Zampini   }
3953202ece2SStefano Zampini   PetscFunctionReturn(0);
3963202ece2SStefano Zampini }
39734a97f8cSStefano Zampini 
39834a97f8cSStefano Zampini #undef __FUNCT__
3991580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
40091af6908SStefano 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)
401b1b3d7a2SStefano Zampini {
402be83ff47SStefano Zampini   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
403be83ff47SStefano Zampini   Mat                    S_all;
4045a94c880SStefano Zampini   Mat                    global_schur_subsets,work_mat,*submats;
4055db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
406b7ab4a40SStefano Zampini   IS                     is_I,is_I_layer;
407dc456d91SStefano Zampini   IS                     all_subsets,all_subsets_mult,all_subsets_n;
408dc456d91SStefano Zampini   PetscInt               *nnz,*all_local_idx_N;
409dc456d91SStefano Zampini   PetscInt               *auxnum1,*auxnum2;
4105a95e1ceSStefano Zampini   PetscInt               i,subset_size,max_subset_size;
411683d3df6SStefano Zampini   PetscInt               n_B,extra,local_size,global_size;
41208122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
41306a4b1faSStefano Zampini   PetscScalar            *Bwork;
4145a95e1ceSStefano Zampini   PetscSubcomm           subcomm;
4155a95e1ceSStefano Zampini   PetscMPIInt            color,rank;
4165a95e1ceSStefano Zampini   MPI_Comm               comm_n;
4179d54b7f4SStefano Zampini   PetscBool              deluxe = PETSC_TRUE, use_getr = PETSC_FALSE;
418b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
419b1b3d7a2SStefano Zampini 
420b1b3d7a2SStefano Zampini   PetscFunctionBegin;
421a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
422a64f4aa4SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
423a64f4aa4SStefano Zampini 
4243301b35fSStefano Zampini   sub_schurs->is_hermitian = PETSC_TRUE;
4253301b35fSStefano Zampini   sub_schurs->is_posdef    = PETSC_TRUE;
426*e62b6521Sstefano_zampini   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
427*e62b6521Sstefano_zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
428*e62b6521Sstefano_zampini   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
429*e62b6521Sstefano_zampini   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
430*e62b6521Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
431*e62b6521Sstefano_zampini 
432*e62b6521Sstefano_zampini   /* convert matrix if needed */
433*e62b6521Sstefano_zampini   if (Ain) {
434*e62b6521Sstefano_zampini     PetscBool isseqaij;
435*e62b6521Sstefano_zampini     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
436a64f4aa4SStefano Zampini     if (isseqaij) {
437a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
438a64f4aa4SStefano Zampini       sub_schurs->A = Ain;
4393301b35fSStefano Zampini     } else {
440a64f4aa4SStefano Zampini       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
441a64f4aa4SStefano Zampini     }
442a64f4aa4SStefano Zampini   }
4433301b35fSStefano Zampini 
444a64f4aa4SStefano Zampini   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
445a64f4aa4SStefano Zampini   sub_schurs->S = Sin;
446df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit) {
447df4d28bfSStefano Zampini     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
448a64f4aa4SStefano Zampini   }
449a64f4aa4SStefano Zampini 
4505a95e1ceSStefano Zampini   /* preliminary checks */
451af25d912SStefano 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");
4525a95e1ceSStefano Zampini 
4535a95e1ceSStefano Zampini   /* restrict work on active processes */
4545a95e1ceSStefano Zampini   color = 0;
4555a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
4565a95e1ceSStefano Zampini   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
4575a95e1ceSStefano Zampini   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
4585a95e1ceSStefano Zampini   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
4595a95e1ceSStefano Zampini   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
4605a95e1ceSStefano Zampini   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
4615a95e1ceSStefano Zampini   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
4625a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4635a95e1ceSStefano Zampini     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
4645a95e1ceSStefano Zampini     PetscFunctionReturn(0);
4655a95e1ceSStefano Zampini   }
46681ea8064SStefano Zampini   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
4675a95e1ceSStefano Zampini 
468b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
469df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) {
470a64f4aa4SStefano Zampini     Mat       tA_IB,tA_BI,tA_BB;
4713301b35fSStefano Zampini     PetscBool isseqsbaij;
472a64f4aa4SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
4733301b35fSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
4743301b35fSStefano Zampini     if (isseqsbaij) {
475a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
476a64f4aa4SStefano Zampini       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
477a64f4aa4SStefano Zampini       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
478a64f4aa4SStefano Zampini     } else {
479a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
480a64f4aa4SStefano Zampini       A_BB = tA_BB;
481a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
482a64f4aa4SStefano Zampini       A_IB = tA_IB;
483a64f4aa4SStefano Zampini       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
484a64f4aa4SStefano Zampini       A_BI = tA_BI;
485f11841e3SStefano Zampini     }
486a58a30b4SStefano Zampini   } else {
4875a95e1ceSStefano Zampini     A_II = NULL;
4885a95e1ceSStefano Zampini     A_IB = NULL;
4895a95e1ceSStefano Zampini     A_BI = NULL;
4905a95e1ceSStefano Zampini     A_BB = NULL;
491b1b3d7a2SStefano Zampini   }
4925a95e1ceSStefano Zampini   S_all = NULL;
493b1b3d7a2SStefano Zampini 
494b1b3d7a2SStefano Zampini   /* determine interior problems */
4953dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
4963dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
497b1b3d7a2SStefano Zampini     PetscBT                touched;
498b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
499b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
500b1b3d7a2SStefano Zampini 
5011633d1f0SStefano Zampini     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
502b1b3d7a2SStefano Zampini     /* get sizes */
503b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
504b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
505b1b3d7a2SStefano Zampini 
506b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
507b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
508b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
509b1b3d7a2SStefano Zampini 
510b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
511b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
512b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
513b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
514b1b3d7a2SStefano Zampini     }
515b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
516b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
517b1b3d7a2SStefano Zampini 
518b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
519b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
520b1b3d7a2SStefano Zampini     n_prev_added = n_B;
521b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
522b1b3d7a2SStefano Zampini       PetscInt n_added;
523b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5246c4ed002SBarry 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);
525b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
526b1b3d7a2SStefano Zampini       n_prev_added = n_added;
527b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
528b1b3d7a2SStefano Zampini       if (!n_added) break;
529b1b3d7a2SStefano Zampini     }
530b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
531b1b3d7a2SStefano Zampini 
532883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
533a9b99552SStefano 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);
534b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
535a9b99552SStefano Zampini     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
536883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
537df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
538b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
539b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
540a9b99552SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
541b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
542b1b3d7a2SStefano Zampini 
543b1b3d7a2SStefano Zampini       /* II block */
544b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
545b1b3d7a2SStefano Zampini     }
546b1b3d7a2SStefano Zampini   } else {
547b1b3d7a2SStefano Zampini     PetscInt n_I;
548b1b3d7a2SStefano Zampini 
549b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
550b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
551a9b99552SStefano Zampini     is_I_layer = sub_schurs->is_I;
552b1b3d7a2SStefano Zampini 
553b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
554df4d28bfSStefano Zampini     if (!sub_schurs->schur_explicit) {
555b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
556b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
557b1b3d7a2SStefano Zampini 
558b1b3d7a2SStefano Zampini       /* II block is the same */
559b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
560b1b3d7a2SStefano Zampini       AE_II = A_II;
561b1b3d7a2SStefano Zampini     }
562b1b3d7a2SStefano Zampini   }
5635a95e1ceSStefano Zampini 
564883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
5655a95e1ceSStefano Zampini   max_subset_size = 0;
566883469d8SStefano Zampini   local_size = 0;
5675a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
5685a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5695a95e1ceSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
570883469d8SStefano Zampini     local_size += subset_size;
571883469d8SStefano Zampini   }
572883469d8SStefano Zampini 
573883469d8SStefano Zampini   /* Work arrays for local indices */
574883469d8SStefano Zampini   extra = 0;
575683d3df6SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
576df4d28bfSStefano Zampini   if (sub_schurs->schur_explicit && is_I_layer) {
577a9b99552SStefano Zampini     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
578883469d8SStefano Zampini   }
579683d3df6SStefano Zampini   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
580883469d8SStefano Zampini   if (extra) {
581883469d8SStefano Zampini     const PetscInt *idxs;
582a9b99552SStefano Zampini     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
583883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
584a9b99552SStefano Zampini     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
585883469d8SStefano Zampini   }
586883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
587dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
588dc456d91SStefano Zampini   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
589883469d8SStefano Zampini 
590883469d8SStefano Zampini   /* Get local indices in local numbering */
591883469d8SStefano Zampini   local_size = 0;
5925a95e1ceSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
593883469d8SStefano Zampini     PetscInt j;
594883469d8SStefano Zampini     const    PetscInt *idxs;
595883469d8SStefano Zampini 
5965a95e1ceSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
5975a95e1ceSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
598eb595f79SStefano Zampini     /* start (smallest in global ordering) and multiplicity */
599eb595f79SStefano Zampini     auxnum1[i] = idxs[0];
600eb595f79SStefano Zampini     auxnum2[i] = subset_size;
601883469d8SStefano Zampini     /* subset indices in local numbering */
602883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
6035a95e1ceSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
604883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
605883469d8SStefano Zampini     local_size += subset_size;
606883469d8SStefano Zampini   }
607883469d8SStefano Zampini 
6085a95e1ceSStefano Zampini   /* allocate extra workspace needed only for GETRI */
60906a4e24aSStefano Zampini   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
61059ac4de7SStefano Zampini     PetscScalar  lwork,dummyscalar = 0.;
61159ac4de7SStefano Zampini     PetscBLASInt dummyint = 0;
612d2627357SStefano Zampini 
61306a4e24aSStefano Zampini     use_getr = PETSC_TRUE;
614d2627357SStefano Zampini     B_lwork = -1;
615d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
616d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
61759ac4de7SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
618d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
619d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
620d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
621d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
622056290a2SStefano Zampini   } else {
623056290a2SStefano Zampini     Bwork = NULL;
624056290a2SStefano Zampini     pivots = NULL;
625d2627357SStefano Zampini   }
626d2627357SStefano Zampini 
627d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
628dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
629dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
630dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
631dc456d91SStefano Zampini   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
6326583bcc1SStefano Zampini   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
633dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
634dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
635dc456d91SStefano Zampini   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
6366c4ed002SBarry Smith   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
637dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
638e176bc59SStefano Zampini   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
639d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
640d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
641d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
642d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6432972d61bSStefano Zampini 
6445a95e1ceSStefano Zampini   /* subset indices in local boundary numbering */
6455a95e1ceSStefano Zampini   if (!sub_schurs->is_Ej_all) {
6465a95e1ceSStefano Zampini     PetscInt *all_local_idx_B;
6475a95e1ceSStefano Zampini 
6485a95e1ceSStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
6495a95e1ceSStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
6506c4ed002SBarry 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);
6515a95e1ceSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
652b1b3d7a2SStefano Zampini   }
653b1b3d7a2SStefano Zampini 
65472b8c272SStefano Zampini   if (change) {
65572b8c272SStefano Zampini     ISLocalToGlobalMapping BtoS;
65672b8c272SStefano Zampini     IS                     change_primal_B;
65772b8c272SStefano Zampini     IS                     change_primal_all;
65872b8c272SStefano Zampini 
659b7ab4a40SStefano Zampini     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
660b7ab4a40SStefano Zampini     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
661b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
66272b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
66372b8c272SStefano Zampini       ISLocalToGlobalMapping NtoS;
66472b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
665b7ab4a40SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
66672b8c272SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
66772b8c272SStefano Zampini     }
66872b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
66972b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
67072b8c272SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
67172b8c272SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
67272b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
673b7ab4a40SStefano Zampini     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
67472b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
67572b8c272SStefano Zampini       Mat change_sub;
67672b8c272SStefano Zampini 
677b7ab4a40SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
67872b8c272SStefano Zampini       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
67972b8c272SStefano Zampini       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
68072b8c272SStefano Zampini       if (!sub_schurs->change_with_qr) {
68172b8c272SStefano Zampini         ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
68272b8c272SStefano Zampini       } else {
68372b8c272SStefano Zampini         Mat change_subt;
68472b8c272SStefano Zampini         ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
68572b8c272SStefano Zampini         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
68672b8c272SStefano Zampini         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
68772b8c272SStefano Zampini       }
68872b8c272SStefano Zampini       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
68972b8c272SStefano Zampini       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
690*e62b6521Sstefano_zampini       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
691*e62b6521Sstefano_zampini       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
69272b8c272SStefano Zampini     }
69372b8c272SStefano Zampini     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
69472b8c272SStefano Zampini   }
69572b8c272SStefano Zampini 
6965a95e1ceSStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
6975a95e1ceSStefano Zampini   if (!sub_schurs->S_Ej_all) {
6985a95e1ceSStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
6995a95e1ceSStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
7005a95e1ceSStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
7015a95e1ceSStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
702aa83b6aeSStefano Zampini   }
703b1b3d7a2SStefano Zampini 
7045a95e1ceSStefano Zampini   /* Compute Schur complements explicitly */
705be83ff47SStefano Zampini   F = NULL;
706df4d28bfSStefano 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 */
7075a95e1ceSStefano Zampini     Mat         S_Ej_expl;
7085a95e1ceSStefano Zampini     PetscScalar *work;
7095a95e1ceSStefano Zampini     PetscInt    j,*dummy_idx;
7105a95e1ceSStefano Zampini     PetscBool   Sdense;
7115a95e1ceSStefano Zampini 
7125a95e1ceSStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
7135a95e1ceSStefano Zampini     local_size = 0;
714b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
7155a95e1ceSStefano Zampini       IS  is_subset_B;
7165a95e1ceSStefano Zampini       Mat AE_EE,AE_IE,AE_EI,S_Ej;
7175a95e1ceSStefano Zampini 
7185a95e1ceSStefano Zampini       /* subsets in original and boundary numbering */
7195a95e1ceSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
7205a95e1ceSStefano Zampini       /* EE block */
7215a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
7225a95e1ceSStefano Zampini       /* IE block */
7235a95e1ceSStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
7245a95e1ceSStefano Zampini       /* EI block */
7255a95e1ceSStefano Zampini       if (sub_schurs->is_hermitian) {
7265a95e1ceSStefano Zampini         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
7275a95e1ceSStefano Zampini       } else {
7285a95e1ceSStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
7295a95e1ceSStefano Zampini       }
730a64f4aa4SStefano Zampini       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
7315a95e1ceSStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
7325a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
7335a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
7345a95e1ceSStefano Zampini       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
735b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
736b1b3d7a2SStefano Zampini         KSP ksp;
737b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
7385a95e1ceSStefano Zampini         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
739b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
740b1b3d7a2SStefano Zampini         KSP       origksp,schurksp;
741b1b3d7a2SStefano Zampini         PC        origpc,schurpc;
742b1b3d7a2SStefano Zampini         KSPType   ksp_type;
743b1b3d7a2SStefano Zampini         PetscInt  n_internal;
7445a95e1ceSStefano Zampini         PetscBool ispcnone;
745b1b3d7a2SStefano Zampini 
746b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
7475a95e1ceSStefano Zampini         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
748b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
749b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
750b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
751b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
7525a95e1ceSStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
7535a95e1ceSStefano Zampini         if (!ispcnone) {
7545a95e1ceSStefano Zampini           PCType pc_type;
755b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
756b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
7575a95e1ceSStefano Zampini         } else {
7585a95e1ceSStefano Zampini           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
7595a95e1ceSStefano Zampini         }
760b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
761b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
762b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
763b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
764b1b3d7a2SStefano Zampini           if (solver) {
765b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
766b1b3d7a2SStefano Zampini           }
767b1b3d7a2SStefano Zampini         }
768b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
769b1b3d7a2SStefano Zampini       }
7705a95e1ceSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
7715a95e1ceSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
7725a95e1ceSStefano Zampini       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
7735a95e1ceSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
7745a95e1ceSStefano Zampini       if (Sdense) {
7755a95e1ceSStefano Zampini         for (j=0;j<subset_size;j++) {
7765a95e1ceSStefano Zampini           dummy_idx[j]=local_size+j;
777b1b3d7a2SStefano Zampini         }
7785a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
7796c4ed002SBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
7805a95e1ceSStefano Zampini       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
781a64f4aa4SStefano Zampini       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
7825a95e1ceSStefano Zampini       local_size += subset_size;
7835a95e1ceSStefano Zampini     }
7845a95e1ceSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
785b1b3d7a2SStefano Zampini     /* free */
786b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
787b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
7885a95e1ceSStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
789883469d8SStefano Zampini   } else {
7905cbda25cSStefano Zampini     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
7919d54b7f4SStefano Zampini     Vec         Dall = NULL;
792ca92afb2SStefano Zampini     IS          is_A_all,*is_p_r = NULL;
7935cbda25cSStefano Zampini     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
7946dba178dSStefano Zampini     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
795df4d28bfSStefano Zampini     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
796d4933d67SStefano Zampini     char        solver[256];
797883469d8SStefano Zampini 
798683d3df6SStefano Zampini     /* get sizes */
79981ea8064SStefano Zampini     n_I = 0;
80081ea8064SStefano Zampini     if (is_I_layer) {
801a9b99552SStefano Zampini       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
80281ea8064SStefano Zampini     }
803683d3df6SStefano Zampini     economic = PETSC_FALSE;
804683d3df6SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
805683d3df6SStefano Zampini     if (cum != n_I) economic = PETSC_TRUE;
806683d3df6SStefano Zampini     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
8079d54b7f4SStefano Zampini     size_active_schur = local_size;
8089d54b7f4SStefano Zampini 
809f17d2ae1SStefano Zampini     /* import scaling vector (wrong formulation if we have 3D edges) */
8109d54b7f4SStefano Zampini     if (scaling && compute_Stilda) {
8119d54b7f4SStefano Zampini       const PetscScalar *array;
8129d54b7f4SStefano Zampini       PetscScalar       *array2;
8139d54b7f4SStefano Zampini       const PetscInt    *idxs;
8149d54b7f4SStefano Zampini       PetscInt          i;
8159d54b7f4SStefano Zampini 
8169d54b7f4SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8179d54b7f4SStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
8189d54b7f4SStefano Zampini       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
8199d54b7f4SStefano Zampini       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
8209d54b7f4SStefano Zampini       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
8219d54b7f4SStefano Zampini       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
8229d54b7f4SStefano Zampini       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
8239d54b7f4SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
8249d54b7f4SStefano Zampini       deluxe = PETSC_FALSE;
8259d54b7f4SStefano Zampini     }
826d62866d3SStefano Zampini 
827683d3df6SStefano Zampini     /* size active schurs does not count any dirichlet or vertex dof on the interface */
828683d3df6SStefano Zampini     cum = n_I+size_active_schur;
829683d3df6SStefano Zampini     if (sub_schurs->is_dir) {
830683d3df6SStefano Zampini       const PetscInt* idxs;
831683d3df6SStefano Zampini       PetscInt        n_dir;
832683d3df6SStefano Zampini 
833683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
834683d3df6SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
835683d3df6SStefano Zampini       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
836683d3df6SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
837683d3df6SStefano Zampini       cum += n_dir;
838d62866d3SStefano Zampini     }
839683d3df6SStefano Zampini     factor_workaround = PETSC_FALSE;
840683d3df6SStefano Zampini     /* include the primal vertices in the Schur complement */
841367aa537SStefano Zampini     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
842683d3df6SStefano Zampini       PetscInt n_v;
843683d3df6SStefano Zampini 
844683d3df6SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
845683d3df6SStefano Zampini       if (n_v) {
846683d3df6SStefano Zampini         const PetscInt* idxs;
847683d3df6SStefano Zampini 
848683d3df6SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
849683d3df6SStefano Zampini         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
850683d3df6SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
851683d3df6SStefano Zampini         cum += n_v;
852683d3df6SStefano Zampini         factor_workaround = PETSC_TRUE;
853683d3df6SStefano Zampini       }
854683d3df6SStefano Zampini     }
855683d3df6SStefano Zampini     size_schur = cum - n_I;
856683d3df6SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
857683d3df6SStefano Zampini     /* get working mat (TODO: factorize without actually permuting)  */
858683d3df6SStefano Zampini     if (cum == n) {
859683d3df6SStefano Zampini       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
860683d3df6SStefano Zampini       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
861683d3df6SStefano Zampini     } else {
8626816873aSStefano Zampini       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
863683d3df6SStefano Zampini     }
864*e62b6521Sstefano_zampini     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
865*e62b6521Sstefano_zampini     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
866ca92afb2SStefano Zampini 
867ca92afb2SStefano Zampini     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
868367aa537SStefano Zampini        n^2 instead of n^1.5 or something. This is a workaround */
869ca92afb2SStefano Zampini     if (benign_n) {
870ca92afb2SStefano Zampini       Vec                    v;
871ca92afb2SStefano Zampini       ISLocalToGlobalMapping N_to_reor;
872ca92afb2SStefano Zampini       IS                     is_p0,is_p0_p;
8735cbda25cSStefano Zampini       PetscScalar            *cs_AIB,*AIIm1_data;
8745cbda25cSStefano Zampini       PetscInt               sizeA;
875ca92afb2SStefano Zampini 
876ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
877ca92afb2SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
878ca92afb2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
879ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
880ca92afb2SStefano Zampini       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
8815cbda25cSStefano Zampini       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
8825cbda25cSStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
883282d6408SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
8845cbda25cSStefano Zampini       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
8855cbda25cSStefano Zampini       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
886ca92afb2SStefano Zampini       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
887ca92afb2SStefano Zampini       /* compute colsum of A_IB restricted to pressures */
888ca92afb2SStefano Zampini       for (i=0;i<benign_n;i++) {
8895cbda25cSStefano Zampini         Vec            benign_AIIm1_ones;
890ca92afb2SStefano Zampini         PetscScalar    *array;
891ca92afb2SStefano Zampini         const PetscInt *idxs;
892ca92afb2SStefano Zampini         PetscInt       j,nz;
893ca92afb2SStefano Zampini 
894ca92afb2SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
895ca92afb2SStefano Zampini         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
896ca92afb2SStefano Zampini         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8975cbda25cSStefano Zampini         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
898ca92afb2SStefano Zampini         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
8995cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
9005cbda25cSStefano Zampini         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
9015cbda25cSStefano Zampini         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
902ca92afb2SStefano Zampini         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
90322db5ddcSStefano Zampini         for (j=0;j<size_schur;j++) {
90422db5ddcSStefano Zampini #if defined(PETSC_USE_COMPLEX)
90522db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
90622db5ddcSStefano Zampini #else
90722db5ddcSStefano Zampini           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
90822db5ddcSStefano Zampini #endif
90922db5ddcSStefano Zampini         }
910ca92afb2SStefano Zampini         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
911ca92afb2SStefano Zampini       }
9125cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9135cbda25cSStefano Zampini       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
914ca92afb2SStefano Zampini       ierr = VecDestroy(&v);CHKERRQ(ierr);
915ca92afb2SStefano Zampini       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
916ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
917ec968b6dSstefano_zampini       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
918ca92afb2SStefano Zampini       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
919ca92afb2SStefano Zampini       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
920ca92afb2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
921ca92afb2SStefano Zampini     }
9223301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
9233301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
9243301b35fSStefano Zampini     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
925df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
926d4933d67SStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
927df4d28bfSStefano Zampini #else
928df4d28bfSStefano Zampini     ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr);
929df4d28bfSStefano Zampini #endif
930*e62b6521Sstefano_zampini     ierr = PetscOptionsGetString(NULL,((PetscObject)A)->prefix,"-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
931883469d8SStefano Zampini 
932683d3df6SStefano Zampini     /* when using the benign subspace trick, the local Schur complements are SPD */
933d47842beSStefano Zampini     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
934d47842beSStefano Zampini 
935683d3df6SStefano Zampini     if (n_I) { /* TODO add ordering from options */
9365a05ddb0SStefano Zampini       IS is_schur;
9375a05ddb0SStefano Zampini 
9389ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
939d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
940883469d8SStefano Zampini       } else {
941d4933d67SStefano Zampini         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
942883469d8SStefano Zampini       }
9438c09ecd8Sstefano_zampini       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
9448c09ecd8Sstefano_zampini 
945883469d8SStefano Zampini       /* subsets ordered last */
9465a05ddb0SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
9475a05ddb0SStefano Zampini       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
9485a05ddb0SStefano Zampini       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
949883469d8SStefano Zampini 
950883469d8SStefano Zampini       /* factorization step */
9519ab7bb16SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
952883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
953be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
954be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
955be83ff47SStefano Zampini #endif
956e8ade678SStefano Zampini         if (sub_schurs->is_posdef) {
957e8ade678SStefano Zampini           ierr = MatFactorSetSchurComplementSolverType(F,1);CHKERRQ(ierr);
958e8ade678SStefano Zampini         } else {
959e8ade678SStefano Zampini           ierr = MatFactorSetSchurComplementSolverType(F,2);CHKERRQ(ierr);
960e8ade678SStefano Zampini         }
961883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
962a0b0af32SStefano Zampini         S_lower_triangular = PETSC_TRUE;
963883469d8SStefano Zampini       } else {
964883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
965be83ff47SStefano Zampini #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
966be83ff47SStefano Zampini         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
967be83ff47SStefano Zampini #endif
968883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
969883469d8SStefano Zampini       }
970ec968b6dSstefano_zampini       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
971883469d8SStefano Zampini 
972883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
9735a05ddb0SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
974d62866d3SStefano Zampini       /* we can reuse the solvers if we are not using the economic version */
975683d3df6SStefano Zampini       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
976683d3df6SStefano Zampini       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
977df4d28bfSStefano Zampini       solver_S = PETSC_TRUE;
978ca92afb2SStefano Zampini 
97972b8c272SStefano Zampini       /* update the Schur complement with the change of basis on the pressures */
980ca92afb2SStefano Zampini       if (benign_n) {
9815cbda25cSStefano Zampini         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
982ca92afb2SStefano Zampini         Vec         v;
9835cbda25cSStefano Zampini         PetscInt    sizeA;
984ca92afb2SStefano Zampini 
985ca92afb2SStefano Zampini         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
9865cbda25cSStefano Zampini         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
9875cbda25cSStefano Zampini         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
988ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
989ca92afb2SStefano Zampini         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
990ca92afb2SStefano Zampini #endif
991ca92afb2SStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
992ca92afb2SStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
993ca92afb2SStefano Zampini #endif
9945cbda25cSStefano Zampini         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
9955cbda25cSStefano Zampini         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
996ca92afb2SStefano Zampini         for (i=0;i<benign_n;i++) {
9975cbda25cSStefano Zampini           Vec            benign_AIIm1_ones;
99847484b83SStefano Zampini           PetscScalar    *array,sum = 0.,one = 1.;
999ca92afb2SStefano Zampini           const PetscInt *idxs;
1000ca92afb2SStefano Zampini           PetscInt       j,nz;
100147484b83SStefano Zampini           PetscBLASInt   B_k,B_n;
1002ca92afb2SStefano Zampini 
10035cbda25cSStefano Zampini           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
10045cbda25cSStefano Zampini           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
10055cbda25cSStefano Zampini           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1006ca92afb2SStefano Zampini           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1007ca92afb2SStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1008282d6408SStefano Zampini           /* p0 dof (eliminated) is excluded from the sum */
1009282d6408SStefano Zampini           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1010ca92afb2SStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10115cbda25cSStefano Zampini           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1012ca92afb2SStefano Zampini           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
101347484b83SStefano Zampini           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
101447484b83SStefano Zampini           /* cs_AIB already scaled by 1./nz */
101547484b83SStefano Zampini           B_k = 1;
101647484b83SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
101747484b83SStefano Zampini           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
10185cbda25cSStefano Zampini           sum  = 1.;
101947484b83SStefano 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));
1020ca92afb2SStefano Zampini           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
10215cbda25cSStefano Zampini           /* set p0 entry of AIIm1_ones to zero */
10225cbda25cSStefano Zampini           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1023282d6408SStefano Zampini           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
10245cbda25cSStefano Zampini           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
10255cbda25cSStefano Zampini           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1026ca92afb2SStefano Zampini         }
10275cbda25cSStefano Zampini   /* restore defaults */
10285cbda25cSStefano Zampini #if defined(PETSC_HAVE_MUMPS)
10295cbda25cSStefano Zampini         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
10305cbda25cSStefano Zampini #endif
10315cbda25cSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
10325cbda25cSStefano Zampini         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
10335cbda25cSStefano Zampini #endif
10345cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
10355cbda25cSStefano Zampini         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1036ca92afb2SStefano Zampini         ierr = VecDestroy(&v);CHKERRQ(ierr);
1037ca92afb2SStefano Zampini         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1038ca92afb2SStefano Zampini       }
1039a3df083aSStefano Zampini       if (!reuse_solvers) {
1040a3df083aSStefano Zampini         for (i=0;i<benign_n;i++) {
1041a3df083aSStefano Zampini           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1042a3df083aSStefano Zampini         }
1043a3df083aSStefano Zampini         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
10445cbda25cSStefano Zampini         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
10455cbda25cSStefano Zampini         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1046a3df083aSStefano Zampini       }
1047df4d28bfSStefano Zampini     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1048be83ff47SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1049683d3df6SStefano Zampini       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1050166598c1SStefano Zampini       factor_workaround = PETSC_FALSE;
1051df4d28bfSStefano Zampini       solver_S = PETSC_FALSE;
1052be83ff47SStefano Zampini     }
1053be83ff47SStefano Zampini 
1054be83ff47SStefano Zampini     if (reuse_solvers) {
1055a00504b5SStefano Zampini       Mat                A_II,Afake;
105653892102SStefano Zampini       Vec                vec1_B;
1057df4d28bfSStefano Zampini       PCBDDCReuseSolvers msolv_ctx;
10583462e049SStefano Zampini       PetscInt           n_R;
1059d5574798SStefano Zampini 
1060df4d28bfSStefano Zampini       if (sub_schurs->reuse_solver) {
1061df4d28bfSStefano Zampini         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1062e28d306cSStefano Zampini       } else {
1063df4d28bfSStefano Zampini         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1064d62866d3SStefano Zampini       }
1065df4d28bfSStefano Zampini       msolv_ctx = sub_schurs->reuse_solver;
1066be83ff47SStefano Zampini       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1067d5574798SStefano Zampini       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1068d5574798SStefano Zampini       msolv_ctx->F = F;
1069683d3df6SStefano Zampini       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1070683d3df6SStefano Zampini       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1071683d3df6SStefano Zampini       {
1072683d3df6SStefano Zampini         PetscScalar *array;
1073683d3df6SStefano Zampini         PetscInt    n;
1074683d3df6SStefano Zampini 
1075683d3df6SStefano Zampini         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1076683d3df6SStefano Zampini         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1077683d3df6SStefano Zampini         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1078683d3df6SStefano Zampini         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1079683d3df6SStefano Zampini       }
1080683d3df6SStefano Zampini       msolv_ctx->has_vertices = factor_workaround;
1081d62866d3SStefano Zampini 
1082d62866d3SStefano Zampini       /* interior solver */
1083683d3df6SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1084d62866d3SStefano Zampini       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1085d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1086d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1087df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1088df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1089d62866d3SStefano Zampini 
1090d62866d3SStefano Zampini       /* correction solver */
10913462e049SStefano Zampini       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1092d62866d3SStefano Zampini       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1093d62866d3SStefano Zampini       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1094df4d28bfSStefano Zampini       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1095df4d28bfSStefano Zampini       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
109653892102SStefano Zampini 
109753892102SStefano Zampini       /* scatter and vecs for Schur complement solver */
109853892102SStefano Zampini       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
109953892102SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1100683d3df6SStefano Zampini       if (!factor_workaround) {
110153892102SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
110253892102SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
110353892102SStefano Zampini         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
110453892102SStefano Zampini         msolv_ctx->is_R = is_A_all;
1105683d3df6SStefano Zampini       } else {
1106683d3df6SStefano Zampini         IS              is_B_all;
1107683d3df6SStefano Zampini         const PetscInt* idxs;
1108683d3df6SStefano Zampini         PetscInt        dual,n_v,n;
1109683d3df6SStefano Zampini 
1110683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1111683d3df6SStefano Zampini         dual = size_schur - n_v;
1112683d3df6SStefano Zampini         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1113683d3df6SStefano Zampini         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1114683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1115683d3df6SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1116683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1117683d3df6SStefano Zampini         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1118683d3df6SStefano Zampini         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1119683d3df6SStefano Zampini         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1120683d3df6SStefano Zampini         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1121683d3df6SStefano Zampini         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1122683d3df6SStefano Zampini       }
11233462e049SStefano Zampini       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
11243462e049SStefano Zampini       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
11253462e049SStefano Zampini       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
11263462e049SStefano Zampini       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1127683d3df6SStefano Zampini       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1128ca92afb2SStefano Zampini 
1129ca92afb2SStefano Zampini       /* communicate benign info to solver context */
1130ca92afb2SStefano Zampini       if (benign_n) {
11315cbda25cSStefano Zampini         PetscScalar *array;
11325cbda25cSStefano Zampini 
1133ca92afb2SStefano Zampini         msolv_ctx->benign_n = benign_n;
1134ca92afb2SStefano Zampini         msolv_ctx->benign_zerodiag_subs = is_p_r;
1135ca92afb2SStefano Zampini         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
11365cbda25cSStefano Zampini         msolv_ctx->benign_csAIB = cs_AIB_mat;
1137282d6408SStefano Zampini         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
11385cbda25cSStefano Zampini         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11395cbda25cSStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
11405cbda25cSStefano Zampini         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
11415cbda25cSStefano Zampini         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1142ca92afb2SStefano Zampini       }
1143d5574798SStefano Zampini     }
114408122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
114553892102SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
11465db18549SStefano Zampini 
1147be83ff47SStefano Zampini     /* Work arrays */
1148be83ff47SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
11495a94c880SStefano Zampini 
11505a94c880SStefano Zampini     /* matrices for adaptive selection */
115112d906b1SStefano Zampini     if (compute_Stilda) {
11525a94c880SStefano Zampini       if (!sub_schurs->sum_S_Ej_tilda_all) {
11535a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
11545a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
11555a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
11565a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
11575a94c880SStefano Zampini       }
11589d54b7f4SStefano Zampini       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
11595a94c880SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
11605a94c880SStefano Zampini         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
11615a94c880SStefano Zampini         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
11625a94c880SStefano Zampini         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
11635a94c880SStefano Zampini       }
116412d906b1SStefano Zampini     }
1165d2627357SStefano Zampini 
1166be83ff47SStefano Zampini     /* S_Ej_all */
1167be83ff47SStefano Zampini     cum = cum2 = 0;
1168be83ff47SStefano Zampini     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
116965d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
117065d8bf0aSStefano Zampini       PetscInt j;
117165d8bf0aSStefano Zampini 
11725a95e1ceSStefano Zampini       /* get S_E */
1173b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1174683d3df6SStefano Zampini       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1175be83ff47SStefano Zampini         PetscInt k;
1176be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1177be83ff47SStefano Zampini           for (j=k;j<subset_size;j++) {
1178be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1179683d3df6SStefano Zampini             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1180be83ff47SStefano Zampini           }
1181be83ff47SStefano Zampini         }
118206a4e24aSStefano Zampini       } else { /* just copy to workspace */
1183be83ff47SStefano Zampini         PetscInt k;
1184be83ff47SStefano Zampini         for (k=0;k<subset_size;k++) {
1185be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) {
1186be83ff47SStefano Zampini             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1187be83ff47SStefano Zampini           }
1188be83ff47SStefano Zampini         }
11899087bf02SStefano Zampini       }
11905a95e1ceSStefano Zampini       /* insert S_E values */
1191be83ff47SStefano Zampini       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1192b7ab4a40SStefano Zampini       if (sub_schurs->change) {
11938760537fSStefano Zampini         Mat            change_sub,SEj,T;
119472b8c272SStefano Zampini 
119572b8c272SStefano Zampini         /* change basis */
119672b8c272SStefano Zampini         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
119772b8c272SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
11988760537fSStefano Zampini         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
11998760537fSStefano Zampini           Mat T2;
12008760537fSStefano Zampini           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
12018760537fSStefano Zampini           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
12028760537fSStefano Zampini           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
120372b8c272SStefano Zampini           ierr = MatDestroy(&T2);CHKERRQ(ierr);
12048760537fSStefano Zampini         } else {
12058760537fSStefano Zampini           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
120672b8c272SStefano Zampini         }
120772b8c272SStefano Zampini         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
120872b8c272SStefano Zampini         ierr = MatDestroy(&T);CHKERRQ(ierr);
1209afa0f562SStefano Zampini         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
121072b8c272SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
121172b8c272SStefano Zampini       }
12129d54b7f4SStefano Zampini       if (deluxe) {
12135a95e1ceSStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1214683d3df6SStefano Zampini         /* if adaptivity is requested, invert S_E blocks */
1215862806e4SStefano Zampini         if (compute_Stilda) {
121608122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
121708122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
121806a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
12196c3e6151SStefano Zampini             PetscInt k;
12205a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
12212972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
12225a95e1ceSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
12232972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
12246c3e6151SStefano Zampini             for (k=0;k<subset_size;k++) {
12256c3e6151SStefano Zampini               for (j=k;j<subset_size;j++) {
12266c3e6151SStefano Zampini                 work[j*subset_size+k] = work[k*subset_size+j];
12276c3e6151SStefano Zampini               }
12286c3e6151SStefano Zampini             }
1229d6462365SStefano Zampini           } else {
1230d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1231d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1232d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1233d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
12342972d61bSStefano Zampini           }
123508122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
12365a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12379087bf02SStefano Zampini         }
12389d54b7f4SStefano Zampini       } else if (compute_Stilda) { /* not using deluxe */
12399d54b7f4SStefano Zampini         Mat         SEj;
12409d54b7f4SStefano Zampini         Vec         D;
12419d54b7f4SStefano Zampini         PetscScalar *array;
12429d54b7f4SStefano Zampini 
12439d54b7f4SStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
12449d54b7f4SStefano Zampini         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
12459d54b7f4SStefano Zampini         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
12469d54b7f4SStefano Zampini         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1247f17d2ae1SStefano Zampini         ierr = VecShift(D,-1.);CHKERRQ(ierr);
12489d54b7f4SStefano Zampini         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
12499d54b7f4SStefano Zampini         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
12509d54b7f4SStefano Zampini         ierr = VecDestroy(&D);CHKERRQ(ierr);
12519d54b7f4SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12529d54b7f4SStefano Zampini       }
1253be83ff47SStefano Zampini       cum += subset_size;
1254be83ff47SStefano Zampini       cum2 += subset_size*(size_schur + 1);
1255be83ff47SStefano Zampini     }
1256be83ff47SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
12574a6c6b0dSStefano Zampini 
1258df4d28bfSStefano Zampini     if (solver_S) {
12595a05ddb0SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
12604a6c6b0dSStefano Zampini     }
1261683d3df6SStefano Zampini 
1262683d3df6SStefano Zampini     schur_factor = NULL;
126345951f25SStefano Zampini     if (compute_Stilda && size_active_schur) {
1264683d3df6SStefano Zampini 
12659d54b7f4SStefano Zampini       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
12664a6c6b0dSStefano Zampini         PetscInt j;
12674a6c6b0dSStefano Zampini         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
12685a94c880SStefano Zampini         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
12694a6c6b0dSStefano Zampini       } else {
1270683d3df6SStefano Zampini         Mat S_all_inv=NULL;
1271df4d28bfSStefano Zampini         if (solver_S) { /* use MatFactor calls to invert S */
127272b8c272SStefano Zampini             /* let MatFactor factor the Schur complement */
12736dba178dSStefano Zampini           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
1274683d3df6SStefano Zampini           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1275683d3df6SStefano 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 */
1276683d3df6SStefano Zampini           if (factor_workaround) {
1277683d3df6SStefano Zampini             PetscScalar *data;
1278683d3df6SStefano Zampini             PetscInt nd = 0;
12796dba178dSStefano Zampini 
1280683d3df6SStefano Zampini             if (!sub_schurs->is_posdef) {
1281683d3df6SStefano Zampini               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1282683d3df6SStefano Zampini             }
12836dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1284683d3df6SStefano Zampini             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1285683d3df6SStefano Zampini             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1286683d3df6SStefano Zampini               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1287683d3df6SStefano Zampini             }
1288683d3df6SStefano Zampini             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
1289683d3df6SStefano Zampini             cum = 0;
1290683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1291683d3df6SStefano Zampini               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1292683d3df6SStefano Zampini               cum += size_active_schur-i;
1293683d3df6SStefano Zampini             }
1294683d3df6SStefano Zampini             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
12956dba178dSStefano Zampini             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1296683d3df6SStefano Zampini             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1297683d3df6SStefano Zampini             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1298683d3df6SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
1299683d3df6SStefano Zampini             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1300683d3df6SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1301683d3df6SStefano Zampini             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1302683d3df6SStefano Zampini           } else {
13036dba178dSStefano Zampini             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
13046dba178dSStefano Zampini             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1305683d3df6SStefano Zampini           }
1306df4d28bfSStefano Zampini         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1307683d3df6SStefano Zampini           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1308683d3df6SStefano Zampini           S_all_inv = S_all;
1309683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1310be83ff47SStefano Zampini           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1311be83ff47SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
131206a4e24aSStefano Zampini           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1313be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1314be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1315be83ff47SStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1316be83ff47SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1317d6462365SStefano Zampini           } else {
1318d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1319d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1320d6462365SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1321d6462365SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1322be83ff47SStefano Zampini           }
1323be83ff47SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1324683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1325be83ff47SStefano Zampini         }
1326be83ff47SStefano Zampini         /* S_Ej_tilda_all */
1327be83ff47SStefano Zampini         cum = cum2 = 0;
1328683d3df6SStefano Zampini         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1329be83ff47SStefano Zampini         for (i=0;i<sub_schurs->n_subs;i++) {
1330be83ff47SStefano Zampini           PetscInt j;
1331862806e4SStefano Zampini 
1332862806e4SStefano Zampini           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1333be83ff47SStefano Zampini           /* get (St^-1)_E */
133472b8c272SStefano Zampini           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
133506a4e24aSStefano Zampini              will be properly accessed later during adaptive selection */
1336a0b0af32SStefano Zampini           if (S_lower_triangular) {
1337be83ff47SStefano Zampini             PetscInt k;
1338b7ab4a40SStefano Zampini             if (sub_schurs->change) {
1339be83ff47SStefano Zampini               for (k=0;k<subset_size;k++) {
1340be83ff47SStefano Zampini                 for (j=k;j<subset_size;j++) {
1341be83ff47SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
13426c3e6151SStefano Zampini                   work[j*subset_size+k] = work[k*subset_size+j];
1343be83ff47SStefano Zampini                 }
1344be83ff47SStefano Zampini               }
134572b8c272SStefano Zampini             } else {
134672b8c272SStefano Zampini               for (k=0;k<subset_size;k++) {
134772b8c272SStefano Zampini                 for (j=k;j<subset_size;j++) {
134872b8c272SStefano Zampini                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
134972b8c272SStefano Zampini                 }
135072b8c272SStefano Zampini               }
135172b8c272SStefano Zampini             }
135272b8c272SStefano Zampini           } else {
1353be83ff47SStefano Zampini             PetscInt k;
1354be83ff47SStefano Zampini             for (k=0;k<subset_size;k++) {
1355be83ff47SStefano Zampini               for (j=0;j<subset_size;j++) {
1356be83ff47SStefano Zampini                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1357be83ff47SStefano Zampini               }
1358be83ff47SStefano Zampini             }
1359be83ff47SStefano Zampini           }
1360b7ab4a40SStefano Zampini           if (sub_schurs->change) {
13618760537fSStefano Zampini             Mat change_sub,SEj,T;
136272b8c272SStefano Zampini 
136372b8c272SStefano Zampini             /* change basis */
136472b8c272SStefano Zampini             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
136572b8c272SStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
13668760537fSStefano Zampini             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
13678760537fSStefano Zampini               Mat T2;
13688760537fSStefano Zampini               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
13698760537fSStefano Zampini               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
137072b8c272SStefano Zampini               ierr = MatDestroy(&T2);CHKERRQ(ierr);
13718760537fSStefano Zampini               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
13728760537fSStefano Zampini             } else {
13738760537fSStefano Zampini               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
137472b8c272SStefano Zampini             }
13758760537fSStefano Zampini             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
13768760537fSStefano Zampini             ierr = MatDestroy(&T);CHKERRQ(ierr);
137772b8c272SStefano Zampini             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1378afa0f562SStefano Zampini             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
137972b8c272SStefano Zampini             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
138072b8c272SStefano Zampini           }
1381be83ff47SStefano Zampini           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
13825a94c880SStefano Zampini           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1383be83ff47SStefano Zampini           cum += subset_size;
1384be83ff47SStefano Zampini           cum2 += subset_size*(size_schur + 1);
1385883469d8SStefano Zampini         }
1386683d3df6SStefano Zampini         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1387df4d28bfSStefano Zampini         if (solver_S) {
13886dba178dSStefano Zampini           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
13895db18549SStefano Zampini         }
1390683d3df6SStefano Zampini         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1391683d3df6SStefano Zampini       }
1392683d3df6SStefano Zampini 
1393683d3df6SStefano Zampini       /* move back factors to Schur data into F */
1394683d3df6SStefano Zampini       if (factor_workaround) {
1395683d3df6SStefano Zampini         Mat S_tmp;
1396683d3df6SStefano Zampini         PetscInt nv,nd=0;
1397683d3df6SStefano Zampini 
1398df4d28bfSStefano Zampini         if (!solver_S) {
1399683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1400683d3df6SStefano Zampini         }
1401683d3df6SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
14026dba178dSStefano Zampini         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1403683d3df6SStefano Zampini         if (sub_schurs->is_posdef) {
1404683d3df6SStefano Zampini           PetscScalar *data;
1405683d3df6SStefano Zampini 
1406683d3df6SStefano Zampini           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1407683d3df6SStefano Zampini           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1408683d3df6SStefano Zampini 
1409683d3df6SStefano Zampini           if (S_lower_triangular) {
1410683d3df6SStefano Zampini             cum = 0;
1411683d3df6SStefano Zampini             for (i=0;i<size_active_schur;i++) {
1412683d3df6SStefano Zampini               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1413683d3df6SStefano Zampini               cum += size_active_schur-i;
1414683d3df6SStefano Zampini 	    }
1415683d3df6SStefano Zampini           } else {
1416683d3df6SStefano Zampini             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1417683d3df6SStefano Zampini           }
1418683d3df6SStefano Zampini           if (sub_schurs->is_dir) {
1419683d3df6SStefano Zampini             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1420683d3df6SStefano Zampini 	    for (i=0;i<nd;i++) {
1421683d3df6SStefano Zampini 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1422683d3df6SStefano Zampini 	    }
1423683d3df6SStefano Zampini           }
14246dba178dSStefano Zampini           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1425683d3df6SStefano Zampini              set the diagonal entry of the Schur factor to a very large value */
1426683d3df6SStefano Zampini           for (i=size_active_schur+nd;i<size_schur;i++) {
14276c3e6151SStefano Zampini             data[i*(size_schur+1)] = infty;
1428683d3df6SStefano Zampini           }
1429683d3df6SStefano Zampini           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1430683d3df6SStefano Zampini         } else {
1431683d3df6SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1432683d3df6SStefano Zampini         }
14336dba178dSStefano Zampini         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
14349087bf02SStefano Zampini       }
1435367aa537SStefano Zampini     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1436367aa537SStefano Zampini       PetscScalar *data;
1437367aa537SStefano Zampini       PetscInt    nd = 0;
1438367aa537SStefano Zampini 
1439367aa537SStefano Zampini       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1440367aa537SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1441367aa537SStefano Zampini       }
1442367aa537SStefano Zampini       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
1443367aa537SStefano Zampini       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1444367aa537SStefano Zampini       for (i=0;i<size_active_schur;i++) {
1445367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1446367aa537SStefano Zampini       }
1447367aa537SStefano Zampini       for (i=size_active_schur+nd;i<size_schur;i++) {
1448367aa537SStefano Zampini         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
14496c3e6151SStefano Zampini         data[i*(size_schur+1)] = infty;
1450367aa537SStefano Zampini       }
1451367aa537SStefano Zampini       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
1452367aa537SStefano Zampini       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
14534a6c6b0dSStefano Zampini     }
14544a6c6b0dSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1455683d3df6SStefano Zampini     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
14569d54b7f4SStefano Zampini     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
14574a6c6b0dSStefano Zampini   }
14585a94c880SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
1459be83ff47SStefano Zampini   ierr = MatDestroy(&F);CHKERRQ(ierr);
1460a9b99552SStefano Zampini   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1461a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1462a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1463a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1464a64f4aa4SStefano Zampini   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
14655db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14665db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14675a95e1ceSStefano Zampini   if (compute_Stilda) {
14685a94c880SStefano Zampini     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14695a94c880SStefano Zampini     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14709d54b7f4SStefano Zampini     if (deluxe) {
14715a94c880SStefano Zampini       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14725a94c880SStefano Zampini       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
147308122e43SStefano Zampini     }
14749d54b7f4SStefano Zampini   }
1475a1337663SStefano Zampini 
14765db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
14775db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
14783927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
14793927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
14805a95e1ceSStefano Zampini 
14815db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
14825a94c880SStefano Zampini   if (!sub_schurs->sum_S_Ej_all) {
1483dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
14845a94c880SStefano Zampini     sub_schurs->sum_S_Ej_all = submats[0];
14855a94c880SStefano Zampini   } else {
14865a94c880SStefano Zampini     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
14875a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_all;
1488dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
14895a94c880SStefano Zampini   }
149008122e43SStefano Zampini 
1491f6f667cfSStefano Zampini   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
14925a95e1ceSStefano Zampini   if (compute_Stilda) {
14935a94c880SStefano Zampini     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1494a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
14955a94c880SStefano Zampini     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1496dc456d91SStefano Zampini     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
14979d54b7f4SStefano Zampini     if (deluxe) {
14985a94c880SStefano Zampini       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
149908122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
15005a94c880SStefano Zampini       submats[0] = sub_schurs->sum_S_Ej_inv_all;
1501dc456d91SStefano Zampini       ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
15029d54b7f4SStefano Zampini     } else {
15039d54b7f4SStefano Zampini       PetscScalar *array;
15049d54b7f4SStefano Zampini       PetscInt    cum;
15059d54b7f4SStefano Zampini 
15069d54b7f4SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15079d54b7f4SStefano Zampini       cum = 0;
15089d54b7f4SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
15099d54b7f4SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
15109d54b7f4SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
15119d54b7f4SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
15129d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
15139d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
15149d54b7f4SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
15159d54b7f4SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
15169d54b7f4SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
15179d54b7f4SStefano Zampini         cum += subset_size*subset_size;
15189d54b7f4SStefano Zampini       }
15199d54b7f4SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
15209d54b7f4SStefano Zampini       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
15219d54b7f4SStefano Zampini       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
15229d54b7f4SStefano Zampini       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
15239d54b7f4SStefano Zampini     }
152408122e43SStefano Zampini   }
15253202ece2SStefano Zampini 
15265a95e1ceSStefano Zampini   /* free workspace */
15275a94c880SStefano Zampini   ierr = PetscFree(submats);CHKERRQ(ierr);
152806a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1529a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
15303202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1531dc456d91SStefano Zampini   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
15325a95e1ceSStefano Zampini   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1533b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1534b1b3d7a2SStefano Zampini }
1535b1b3d7a2SStefano Zampini 
1536b1b3d7a2SStefano Zampini #undef __FUNCT__
1537b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
15388b6046baSStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1539b1b3d7a2SStefano Zampini {
15409bb4a8caSStefano Zampini   IS              *faces,*edges,*all_cc,vertices;
15415a95e1ceSStefano Zampini   PetscInt        i,n_faces,n_edges,n_all_cc;
1542b1b3d7a2SStefano Zampini   PetscBool       is_sorted;
1543b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
1544b1b3d7a2SStefano Zampini 
1545b1b3d7a2SStefano Zampini   PetscFunctionBegin;
1546b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
15476c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1548b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
15496c4ed002SBarry Smith   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1550b1b3d7a2SStefano Zampini 
1551b1b3d7a2SStefano Zampini   /* reset any previous data */
1552b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1553b1b3d7a2SStefano Zampini 
15545a95e1ceSStefano Zampini   /* get index sets for faces and edges (already sorted by global ordering) */
15559bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1556b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
155708122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1558b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1559b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
15608b6046baSStefano Zampini     if (copycc) {
15618b6046baSStefano Zampini       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
15628b6046baSStefano Zampini     } else {
1563c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1564b1b3d7a2SStefano Zampini       all_cc[i] = faces[i];
1565b1b3d7a2SStefano Zampini     }
15668b6046baSStefano Zampini   }
1567b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
15688b6046baSStefano Zampini     if (copycc) {
15698b6046baSStefano Zampini       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
15708b6046baSStefano Zampini     } else {
1571c8272957SStefano Zampini       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1572b1b3d7a2SStefano Zampini       all_cc[n_faces+i] = edges[i];
15738b6046baSStefano Zampini     }
157408122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1575b1b3d7a2SStefano Zampini   }
1576c8272957SStefano Zampini   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1577c8272957SStefano Zampini   sub_schurs->is_vertices = vertices;
1578c8272957SStefano Zampini   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1579d62866d3SStefano Zampini   sub_schurs->is_dir = NULL;
1580d62866d3SStefano Zampini   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1581b1b3d7a2SStefano Zampini 
1582df4d28bfSStefano Zampini   /* Determine if MatFactor can be used */
1583df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_FALSE;
1584883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
1585df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1586df4d28bfSStefano Zampini #endif
1587df4d28bfSStefano Zampini #if defined(PETSC_HAVE_MKL_PARDISO)
1588df4d28bfSStefano Zampini   sub_schurs->schur_explicit = PETSC_TRUE;
1589883469d8SStefano Zampini #endif
1590b1b3d7a2SStefano Zampini 
1591b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1592b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
1593b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1594b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
15955db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
15965db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
15975db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
15985db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
15995a95e1ceSStefano Zampini   sub_schurs->n_subs = n_all_cc;
1600b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
1601df4d28bfSStefano Zampini   if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1602b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1603b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1604b96c3477SStefano Zampini     }
16059bb4a8caSStefano Zampini   }
1606b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
1607b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
160808122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
1609b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
1610b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
1611b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
1612b1b3d7a2SStefano Zampini }
1613b1b3d7a2SStefano Zampini 
1614b1b3d7a2SStefano Zampini #undef __FUNCT__
161534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
161634a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
161734a97f8cSStefano Zampini {
161834a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
161934a97f8cSStefano Zampini   PetscErrorCode  ierr;
162034a97f8cSStefano Zampini 
162134a97f8cSStefano Zampini   PetscFunctionBegin;
162234a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
16235ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
162434a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
162534a97f8cSStefano Zampini   PetscFunctionReturn(0);
162634a97f8cSStefano Zampini }
162734a97f8cSStefano Zampini 
162834a97f8cSStefano Zampini #undef __FUNCT__
162934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
163034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
163134a97f8cSStefano Zampini {
163234a97f8cSStefano Zampini   PetscInt       i;
163334a97f8cSStefano Zampini   PetscErrorCode ierr;
163434a97f8cSStefano Zampini 
163534a97f8cSStefano Zampini   PetscFunctionBegin;
1636aea80f77Sstefano_zampini   if (!sub_schurs) PetscFunctionReturn(0);
16371e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1638b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1639b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1640b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
16415db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
16425db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
164341c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
164441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
164508122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1646a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
16475db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1648d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1649d62866d3SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
165008122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
165134a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
1652b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
165334a97f8cSStefano Zampini   }
16545ff63025SStefano Zampini   if (sub_schurs->n_subs) {
1655b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
16563dc780c3SStefano Zampini   }
1657df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
1658df4d28bfSStefano Zampini     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1659d62866d3SStefano Zampini   }
1660df4d28bfSStefano Zampini   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
166172b8c272SStefano Zampini   if (sub_schurs->change) {
166272b8c272SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
166372b8c272SStefano Zampini       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1664b7ab4a40SStefano Zampini       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
166572b8c272SStefano Zampini     }
166672b8c272SStefano Zampini   }
166772b8c272SStefano Zampini   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1668b7ab4a40SStefano Zampini   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
166934a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
167034a97f8cSStefano Zampini   PetscFunctionReturn(0);
167134a97f8cSStefano Zampini }
167234a97f8cSStefano Zampini 
167334a97f8cSStefano Zampini #undef __FUNCT__
1674aea80f77Sstefano_zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
1675aea80f77Sstefano_zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1676aea80f77Sstefano_zampini {
1677aea80f77Sstefano_zampini   PetscErrorCode ierr;
1678aea80f77Sstefano_zampini 
1679aea80f77Sstefano_zampini   PetscFunctionBegin;
1680aea80f77Sstefano_zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1681aea80f77Sstefano_zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1682aea80f77Sstefano_zampini   PetscFunctionReturn(0);
1683aea80f77Sstefano_zampini }
1684aea80f77Sstefano_zampini 
1685aea80f77Sstefano_zampini #undef __FUNCT__
168634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
16872a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
168834a97f8cSStefano Zampini {
168934a97f8cSStefano Zampini   PetscInt       i,j,n;
169034a97f8cSStefano Zampini   PetscErrorCode ierr;
169134a97f8cSStefano Zampini 
169234a97f8cSStefano Zampini   PetscFunctionBegin;
169334a97f8cSStefano Zampini   n = 0;
169434a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
169534a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
169634a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
169734a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
169834a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
169934a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
170034a97f8cSStefano Zampini         queue_tip[n] = dof;
170134a97f8cSStefano Zampini         n++;
170234a97f8cSStefano Zampini       }
170334a97f8cSStefano Zampini     }
170434a97f8cSStefano Zampini   }
170534a97f8cSStefano Zampini   *n_added = n;
170634a97f8cSStefano Zampini   PetscFunctionReturn(0);
170734a97f8cSStefano Zampini }
1708